Abstract interface to a mesh. More...
Public Member Functions | |
Constructors/Destructor | |
| MeshBase (int dim, const MPIComm &comm, const MeshEntityOrder &meshOrder) | |
| | |
| virtual | ~MeshBase () |
| | |
Topological Information | |
| const MeshEntityOrder & | meshOrder () const |
| Get the ordering convention used by this mesh. | |
| int | spatialDim () const |
| Get the spatial dimension of the mesh. | |
| virtual int | numCells (int cellDim) const =0 |
| Get the number of local cells having dimension cellDim. | |
| virtual int | numFacets (int cellDim, int cellLID, int facetDim) const =0 |
| Return the number of facets of the given relative cell. | |
| virtual int | facetLID (int cellDim, int cellLID, int facetDim, int facetIndex, int &facetOrientation) const =0 |
| Return the LID of a single facet cell (and optionally its orientation) with respect to a relative parent cell. | |
| virtual void | getFacetLIDs (int cellDim, const Array< int > &cellLIDs, int facetDim, Array< int > &facetLIDs, Array< int > &facetOrientations) const =0 |
| Return an array containing the LIDs of facets of dimension facetDim for the given batch of relative parent cells. | |
| void | getFacetArray (int cellDim, int cellLID, int facetDim, Array< int > &facetLIDs, Array< int > &facetOrientations) const |
| Return an array containing the LIDs of the facets of dimension facetDim of the single relative parent cell (optionally also returning the facet orientations). | |
| virtual const int * | elemZeroFacetView (int maxCellLID) const |
| Return a view of an array of LIDs for a maximum-dimensional cell's zero-dimensional facets (i.e. vertexes). | |
| virtual int | numMaxCofacets (int cellDim, int cellLID) const =0 |
| Return the number of maximal co-facets of the given cell. | |
| virtual int | maxCofacetLID (int cellDim, int cellLID, int cofacetIndex, int &facetIndex) const =0 |
| Return the LID of a maximal co-facet of a cell. | |
| virtual void | getMaxCofacetLIDs (const Array< int > &cellLIDs, MaximalCofacetBatch &cofacets) const |
| virtual void | getCofacets (int cellDim, int cellLID, int cofacetDim, Array< int > &cofacetLIDs) const =0 |
| Return an array of the LIDs of all of the co-facets for a given relative cell. | |
| virtual CellType | cellType (int cellDim) const =0 |
| Get the cell type of the given cell dimension. | |
| virtual int | label (int cellDim, int cellLID) const =0 |
| Get the label of the given cell. | |
| void | getLabels (int cellDim, const Array< int > &cellLIDs, Array< int > &labels) const |
| Get the labels for a batch of cells. | |
| virtual void | setLabel (int cellDim, int cellLID, int label)=0 |
| Set the label for the given cell. | |
| virtual Set< int > | getAllLabelsForDimension (int cellDim) const |
| virtual void | getLIDsForLabel (int cellDim, int label, Array< int > &cellLIDs) const |
Geometric Information | |
| virtual Point | nodePosition (int vertexLID) const =0 |
| Return the position of a local vertex. | |
| virtual const double * | nodePositionView (int vertexLID) const =0 |
| Return a const view into an raw array for the position of a local vertex. | |
| virtual Point | centroid (int cellDim, int cellLID) const |
| Return the centroid of a cell. | |
| virtual void | getJacobians (int cellDim, const Array< int > &cellLID, CellJacobianBatch &jBatch) const |
| Compute (or get) the Jacobians for a batch of cells. | |
| virtual void | getCellDiameters (int cellDim, const Array< int > &cellLIDs, Array< double > &diameters) const |
| Compute the diameters of a batch of cells. | |
| virtual void | outwardNormals (const Array< int > &cellLIDs, Array< Point > &outwardNormals) const |
| virtual void | tangentsToEdges (const Array< int > &cellLIDs, Array< Point > &tangenVectors) const |
| virtual void | pushForward (int cellDim, const Array< int > &cellLIDs, const Array< Point > &refPts, Array< Point > &physPts) const |
| Map points from a reference cell to physical points for a batch of cells. | |
Parallel Information | |
| const MPIComm & | comm () const |
| Return the MPI communicator over which this mesh is distributed. | |
| virtual int | ownerProcID (int cellDim, int cellLID) const =0 |
| Return the rank of the processor that owns the given cell. | |
| virtual bool | hasGID (int cellDim, int cellGID) const =0 |
| Determine whether a given cell GID exists on this processor. | |
| virtual int | mapGIDToLID (int cellDim, int cellGID) const =0 |
| Find the LID of a cell given its GID. | |
| virtual int | mapLIDToGID (int cellDim, int cellLID) const =0 |
| Find the global ID of a cell given its LID. | |
| virtual bool | hasIntermediateGIDs (int cellDim) const =0 |
| Return if cells of dimension cellDim have been assigned global IDs or not. | |
| virtual void | assignIntermediateCellGIDs (int cellDim)=0 |
| Assign global IDs to cells of dimension cellDim. | |
Reordering | |
| void | setReorderer (const CellReorderer &reorderer) |
| Set the reordering strategy to be used with this mesh. | |
| const CellReordererImplemBase * | reorderer () const |
| Get the reordering strategy to be used with this mesh. | |
Functions for Mesh with hanging nodes | |
| virtual bool | allowsHangingHodes () const |
| virtual bool | isElementHangingNode (int cellDim, int cellLID) const |
| virtual int | indexInParent (int maxCellLID) const |
| virtual int | maxChildren () const |
| virtual void | returnParentFacets (int childCellLID, int dimFacets, Array< int > &facetsLIDs, int &parentCellLIDs) const |
Store special weights in the mesh (for Adaptive Cell Integration) | |
| virtual bool | IsSpecialWeightValid () const |
| virtual void | setSpecialWeightValid (bool val) const |
| virtual void | flushSpecialWeights () const |
| virtual bool | hasSpecialWeight (int dim, int cellLID) const |
| virtual void | setSpecialWeight (int dim, int cellLID, Array< double > &w) const |
| virtual void | getSpecialWeight (int dim, int cellLID, Array< double > &w) const |
Static Public Member Functions | |
Output | |
| static bool & | staggerOutput () |
| Set whether to stagger output in parallel. Set to true for readable output in parallel debugging sessions. | |
Private Attributes | |
| int | dim_ |
| MPIComm | comm_ |
| MeshEntityOrder | order_ |
| RCP< CellReordererImplemBase > | reorderer_ |
| bool | validWeights_ |
| Array< Sundance::Map< int, Array< double > > > | specialWeights_ |
| bool | curvePoints_Are_Valid_ |
| int | nrCurvesForIntegral_ |
| Array< Sundance::Map< int, Array< Point > > > | curvePoints_ |
| Array< Sundance::Map< int, Array< Point > > > | curveDerivative_ |
| Array< Sundance::Map< int, Array< Point > > > | curveNormal_ |
| Sundance::Map< int, int > | curveID_to_ArrayIndex_ |
Store the intersection/quadrature points for the curve/surf integral <br> | |
for a curve or surf integral we need some quadrature points along the curve in one curve | |
| virtual bool | IsCurvePointsValid () const |
| virtual void | setCurvePointsValid (bool val) const |
| virtual void | flushCurvePoints () const |
| virtual bool | hasCurvePoints (int maxCellLID, int curveID) const |
| virtual void | setCurvePoints (int maxCellLID, int curveID, Array< Point > &points, Array< Point > &derivs, Array< Point > &normals) const |
| virtual void | getCurvePoints (int maxCellLID, int curveID, Array< Point > &points, Array< Point > &derivs, Array< Point > &normals) const |
| virtual int | mapCurveID_to_Index (int curveID) const |
Abstract interface to a mesh.
A mesh is a collection of connected cells (also called elements). This mesh interface is designed to allow for the representation of 1D, 2D, and 3D (and even higher dimensional) meshes. A mesh contains both topological information and geometric information about the cells and their subparts (i.e. facets). Topological information refers to information about how cells of various types and dimensions are connected to each other within the mesh. Geometric information refers to information about the position, size, and shape of cells within the mesh.
Currently, this interface only handles meshes where all cells of a given cell dimension have the same cell type. For example, this only allows meshes with all triangles or all quads in 2D or all tetrahedrals in 3D.
The following definitions are used in the specification of this mesh interface.
Topological information
Spatial Dimension (spatialDim): The maximum cell dimension of any cell type in the mesh. For example, a 3D mesh will have spatialDim=3, a 2D mesh will have spatialDim=2 and a 1D mesh will have spatialDim=1. Typically, when we use the short hand 1D, 2D or 3D, we are referring to the spatial dimension of the mesh. Note that 4D and higher-dimensional meshes are also representable by this interface in principle with some small changes (e.g. like making the size of a point more flexible)
Cell: A mesh object of some dimension. For example, in a 3D mesh, cells are things like elements, faces, edges, and nodes.
Cell Type (cellType): An enumeration of type CellType that lists some common cell types. For example, common cell types are things like vertexes (0D), lines or edges (1D), triangles (2D), and quadrilaterals (2D).
Cell Dimension (cellDim): The dimension of a cell. More formally, the cell dimension is equal to the number of lower level facet types in the cell. For example, A 3D tetrahedral cell (or element) has cellDim=3 and it has the three lower-level facet types faces (cellDim=2), edges (cellDim=1), and vertexes (cellDim=0).
Relative Cell: The term "relative cell" is used below to refer to the cell object for which we consider other cells called "facets" and "co-facets".
Facet: A lower dimensional cell directly connected to and contained within a relative parent cell. For example, a 3D cell (element) has three different facet types: faces (facetDim=2), edges (facetDim=1), and vertexes (facetDim=0). A face (cellDim=2) relative cell (or a 2D element) has two different facet types: edges (facetDim=1), and vertexes (facetDim=0). A particular relative cell (cellDim,cellLID) has nf=numFacets(cellDim,cellLID,facetDim) facets of a particular dimension facetDim.
Facet Dim (facetDim): The cell dimension of a given facet. For example, the dimension of a face facet of a 3D element is facetDim=2.
Facet Index (facetIndex): This is the local index of the facet with respect to its relative parent cell. For example, a 3D brick relative cell/element will have 6 faces which can be accessed using the facet indexes 0, 1, 2, 3, 4, and 5.
Co-facet: A higher dimensional cell of which the given relative cell is a facet. For example, the co-facets of a relative face cell in a 3D mesh are the one or two 3D element cells that share this face. In the case of a boundary face in a 3D mesh, only one co-facet is connected to the boundary cell and that is its boundary element. Note that in 3D, edges can have many different face and element co-facets. There are some important properties of co-facets. Relative cells of dimension cellDim=spatialDim-1 will have either one or two co-facets of dimension cofacetDim=spatialDim. Therefore, relative cells of the dimension cellDim=spatialDim-1 that have only one co-facet must be a boundary cell (i.e. a boundary face in 3D, a boundary edge in 2D, or a boundary vertex in 1D). Note that relative cells with dimension cellDim=spatialDim-n, where n >=2, can have more than two co-facets in general.
Co-facet Dimension (cofacetDim): The cell dimension of a co-facet cell.
Maximal co-facet: A co-facet of a relative cell whose dimension is equal to the spatial dimension of the mesh.
Assertion: maxCofacetDim == spatialDim.
Vertex: A zero dimensional cell (cellDim=0).
Assertion: vertexCellDim == 0.
Element: The term "element" is used to refer to cells whose dimension is equal to the spatial dimension of the mesh. Warning, the term "element" is an overloaded term and is also used to refer to an entry in an array of objects.
Facet orientation (facetOrientation): The orientation of a facet with respect to its parent cell. This is given as an integer value that requires agreement of interpretation. In 2D, an edge for instance will have just a positive and a negative orientation. In 3D, a face can have more than two orientations and each orientation is a permutation of the different nodal listings on the face. The facet orientation integer value is typically interpreted as an index into an array of possible orientations. The mesh reader/creator and the assembly code have to agree on how these orientations are interpreted.
Co-facet index (cofacetIndex): The relative index of a co-facet of some relative cell. For example, an internal face in a 3D mesh will have two co-facet element objects which can be accessed using the co-facet indexes 0 and 1.
Label: Every cell can be given a single integer value that can be used in a variety of ways. For example, boundary faces in 3D can be given an integer label to specify different types of boundary conditions. This integer could also be used with bit manipulation to allow for multiple classifications for a single cell.
Intermediate Cell: An intermediate cell is a cell that has a dimension cellDim in the range 0 < cellDim < spatialDim. For example, in a 3D mesh, faces and edges are intermediate cell types.
Geometric information
Point: A point is an array of coordinate positions. Currently this is stored as a array with three doubles (i.e. x=0, y=1, and z=2). For 2D meshes only the x=0 and y=1 positions are used. For 1D meshes, only the x=0 position is used.
Cell Centroid: A point that represents the geometric center of a cell.
Cell Diameter: The cell diameter is defined as the maximum line segment that can be drawn in the cell. This gives a measure of how big the cell is and, for instance, can be used in some types of finite-element stabilization methods.
Cell Curvature: The curvature of cell refers to the shape of a cell. This is typically used to represent the shapes of faces in 3D or edges in 2D which is critical in the implementation of higher-order higher-order discrimination methods. Cell curvature is currently not supported by this mesh interface but could be in the future.
Reference Cell: The term "reference cell" is used below to refer that a cell type's reference cell which is a geometric description of the cell in a simple coordinate system. By using a standard definition of a "reference cell", we can decouple different types of mathematical operations defined of a cell of a certain basic type (e.g. triangles and quadrilaterals in 2D, and tetrahedral and bricks in 3D). Typically, the nodes of a reference cell will have coordinates that integer multiples of one (including zero).
Cell Jacobian: The cell Jacobian is a transformation matrix that maps cell vertex coordinate values (i.e. points) from the reference cell to the physical cell. The inverse of the cell Jacobian maps points from the physical cell to the reference cell. A cell with the cell dimension cellDim has a cell Jacobian matrix of dimension cellDim x cellDim.
Parallel Information
Cell Ownership: The ownership of a cell refers to which process owns the cell. All other processes that need to access that cell get a "ghost" copy of the cell locally.
Ghost Cell: A ghost cell is a cell that is not owned by the local process that is copied to a process in order to perform certain computations.
Cell Local ID (cellLID): The LID of a cell is an index to the cell within the local process and these IDs are typically ordered from 0 to numLocalCells-1. Local cells can include both owned cells and ghosted cells.
Cell Global ID (cellGID): The global ID of a cell is an index that is unique across all processes and is used as a process-independent identifier for the cell.
Below, a common set of Mesh arguments are described along with relevant pre- and post conditions.
cellDim [in] The dimension of a cell
Precondition: 0 <= cellDim <= spatialDim
cellLID(s) [in] Local ID (LID) of a cell(s).
Precondition: 0 <= cellLID < numCells(cellDim), where cellDim is the dimension of the given cell.
cellGID [in] Global ID (GID) of a cell which is unique across processes.
Precondition: ??? <= cellGID < ??? (How do we write this?).
facetDim [in] The dimension of a given facet cell.
Precondition: 0 <= facetDim < cellDim, where cellDim is the dimension of the relative cell
facetIndex [in] Local index of the facet with respect to its parent relative cell.
Precondition: 0 <= facetIndex < numFacets(cellDim,cellLID,facetDim), where facetDim is the dimension of the facet of the relative parent cell (cellDim,cellLID)
facetOrientation(s) [out] The orientation of a facet with respect to its parent cell. This is given as an integer value that requires agreement of interpretation (see above).
cofacetDim [in] Cell dimension of a co-facet cell.
Precondition: 0 <= cellDim < cofacetDim <= spatialDim, where cellDim is the dimension of the relative cell
cofacetIndex [in] Relative index of a co-facet of some relative cell.
Precondition: 0 <= cofacetIndex < numMaxCofacets(cellDim,cellLID), where (cellDim,cellLID) is the relative cell.
Copy ObjectWithInstanceID base class into a special Nihilo tools package or into Teuchos.
Inherit MeshBase from Teuchos::VerboseObject instead of ObjectWithVerbosity.
Remove any concrete implementations and data from this base interface to provide a clean abstract interface. This would include removing the constructor.
Either remove the non-batch access functions all together, or make the namings of batch and non-batch access functions more consistent.
Have all of the functions return just a single return value or when more than one object is returned, return all values from the argument list. Or, when the value that is being returned from the argument list is only occasionally returned, change it to a pointer type and give it a default value of NULL so that it can be safely ignored! This will actually speed up performance in some cases. For example, facetOrientation is a return argument that is often ignored and would benefit from making it an optional argument.
Remove the staggered output function since the Teuchos::FancyOStream should handle parallel output from every process well.
Add a new cell type called a "general" cell for very irregular cells that defy fixed classification. This could be used to represent many of the irregular cells that remain after non-conforming mesh refinement.
Make sure to tell the DOF map how to interpret the facet orientation integers. Currently the interpretation is hard coded in Sundance. These permutations need to be enumerated in a C++ enum and the DOF map must be given the array of these enums to interpret the orientation integer values. To be 100 percent safe, we need a way of allowing the mesh to directly return the interpretation of these orientations. What you can do is to require that the int value be an index into an an array of enumeration values that corresponds to the cell type. This would require an enum type and an array of enum values for every facet cell type. For the current cell types for regular elements, this would be no problem. However, for the facets of a "general" cell type, the facets may also be classified as "general" cell types and therefore may not allow such an orientation array to be defined.
Clearly define the coordinate system of the reference cell of each of the standard cell types in order to clearly define interpretation of the reference points passed into the pushForward() function. Without a clear definition of the reference cell, these reference points are meaningless. Actually, for maximal flexibility, each cell type should define, at runtime, the coordinate system for the reference cell. Or, we must just all agree on a single set of standard definitions of reference cells for different cell types.
Definition at line 392 of file SundanceMeshBase.hpp.
| MeshBase::MeshBase | ( | int | dim, |
| const MPIComm & | comm, | ||
| const MeshEntityOrder & | meshOrder | ||
| ) |
Definition at line 56 of file SundanceMeshBase.cpp.
References curveDerivative_, curveNormal_, curvePoints_, dim_, and specialWeights_.
| virtual Sundance::MeshBase::~MeshBase | ( | ) | [inline, virtual] |
Definition at line 406 of file SundanceMeshBase.hpp.
| virtual bool Sundance::MeshBase::allowsHangingHodes | ( | ) | const [inline, virtual] |
Function returns true if the mesh allows hanging nodes (by refinement), false otherwise
Reimplemented in Sundance::HNMesh3D, and Sundance::HNMesh2D.
Definition at line 941 of file SundanceMeshBase.hpp.
| virtual void Sundance::MeshBase::assignIntermediateCellGIDs | ( | int | cellDim | ) | [pure virtual] |
Assign global IDs to cells of dimension cellDim.
| cellDim | [in] Dimension of the cell |
Postconditions:
hasIntermediateGIDs(cellDim)==true Implemented in Sundance::BasicSimplicialMesh, Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual CellType Sundance::MeshBase::cellType | ( | int | cellDim | ) | const [pure virtual] |
Get the cell type of the given cell dimension.
Note: This function by its very definition assumes that all cells of a given dimension have the same cell type within a mesh!
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| Point MeshBase::centroid | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [virtual] |
Return the centroid of a cell.
The default implementation just averages the positions of the zero-dimensional facets (i.e. vertexes).
Definition at line 84 of file SundanceMeshBase.cpp.
References facetLID(), nodePosition(), and numFacets().
Referenced by outwardNormals().
| const MPIComm& Sundance::MeshBase::comm | ( | ) | const [inline] |
Return the MPI communicator over which this mesh is distributed.
Definition at line 859 of file SundanceMeshBase.hpp.
References comm_.
Referenced by Sundance::BasicSimplicialMesh::addEdge(), Sundance::BasicSimplicialMesh::addElement(), Sundance::BasicSimplicialMesh::addFace(), Sundance::BasicSimplicialMesh::addVertex(), Sundance::BasicSimplicialMesh::assignIntermediateCellGIDs(), Sundance::BasicSimplicialMesh::printCells(), Sundance::BasicSimplicialMesh::resolveEdgeOwnership(), Sundance::BasicSimplicialMesh::synchronizeGIDNumbering(), and Sundance::BasicSimplicialMesh::synchronizeNeighborLists().
| const int * MeshBase::elemZeroFacetView | ( | int | maxCellLID | ) | const [virtual] |
Return a view of an array of LIDs for a maximum-dimensional cell's zero-dimensional facets (i.e. vertexes).
| maxCellLID | [in] Local ID of the maximum dimensional cell (i.e. element). |
Preconditions:
maxCellLID must be a maximum-dimensional cell! cellType(maxCellLID).Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 76 of file SundanceMeshBase.cpp.
| virtual int Sundance::MeshBase::facetLID | ( | int | cellDim, |
| int | cellLID, | ||
| int | facetDim, | ||
| int | facetIndex, | ||
| int & | facetOrientation | ||
| ) | const [pure virtual] |
Return the LID of a single facet cell (and optionally its orientation) with respect to a relative parent cell.
| cellDim | [in] Dimension of the parent cell whose facets are being obtained. |
| cellLID | [in] Local ID of the parent cell whose facets are being obtained |
| facetDim | [in] Dimension of the desired facet. |
| facetIndex | [in] The relative index into the list of the cell's facets. |
| facetOrientation | [out] The orientation of the facet w.r.t the parent cell. |
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by centroid(), getFacetArray(), outwardNormals(), and tangentsToEdges().
| void MeshBase::flushCurvePoints | ( | ) | const [virtual] |
detletes all the points and its normals which have been stored
Definition at line 266 of file SundanceMeshBase.cpp.
References curveDerivative_, curveID_to_ArrayIndex_, curveNormal_, curvePoints_, curvePoints_Are_Valid_, and nrCurvesForIntegral_.
| void MeshBase::flushSpecialWeights | ( | ) | const [virtual] |
deletes all special weights so those have to be recreated
Definition at line 255 of file SundanceMeshBase.cpp.
References dim_, specialWeights_, and validWeights_.
| Set< int > MeshBase::getAllLabelsForDimension | ( | int | cellDim | ) | const [virtual] |
Get the list of all labels defined for cells of the given dimension
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 225 of file SundanceMeshBase.cpp.
References label(), numCells(), and Sundance::Set< Key, Compare >::put().
| virtual void Sundance::MeshBase::getCellDiameters | ( | int | cellDim, |
| const Array< int > & | cellLIDs, | ||
| Array< double > & | diameters | ||
| ) | const [inline, virtual] |
Compute the diameters of a batch of cells.
| cellDim | [in] Dimension of the cells whose diameters are to be computed |
| cellLIDs | [in] Local IDs of the cells for which diameters are to be computed |
| diameters | [out] Array (size = cellLIDs.size()) of cell diameters. |
Warning! The default implementation returns an empty array of cell diameters! ToDo: Change the default implementation to compute diameters based on calls to the node position accessors. Going through the Mesh interface in that way will be less efficient than a low-level implementation, but would be a reasonable intermediate step for mesh developers.
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 783 of file SundanceMeshBase.hpp.
| virtual void Sundance::MeshBase::getCofacets | ( | int | cellDim, |
| int | cellLID, | ||
| int | cofacetDim, | ||
| Array< int > & | cofacetLIDs | ||
| ) | const [pure virtual] |
Return an array of the LIDs of all of the co-facets for a given relative cell.
| cellDim | [in] Dimension of the relative cell whose co-facets are being obtained |
| cellLID | [in] Local index of the relative cell whose co-facets are being obtained |
| cofacetDim | [in] Dimension of the co-facets to get |
| cofacetLIDs | [out] Array containing the LIDs of the co-facets for the given relative cell (cellDim,cellLID). |
getCofacetArray() to be consistent with getFacetArray()! Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| void MeshBase::getCurvePoints | ( | int | maxCellLID, |
| int | curveID, | ||
| Array< Point > & | points, | ||
| Array< Point > & | derivs, | ||
| Array< Point > & | normals | ||
| ) | const [virtual] |
Gets the points, curve derivatives and curve normals for one maxCell needed for curve/surf integral
Definition at line 299 of file SundanceMeshBase.cpp.
References curveDerivative_, curveNormal_, curvePoints_, mapCurveID_to_Index(), nrCurvesForIntegral_, and SUNDANCE_MSG3.
| void MeshBase::getFacetArray | ( | int | cellDim, |
| int | cellLID, | ||
| int | facetDim, | ||
| Array< int > & | facetLIDs, | ||
| Array< int > & | facetOrientations | ||
| ) | const |
Return an array containing the LIDs of the facets of dimension facetDim of the single relative parent cell (optionally also returning the facet orientations).
| cellDim | [in] Dimension of the parent cell whose facets are being obtained. |
| cellLID | [in] Local ID of the parent cell whose facets are being obtained |
| facetDim | [in] Dimension of the desired facets |
| facetLIDs | [out] On exit this array gives the local facet IDs for the parent cell with size = nf, where nf=numFacets(cellDim,cellLID,facetDim). |
| facetOrientations | [out] On exist, if facetDim > 0, this array gives the integer orientation of the facet with respect to its parent cell (see above definition of "Facet Orientation"). On exist this array will be resized to size = nf, where nf=numFacets(cellDim,cellLID,facetDim)). If facetDim==0 this this array argument is ignored! |
The default implementation loops over calls to facetLID(). Subclasses can provide a more efficient implementation if desired.
Array<int>* type so that it can be NULL and therefore ignored (which is a common use case). Definition at line 165 of file SundanceMeshBase.cpp.
References facetLID(), and numFacets().
| virtual void Sundance::MeshBase::getFacetLIDs | ( | int | cellDim, |
| const Array< int > & | cellLIDs, | ||
| int | facetDim, | ||
| Array< int > & | facetLIDs, | ||
| Array< int > & | facetOrientations | ||
| ) | const [pure virtual] |
Return an array containing the LIDs of facets of dimension facetDim for the given batch of relative parent cells.
| cellDim | [in] Dimension of the relative parent cells whose facets are being obtained |
| cellLIDs | [in] Array of LIDs for the relative parent cells |
| facetDim | [in] Dimension of the desired facets |
| facetLIDs | [out] On exit this array gives the local facet IDs for all of the cells given in cellLIDs in one flat array with size = cellLIDs.size()*nf, where nf=numFacets(cellDim,cellLIDs[j],facetDim)) where j can be any index 0 <= j < numCells(cellDim). Specifically, the local facet IDs for cellLIDs[k], where k=0...cellLIDs.size()-1, is given by facetLIDs[k*nf+j], where j=0...nf-1. |
| facetOrientations | [out] On exist, if facetDim > 0, this array gives the integer orientation of the facet with respect to its parent cell (see above definition of "Facet Orientation"). Oon exist this array will be resized to size = cellLIDs.size()*nf, where nf=numFacets(cellDim,cellLID,facetDim)). Specifically, the local facet orientation for the cell cellLIDs[k], where k=0...cellLIDs.size()-1, is given by facetOrientations[k*nf+j], where j=0...nf-1. If facetDim==0 then this array is ignored. |
Warning! This function will only work for homogeneous cell types!
Array<int>* type so that it can be NULL and therefore ignored (which is a common use case). Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual void Sundance::MeshBase::getJacobians | ( | int | cellDim, |
| const Array< int > & | cellLID, | ||
| CellJacobianBatch & | jBatch | ||
| ) | const [inline, virtual] |
Compute (or get) the Jacobians for a batch of cells.
| cellDim | [in] Dimension of the cells whose Jacobians are to be computed |
| cellLID | [in] Local IDs of the cells for which Jacobians are to be computed |
| jBatch | [out] Batch of cell Jacobians. |
Warning! The default implementation returns an empty batch of cell Jacobians!
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 757 of file SundanceMeshBase.hpp.
| void MeshBase::getLabels | ( | int | cellDim, |
| const Array< int > & | cellLIDs, | ||
| Array< int > & | labels | ||
| ) | const |
Get the labels for a batch of cells.
| cellDim | [in] Dimension of the parent cell whose facets are being obtained |
| cellLIDs | [in] Array of cell LIDs |
| labels | [out] On exit, contains an array (size=cellLIDs.size()) of the labels for each of the given cells. |
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 180 of file SundanceMeshBase.cpp.
References label().
| void MeshBase::getLIDsForLabel | ( | int | cellDim, |
| int | label, | ||
| Array< int > & | cellLIDs | ||
| ) | const [virtual] |
Get the cells associated with a specified label. The array cellLID will be filled with those cells of dimension cellDim having the given label.
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 212 of file SundanceMeshBase.cpp.
References label(), and numCells().
| void MeshBase::getMaxCofacetLIDs | ( | const Array< int > & | cellLIDs, |
| MaximalCofacetBatch & | cofacets | ||
| ) | const [virtual] |
Get the LIDs of the maximal cofacets for a batch of codimension-one cells. The default implementation simply loops over the cells in the cellLID array, taking no advantage of any internal data structures.
| cellLIDs | [in] array of LIDs of the cells whose cofacets are being obtained |
| cofacets | [out] |
| facetIndex | [out] index of each calling cell into the list of its maximal cofacet's facets |
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 188 of file SundanceMeshBase.cpp.
References Sundance::MaximalCofacetBatch::addSingleCofacet(), Sundance::MaximalCofacetBatch::addTwoCofacets(), maxCofacetLID(), numMaxCofacets(), Sundance::MaximalCofacetBatch::reset(), and spatialDim().
| void MeshBase::getSpecialWeight | ( | int | dim, |
| int | cellLID, | ||
| Array< double > & | w | ||
| ) | const [virtual] |
Returns the special weights
Definition at line 250 of file SundanceMeshBase.cpp.
References specialWeights_.
| bool MeshBase::hasCurvePoints | ( | int | maxCellLID, |
| int | curveID | ||
| ) | const [virtual] |
verifies if the specified maxCell has already precalculated quadrature point for one curve
Definition at line 280 of file SundanceMeshBase.cpp.
References curvePoints_, and mapCurveID_to_Index().
| virtual bool Sundance::MeshBase::hasGID | ( | int | cellDim, |
| int | cellGID | ||
| ) | const [pure virtual] |
Determine whether a given cell GID exists on this processor.
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual bool Sundance::MeshBase::hasIntermediateGIDs | ( | int | cellDim | ) | const [pure virtual] |
Return if cells of dimension cellDim have been assigned global IDs or not.
| cellDim | [in] Dimension of the cell |
Implemented in Sundance::BasicSimplicialMesh, Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| bool MeshBase::hasSpecialWeight | ( | int | dim, |
| int | cellLID | ||
| ) | const [virtual] |
verifies if the specified cell with the given dimension has special weights
Definition at line 239 of file SundanceMeshBase.cpp.
References specialWeights_.
| virtual int Sundance::MeshBase::indexInParent | ( | int | maxCellLID | ) | const [inline, virtual] |
Returns the index in the parent maxdim Cell of the refinement tree
| maxCellLID | [in] the LID of the cell |
Reimplemented in Sundance::HNMesh3D, and Sundance::HNMesh2D.
Definition at line 951 of file SundanceMeshBase.hpp.
| virtual bool Sundance::MeshBase::IsCurvePointsValid | ( | ) | const [inline, virtual] |
Definition at line 999 of file SundanceMeshBase.hpp.
References curvePoints_Are_Valid_.
| virtual bool Sundance::MeshBase::isElementHangingNode | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [inline, virtual] |
Function returns true if the specified element is a "hanging" element false otherwise
| cellDim | [in] should be between 0 , D-1 |
| cellLID | [in] the local ID of the element |
Reimplemented in Sundance::HNMesh3D, and Sundance::HNMesh2D.
Definition at line 947 of file SundanceMeshBase.hpp.
| virtual bool Sundance::MeshBase::IsSpecialWeightValid | ( | ) | const [inline, virtual] |
returns the status of the special weights if they are valid
These weights are usually computed for one setting of the curve (Adaptive Cell Integration)
Definition at line 973 of file SundanceMeshBase.hpp.
References validWeights_.
| virtual int Sundance::MeshBase::label | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [pure virtual] |
Get the label of the given cell.
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by getAllLabelsForDimension(), getLabels(), and getLIDsForLabel().
| int MeshBase::mapCurveID_to_Index | ( | int | curveID | ) | const [private, virtual] |
Definition at line 311 of file SundanceMeshBase.cpp.
References Sundance::Map< Key, Value, Compare >::containsKey(), curveDerivative_, curveID_to_ArrayIndex_, curveNormal_, curvePoints_, Sundance::Map< Key, Value, Compare >::get(), nrCurvesForIntegral_, Sundance::Map< Key, Value, Compare >::put(), and SUNDANCE_MSG3.
Referenced by getCurvePoints(), hasCurvePoints(), and setCurvePoints().
| virtual int Sundance::MeshBase::mapGIDToLID | ( | int | cellDim, |
| int | cellGID | ||
| ) | const [pure virtual] |
Find the LID of a cell given its GID.
| cellDim | [in] Dimension of the cell |
| cellGID | [in] Global ID of the cell |
Preconditions:
Postconditions:
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual int Sundance::MeshBase::mapLIDToGID | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [pure virtual] |
Find the global ID of a cell given its LID.
| cellDim | [in] Dimension of the cell |
| cellLID | [in] Local ID of the cell |
Postconditions:
returnVal >= 0 Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual int Sundance::MeshBase::maxChildren | ( | ) | const [inline, virtual] |
How many children has a refined element.
This function provides information of either we have bi or trisection
Reimplemented in Sundance::HNMesh3D, and Sundance::HNMesh2D.
Definition at line 955 of file SundanceMeshBase.hpp.
| virtual int Sundance::MeshBase::maxCofacetLID | ( | int | cellDim, |
| int | cellLID, | ||
| int | cofacetIndex, | ||
| int & | facetIndex | ||
| ) | const [pure virtual] |
Return the LID of a maximal co-facet of a cell.
| cellDim | [in] Dimension of the cell whose co-facets are being obtained |
| cellLID | [in] Local ID of the cell whose co-facets are being obtained |
| cofacetIndex | [in] Index into the list of the cell's co-facets. |
| facetIndex | [out] Returns the local index into the facet array for the relative cell (cellDim,cellLID) with respect to the maximal co-facet. In other words cellLID==facetLID(spatialDim(),returnVal,cellDim,facetIndex). |
Preconditions:
0 <= cellDim < spatialDim() facetIndex an int* argument and give it a default value of NUL so that it can be easily ignored! Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by getMaxCofacetLIDs(), and outwardNormals().
| const MeshEntityOrder& Sundance::MeshBase::meshOrder | ( | ) | const [inline] |
Get the ordering convention used by this mesh.
Definition at line 414 of file SundanceMeshBase.hpp.
References order_.
| virtual Point Sundance::MeshBase::nodePosition | ( | int | vertexLID | ) | const [pure virtual] |
Return the position of a local vertex.
| vertexLID | [in] The LID of the vertex. |
Preconditions:
0 <= vertexLID < this->numCells(0) getVertexPosition(...). Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by centroid(), outwardNormals(), and tangentsToEdges().
| virtual const double* Sundance::MeshBase::nodePositionView | ( | int | vertexLID | ) | const [pure virtual] |
Return a const view into an raw array for the position of a local vertex.
| vertexLID | [in] The LID of the vertex |
Preconditions:
0 <= vertexLID < this->numCells(0) spatialDim.getVertexPositionView(). Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual int Sundance::MeshBase::numCells | ( | int | cellDim | ) | const [pure virtual] |
Get the number of local cells having dimension cellDim.
numLocalCells(cellDim). Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by Sundance::CellReordererImplemBase::end(), getAllLabelsForDimension(), and getLIDsForLabel().
| virtual int Sundance::MeshBase::numFacets | ( | int | cellDim, |
| int | cellLID, | ||
| int | facetDim | ||
| ) | const [pure virtual] |
Return the number of facets of the given relative cell.
| cellDim | [in] The dimension of the relative cell |
| cellLID | [in] The LID of the relative cell |
| facetDim | [in] The dimension of the facets for the relative cell |
Postconditions:
returnVal >= 2: Every cell has two or more facets of a given dimension. Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by centroid(), and getFacetArray().
| virtual int Sundance::MeshBase::numMaxCofacets | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [pure virtual] |
Return the number of maximal co-facets of the given cell.
Preconditions:
0 <= cellDim < spatialDim() Note that if cellDim==spatialDim()-1 and returnVal==1 then this cell must be a boundary cell (i.e. a boundary face in 3D, a boundary edge in 1D, or a boundary node in 1D)!
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Referenced by getMaxCofacetLIDs(), and outwardNormals().
| void MeshBase::outwardNormals | ( | const Array< int > & | cellLIDs, |
| Array< Point > & | outwardNormals | ||
| ) | const [virtual] |
Get the outward normals for the batch of cells of dimension spatialDim()-1. If any cell in the batch is not on the boundary, an exception is thrown.
| cellLIDs | [in] LIDs for the cells whose normals are to be computed. |
| outwardNormals | [out] Outward normal unit vectors for each cell in the batch. |
Definition at line 95 of file SundanceMeshBase.cpp.
References centroid(), Sundance::cross(), facetLID(), maxCofacetLID(), nodePosition(), numMaxCofacets(), and spatialDim().
| virtual int Sundance::MeshBase::ownerProcID | ( | int | cellDim, |
| int | cellLID | ||
| ) | const [pure virtual] |
Return the rank of the processor that owns the given cell.
If returnVal==comm().getRank() then this cell is owned by this process.
Implemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| virtual void Sundance::MeshBase::pushForward | ( | int | cellDim, |
| const Array< int > & | cellLIDs, | ||
| const Array< Point > & | refPts, | ||
| Array< Point > & | physPts | ||
| ) | const [inline, virtual] |
Map points from a reference cell to physical points for a batch of cells.
| cellDim | [in] Dimension of the cells |
| cellLIDs | [in] Local IDs of a batch of cells |
| refPts | [in] Array of points on the single reference cell with respect to the reference cell's coordinate system. Note that the interpretation of these reference points is strictly determined by the coordinate system of the cell type cellType(cellDim) and must be clearly defined by this interface. |
| physPts | [out] Array (size = cellLIDs.size()*refPts.size()) of the physical points given in a flat array for the batch of cells. Specifically, the physical points for each cell cellLIDs[k], for k=0...cellLIDs.size()-1, is given by physPts[k*nrp+j], for j=0...nrp-1 and nrp=refPts.size(). |
Warning! The default implementation returns an empty array of physical points!
Reimplemented in Sundance::HNMesh3D, Sundance::HNMesh2D, Sundance::BasicSimplicialMesh, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
Definition at line 846 of file SundanceMeshBase.hpp.
| const CellReordererImplemBase* Sundance::MeshBase::reorderer | ( | ) | const [inline] |
Get the reordering strategy to be used with this mesh.
Definition at line 930 of file SundanceMeshBase.hpp.
References reorderer_.
| virtual void Sundance::MeshBase::returnParentFacets | ( | int | childCellLID, |
| int | dimFacets, | ||
| Array< int > & | facetsLIDs, | ||
| int & | parentCellLIDs | ||
| ) | const [inline, virtual] |
Function returns the facets of the parent cell (needed for HN treatment)
| childCellLID | [in] the LID of the maxdim cell, whos parents facets we want |
| dimFacets | [in] the dimension of the facets which we want to have |
| facetsLIDs | [out] the LID of the parents facets (all) in the defined order |
| parentCellLIDs | [out] the maxdim parent cell LID |
Reimplemented in Sundance::HNMesh3D, and Sundance::HNMesh2D.
Definition at line 962 of file SundanceMeshBase.hpp.
| void MeshBase::setCurvePoints | ( | int | maxCellLID, |
| int | curveID, | ||
| Array< Point > & | points, | ||
| Array< Point > & | derivs, | ||
| Array< Point > & | normals | ||
| ) | const [virtual] |
Sets the points, curve derivatives and curve normals for one maxCell needed for curve/surf integral
Definition at line 287 of file SundanceMeshBase.cpp.
References curveDerivative_, curveNormal_, curvePoints_, mapCurveID_to_Index(), nrCurvesForIntegral_, and SUNDANCE_MSG3.
| virtual void Sundance::MeshBase::setCurvePointsValid | ( | bool | val | ) | const [inline, virtual] |
Definition at line 1002 of file SundanceMeshBase.hpp.
References curvePoints_Are_Valid_.
| virtual void Sundance::MeshBase::setLabel | ( | int | cellDim, |
| int | cellLID, | ||
| int | label | ||
| ) | [pure virtual] |
Set the label for the given cell.
Implemented in Sundance::HNMesh3D, Sundance::BasicSimplicialMesh, Sundance::HNMesh2D, Sundance::PeriodicMesh1D, and Sundance::PeriodicSingleCellMesh1D.
| void Sundance::MeshBase::setReorderer | ( | const CellReorderer & | reorderer | ) | [inline] |
Set the reordering strategy to be used with this mesh.
Definition at line 926 of file SundanceMeshBase.hpp.
References Sundance::CellReorderer::createInstance(), and reorderer_.
| void MeshBase::setSpecialWeight | ( | int | dim, |
| int | cellLID, | ||
| Array< double > & | w | ||
| ) | const [virtual] |
Sets the special weights
Definition at line 246 of file SundanceMeshBase.cpp.
References specialWeights_.
| virtual void Sundance::MeshBase::setSpecialWeightValid | ( | bool | val | ) | const [inline, virtual] |
specifies if the special weights are valid
if this is false then usually the special weights have to be recomputed
Definition at line 977 of file SundanceMeshBase.hpp.
References validWeights_.
| int Sundance::MeshBase::spatialDim | ( | ) | const [inline] |
Get the spatial dimension of the mesh.
Postconditions:
0 < returnVal <= 3 Definition at line 422 of file SundanceMeshBase.hpp.
References dim_.
Referenced by Sundance::BasicSimplicialMesh::addEdge(), Sundance::BasicSimplicialMesh::addElement(), Sundance::BasicSimplicialMesh::addFace(), Sundance::BasicSimplicialMesh::addVertex(), Sundance::BasicSimplicialMesh::BasicSimplicialMesh(), Sundance::BasicSimplicialMesh::cellToStr(), Sundance::CellReordererImplemBase::end(), Sundance::BasicSimplicialMesh::estimateNumElements(), Sundance::BasicSimplicialMesh::facetLID(), Sundance::BasicSimplicialMesh::getCellDiameters(), Sundance::HNMesh2D::getCellDiameters(), Sundance::HNMesh3D::getCellDiameters(), Sundance::BasicSimplicialMesh::getCofacets(), Sundance::BasicSimplicialMesh::getFacetLIDs(), Sundance::PeriodicSingleCellMesh1D::getJacobians(), Sundance::PeriodicMesh1D::getJacobians(), Sundance::BasicSimplicialMesh::getJacobians(), Sundance::HNMesh2D::getJacobians(), Sundance::HNMesh3D::getJacobians(), getMaxCofacetLIDs(), Sundance::BasicSimplicialMesh::maxCofacetLID(), outwardNormals(), Sundance::BasicSimplicialMesh::printCells(), Sundance::BasicSimplicialMesh::pushForward(), Sundance::HNMesh2D::pushForward(), Sundance::HNMesh3D::pushForward(), and tangentsToEdges().
| static bool& Sundance::MeshBase::staggerOutput | ( | ) | [inline, static] |
Set whether to stagger output in parallel. Set to true for readable output in parallel debugging sessions.
This should be normally, as it causes one synchronization point per process.
Definition at line 1033 of file SundanceMeshBase.hpp.
Referenced by Sundance::BasicSimplicialMesh::assignIntermediateCellGIDs(), and Sundance::BasicSimplicialMesh::resolveEdgeOwnership().
| void MeshBase::tangentsToEdges | ( | const Array< int > & | cellLIDs, |
| Array< Point > & | tangenVectors | ||
| ) | const [virtual] |
Get tangent vectors for a batch of edges
| cellLIDs | [in] LIDs for the cells whose tangents are to be computed. |
| tangentVectors | [out] Unit tangents for each cell |
Definition at line 143 of file SundanceMeshBase.cpp.
References facetLID(), nodePosition(), and spatialDim().
MPIComm Sundance::MeshBase::comm_ [private] |
Definition at line 1041 of file SundanceMeshBase.hpp.
Referenced by comm().
Array< Sundance::Map< int , Array<Point> > > Sundance::MeshBase::curveDerivative_ [mutable, private] |
store directional derivative informations to calculate the curve integral
Definition at line 1067 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), getCurvePoints(), mapCurveID_to_Index(), MeshBase(), and setCurvePoints().
Sundance::Map< int , int > Sundance::MeshBase::curveID_to_ArrayIndex_ [mutable, private] |
map the curve ID to index in the previous arrays
Definition at line 1074 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), and mapCurveID_to_Index().
Array< Sundance::Map< int , Array<Point> > > Sundance::MeshBase::curveNormal_ [mutable, private] |
store normal directional used in the curve or in the surf integral
in case of the surf integral it is the cross product from the integral
Definition at line 1071 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), getCurvePoints(), mapCurveID_to_Index(), MeshBase(), and setCurvePoints().
Array< Sundance::Map< int , Array<Point> > > Sundance::MeshBase::curvePoints_ [mutable, private] |
store intersection informations to calculate the curve integral
Definition at line 1064 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), getCurvePoints(), hasCurvePoints(), mapCurveID_to_Index(), MeshBase(), and setCurvePoints().
bool Sundance::MeshBase::curvePoints_Are_Valid_ [mutable, private] |
true if the curve did not moved, false if those points are not reusable
Definition at line 1058 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), IsCurvePointsValid(), and setCurvePointsValid().
int Sundance::MeshBase::dim_ [private] |
Definition at line 1039 of file SundanceMeshBase.hpp.
Referenced by flushSpecialWeights(), MeshBase(), and spatialDim().
int Sundance::MeshBase::nrCurvesForIntegral_ [mutable, private] |
how many curves are participating in curve integrals
Definition at line 1061 of file SundanceMeshBase.hpp.
Referenced by flushCurvePoints(), getCurvePoints(), mapCurveID_to_Index(), and setCurvePoints().
MeshEntityOrder Sundance::MeshBase::order_ [private] |
Definition at line 1043 of file SundanceMeshBase.hpp.
Referenced by meshOrder().
RCP<CellReordererImplemBase> Sundance::MeshBase::reorderer_ [private] |
Definition at line 1045 of file SundanceMeshBase.hpp.
Referenced by reorderer(), and setReorderer().
Array< Sundance::Map< int , Array<double> > > Sundance::MeshBase::specialWeights_ [mutable, private] |
Object to store the special weights for integration
Definition at line 1053 of file SundanceMeshBase.hpp.
Referenced by flushSpecialWeights(), getSpecialWeight(), hasSpecialWeight(), MeshBase(), and setSpecialWeight().
bool Sundance::MeshBase::validWeights_ [mutable, private] |
flag to indicate if the weights stored are valid
Definition at line 1050 of file SundanceMeshBase.hpp.
Referenced by flushSpecialWeights(), IsSpecialWeightValid(), and setSpecialWeightValid().