00001 /* @HEADER@ */ 00002 // ************************************************************************ 00003 // 00004 // Sundance 00005 // Copyright 2011 Sandia Corporation 00006 // 00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00008 // the U.S. Government retains certain rights in this software. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu) 00038 // 00039 00040 /* @HEADER@ */ 00041 00042 #ifndef SUNDANCE_MESHBASE_H 00043 #define SUNDANCE_MESHBASE_H 00044 00045 #include "SundanceDefs.hpp" 00046 #include "SundancePoint.hpp" 00047 #include "SundanceSet.hpp" 00048 #include "SundanceMap.hpp" 00049 #include "SundanceCellType.hpp" 00050 #include "SundanceObjectWithVerbosity.hpp" 00051 #include "SundanceObjectWithInstanceID.hpp" 00052 #include "PlayaMPIComm.hpp" 00053 #include "Teuchos_Array.hpp" 00054 #include "SundanceCellReorderer.hpp" 00055 #include "SundanceCellReordererImplemBase.hpp" 00056 00057 namespace Sundance { 00058 00059 00060 /** Identifier for ordering convention */ 00061 enum MeshEntityOrder {UFCMeshOrder, ExodusMeshOrder}; 00062 00063 class MaximalCofacetBatch; 00064 00065 using namespace Teuchos; 00066 using Playa::MPIComm; 00067 00068 class CellJacobianBatch; 00069 00070 00071 /** \brief Abstract interface to a mesh. 00072 * 00073 * \section Sundance_MeshBase_Outline_sec Outline 00074 * 00075 * <ul> 00076 * <li>\ref Sundance_MeshBase_Introduction_sec 00077 * <li>\ref Sundance_MeshBase_Definitions_sec 00078 * <li>\ref Sundance_MeshBase_Common_Function_Arguments_sec 00079 * <li>\ref Sundance_MeshBase_Refactorings_sec 00080 * </ul> 00081 * 00082 * \section Sundance_MeshBase_Introduction_sec Introduction 00083 * 00084 * A mesh is a collection of connected cells (also called elements). This 00085 * mesh interface is designed to allow for the representation of 1D, 2D, and 3D (and 00086 * even higher dimensional) meshes. A mesh contains both topological 00087 * information and geometric information about the cells and their subparts 00088 * (i.e. facets). Topological information refers to information about how 00089 * cells of various types and dimensions are connected to each other within 00090 * the mesh. Geometric information refers to information about the position, 00091 * size, and shape of cells within the mesh. 00092 * 00093 * Currently, this interface only handles meshes where all cells of a given 00094 * cell dimension have the same cell type. For example, this only allows 00095 * meshes with all triangles or all quads in 2D or all tetrahedrals in 3D. 00096 * 00097 * \todo Add some diagrams to show different cell types to help define cells, 00098 * facets, co-facets, vertices, etc. 00099 * 00100 * \section Sundance_MeshBase_Definitions_sec Definitions 00101 * 00102 * The following definitions are used in the specification of this mesh 00103 * interface. 00104 * 00105 *<ul> 00106 * 00107 * <li><b>Topological information</b> 00108 * 00109 * <ul> 00110 * 00111 * <li><b>Spatial Dimension (spatialDim)</b>: The maximum cell dimension of 00112 * any cell type in the mesh. For example, a 3D mesh will have 00113 * spatialDim=3, a 2D mesh will have spatialDim=2 and a 1D mesh will have 00114 * spatialDim=1. Typically, when we use the short hand 1D, 2D or 3D, we are 00115 * referring to the spatial dimension of the mesh. Note that 4D and 00116 * higher-dimensional meshes are also representable by this interface in 00117 * principle with some small changes (e.g. like making the size of a point 00118 * more flexible) 00119 * 00120 * <li><b>Cell</b>: A mesh object of some dimension. For example, in a 3D 00121 * mesh, cells are things like elements, faces, edges, and nodes. 00122 * 00123 * <li><b>Cell Type (cellType)</b>: An enumeration of type <tt>CellType</tt> 00124 * that lists some common cell types. For example, common cell types are 00125 * things like vertexes (0D), lines or edges (1D), triangles (2D), and 00126 * quadrilaterals (2D). 00127 * 00128 * <li><b>Cell Dimension (cellDim)</b>: The dimension of a cell. More 00129 * formally, the cell dimension is equal to the number of lower level facet 00130 * types in the cell. For example, A 3D tetrahedral cell (or element) has 00131 * cellDim=3 and it has the three lower-level facet types faces (cellDim=2), 00132 * edges (cellDim=1), and vertexes (cellDim=0). 00133 * 00134 * <li><b>Relative Cell</b>: The term "relative cell" is used below to refer 00135 * to the cell object for which we consider other cells called "facets" and 00136 * "co-facets". 00137 * 00138 * <li><b>Facet</b>: A lower dimensional cell directly connected to and 00139 * contained within a relative <em>parent</em> cell. For example, a 3D cell 00140 * (element) has three different facet types: faces (facetDim=2), edges 00141 * (facetDim=1), and vertexes (facetDim=0). A face (cellDim=2) relative 00142 * cell (or a 2D element) has two different facet types: edges (facetDim=1), 00143 * and vertexes (facetDim=0). A particular relative cell 00144 * <tt>(cellDim,cellLID)</tt> has 00145 * <tt>nf=numFacets(cellDim,cellLID,facetDim)</tt> facets of a particular 00146 * dimension <tt>facetDim</tt>. 00147 * 00148 * <li><b>Facet Dim (facetDim)</b>: The cell dimension of a given facet. 00149 * For example, the dimension of a face facet of a 3D element is facetDim=2. 00150 * 00151 * <li><b>Facet Index (facetIndex)</b>: This is the local index of the facet 00152 * with respect to its relative parent cell. For example, a 3D brick 00153 * relative cell/element will have 6 faces which can be accessed using the 00154 * facet indexes 0, 1, 2, 3, 4, and 5. 00155 * 00156 * <li><b>Co-facet</b>: A higher dimensional cell of which the given 00157 * relative cell is a facet. For example, the co-facets of a relative 00158 * face cell in a 3D mesh are the one or two 3D element cells that share 00159 * this face. In the case of a boundary face in a 3D mesh, only one 00160 * co-facet is connected to the boundary cell and that is its boundary 00161 * element. Note that in 3D, edges can have many different face and element 00162 * co-facets. There are some important properties of co-facets. Relative 00163 * cells of dimension cellDim=spatialDim-1 will have either one or two 00164 * co-facets of dimension cofacetDim=spatialDim. Therefore, relative cells 00165 * of the dimension cellDim=spatialDim-1 that have only one co-facet must be 00166 * a boundary cell (i.e. a boundary face in 3D, a boundary edge in 2D, or a 00167 * boundary vertex in 1D). Note that relative cells with dimension 00168 * cellDim=spatialDim-n, where n >=2, can have more than two co-facets in 00169 * general. 00170 * 00171 * <li><b>Co-facet Dimension (cofacetDim)</b>: The cell dimension of a 00172 * co-facet cell. 00173 * 00174 * <li><b>Maximal co-facet</b>: A co-facet of a relative cell whose 00175 * dimension is equal to the spatial dimension of the mesh. 00176 * 00177 * <em>Assertion</em>: <tt>maxCofacetDim == spatialDim</tt>. 00178 * 00179 * <li><b>Vertex</b>: A zero dimensional cell (cellDim=0). 00180 * 00181 * <em>Assertion</em>: <tt>vertexCellDim == 0</tt>. 00182 * 00183 * <li><b>Element</b>: The term "element" is used to refer to cells whose 00184 * dimension is equal to the spatial dimension of the mesh. Warning, the 00185 * term "element" is an overloaded term and is also used to refer to an 00186 * entry in an array of objects. 00187 * 00188 * <li><b>Facet orientation (facetOrientation)</b>: The orientation of a 00189 * facet with respect to its parent cell. This is given as an integer value 00190 * that requires agreement of interpretation. In 2D, an edge for instance 00191 * will have just a positive and a negative orientation. In 3D, a face can 00192 * have more than two orientations and each orientation is a permutation of 00193 * the different nodal listings on the face. The facet orientation integer 00194 * value is typically interpreted as an index into an array of possible 00195 * orientations. The mesh reader/creator and the assembly code have to 00196 * agree on how these orientations are interpreted. 00197 * 00198 * <li><b>Co-facet index (cofacetIndex)</b>: The relative index of a 00199 * co-facet of some relative cell. For example, an internal face in a 3D 00200 * mesh will have two co-facet element objects which can be accessed using 00201 * the co-facet indexes 0 and 1. 00202 * 00203 * <li><b>Label</b>: Every cell can be given a single integer value that 00204 * can be used in a variety of ways. For example, boundary faces in 3D can 00205 * be given an integer label to specify different types of boundary 00206 * conditions. This integer could also be used with bit manipulation to 00207 * allow for multiple classifications for a single cell. 00208 * 00209 * <li><b>Intermediate Cell</b>: An intermediate cell is a cell that has a 00210 * dimension cellDim in the range 0 < cellDim < spatialDim. For example, in 00211 * a 3D mesh, faces and edges are intermediate cell types. 00212 * 00213 * </ul> 00214 * 00215 * <li><b>Geometric information</b> 00216 * 00217 * <ul> 00218 * 00219 * <li><b>Point</b>: A point is an array of coordinate positions. Currently 00220 * this is stored as a array with three doubles (i.e. x=0, y=1, and z=2). 00221 * For 2D meshes only the x=0 and y=1 positions are used. For 1D meshes, 00222 * only the x=0 position is used. 00223 * 00224 * <li><b>Cell Centroid</b>: A point that represents the geometric center of 00225 * a cell. 00226 * 00227 * <li> <b>Cell Diameter</b>: The cell diameter is defined as the maximum 00228 * line segment that can be drawn in the cell. This gives a measure of how 00229 * big the cell is and, for instance, can be used in some types of 00230 * finite-element stabilization methods. 00231 * 00232 * <li><b>Cell Curvature</b>: The curvature of cell refers to the shape of a 00233 * cell. This is typically used to represent the shapes of faces in 3D or 00234 * edges in 2D which is critical in the implementation of higher-order 00235 * higher-order discrimination methods. Cell curvature is currently not 00236 * supported by this mesh interface but could be in the future. 00237 * 00238 * <li><b>Reference Cell</b>: The term "reference cell" is used below to 00239 * refer that a cell type's reference cell which is a geometric description 00240 * of the cell in a simple coordinate system. By using a standard 00241 * definition of a "reference cell", we can decouple different types of 00242 * mathematical operations defined of a cell of a certain basic type (e.g. 00243 * triangles and quadrilaterals in 2D, and tetrahedral and bricks in 3D). 00244 * Typically, the nodes of a reference cell will have coordinates that 00245 * integer multiples of one (including zero). 00246 * 00247 * <li><b>Cell Jacobian</b>: The cell Jacobian is a transformation matrix 00248 * that maps cell vertex coordinate values (i.e. points) from the reference 00249 * cell to the physical cell. The inverse of the cell Jacobian maps points 00250 * from the physical cell to the reference cell. A cell with the cell 00251 * dimension cellDim has a cell Jacobian matrix of dimension cellDim x 00252 * cellDim. 00253 * 00254 * </ul> 00255 * 00256 * <li><b>Parallel Information</b> 00257 * 00258 * <ul> 00259 * 00260 * <li><b>Cell Ownership</b>: The ownership of a cell refers to which 00261 * process owns the cell. All other processes that need to access that cell 00262 * get a "ghost" copy of the cell locally. 00263 * 00264 * <li><b>Ghost Cell</b>: A ghost cell is a cell that is not owned by the 00265 * local process that is copied to a process in order to perform certain 00266 * computations. 00267 * 00268 * <li><b>Cell Local ID (cellLID)</b>: The LID of a cell is an index to 00269 * the cell within the local process and these IDs are typically ordered 00270 * from 0 to numLocalCells-1. Local cells can include both owned cells and 00271 * ghosted cells. 00272 * 00273 * <li><b>Cell Global ID (cellGID)</b>: The global ID of a cell is an index 00274 * that is unique across all processes and is used as a process-independent 00275 * identifier for the cell. 00276 * 00277 * </ul> 00278 * 00279 * </ul> 00280 * 00281 * \section Sundance_MeshBase_Common_Function_Arguments_sec Common Function Arguments and Pre/Post-Conditions 00282 * 00283 * Below, a common set of Mesh arguments are described along with relevant 00284 * pre- and post conditions. 00285 * 00286 * <ul> 00287 * 00288 * <li><b>cellDim</b> [in] The dimension of a cell 00289 * 00290 * <em>Precondition</em>: <tt>0 <= cellDim <= spatialDim</tt> 00291 * 00292 * <li><b>cellLID(s)</b> [in] Local ID (LID) of a cell(s). 00293 * 00294 * <em>Precondition</em>: <tt>0 <= cellLID < numCells(cellDim)</tt>, where 00295 * <tt>cellDim</tt> is the dimension of the given cell. 00296 * 00297 * <li><b>cellGID</b> [in] Global ID (GID) of a cell which is unique across 00298 * processes. 00299 * 00300 * <em>Precondition</em>: <tt>??? <= cellGID < ???</tt> (How do we write 00301 * this?). 00302 * 00303 * <li><b>facetDim</b> [in] The dimension of a given facet cell. 00304 * 00305 * <em>Precondition</em>: <tt>0 <= facetDim < cellDim</tt>, where 00306 * <tt>cellDim</tt> is the dimension of the relative cell 00307 * 00308 * <li><b>facetIndex</b> [in] Local index of the facet with respect to its 00309 * parent relative cell. 00310 * 00311 * <em>Precondition</em>: <tt>0 <= facetIndex < 00312 * numFacets(cellDim,cellLID,facetDim)</tt>, where <tt>facetDim</tt> is the 00313 * dimension of the facet of the relative parent cell 00314 * <tt>(cellDim,cellLID)</tt> 00315 * 00316 * <li><b>facetOrientation(s)</b> [out] The orientation of a facet with 00317 * respect to its parent cell. This is given as an integer value that 00318 * requires agreement of interpretation (see above). 00319 * 00320 * <li><b>cofacetDim</b> [in] Cell dimension of a co-facet cell. 00321 * 00322 * <em>Precondition</em>: <tt>0 <= cellDim < cofacetDim <= spatialDim</tt>, 00323 * where <tt>cellDim</tt> is the dimension of the relative cell 00324 * 00325 * <li><b>cofacetIndex</b> [in] Relative index of a co-facet of some relative 00326 * cell. 00327 * 00328 * <em>Precondition</em>: <tt>0 <= cofacetIndex < 00329 * numMaxCofacets(cellDim,cellLID)</tt>, where <tt>(cellDim,cellLID)</tt> is 00330 * the relative cell. 00331 * 00332 * </ul> 00333 * 00334 * \section Sundance_MeshBase_Refactorings_sec Possible Refactorings 00335 * 00336 * <ul> 00337 * 00338 * <li>Copy ObjectWithInstanceID base class into a special Nihilo tools 00339 * package or into Teuchos. 00340 * 00341 * <li>Inherit MeshBase from Teuchos::VerboseObject instead of 00342 * ObjectWithVerbosity. 00343 * 00344 * <li>Remove any concrete implementations and data from this base interface 00345 * to provide a clean abstract interface. This would include removing the 00346 * constructor. 00347 * 00348 * <li>Either remove the non-batch access functions all together, or make the 00349 * namings of batch and non-batch access functions more consistent. 00350 * 00351 * <li>Have all of the functions return just a single return value or when 00352 * more than one object is returned, return all values from the argument list. 00353 * Or, when the value that is being returned from the argument list is only 00354 * occasionally returned, change it to a pointer type and give it a default 00355 * value of NULL so that it can be safely ignored! This will actually speed 00356 * up performance in some cases. For example, <tt>facetOrientation</tt> is a 00357 * return argument that is often ignored and would benefit from making it an 00358 * optional argument. 00359 * 00360 * <li>Remove the staggered output function since the Teuchos::FancyOStream 00361 * should handle parallel output from every process well. 00362 * 00363 * <li>Add a new cell type called a "general" cell for very irregular cells 00364 * that defy fixed classification. This could be used to represent many of 00365 * the irregular cells that remain after non-conforming mesh refinement. 00366 * 00367 * <li>Make sure to tell the DOF map how to interpret the facet orientation 00368 * integers. Currently the interpretation is hard coded in Sundance. These 00369 * permutations need to be enumerated in a C++ enum and the DOF map must be 00370 * given the array of these enums to interpret the orientation integer values. 00371 * To be 100 percent safe, we need a way of allowing the mesh to directly 00372 * return the interpretation of these orientations. What you can do is to 00373 * require that the int value be an index into an an array of enumeration 00374 * values that corresponds to the cell type. This would require an enum type 00375 * and an array of enum values for every facet cell type. For the current 00376 * cell types for regular elements, this would be no problem. However, for 00377 * the facets of a "general" cell type, the facets may also be classified as 00378 * "general" cell types and therefore may not allow such an orientation array 00379 * to be defined. 00380 * 00381 * <li>Clearly define the coordinate system of the reference cell of each of 00382 * the standard cell types in order to clearly define interpretation of the 00383 * reference points passed into the <tt>pushForward()</tt> function. Without 00384 * a clear definition of the reference cell, these reference points are 00385 * meaningless. Actually, for maximal flexibility, each cell type should 00386 * define, at runtime, the coordinate system for the reference cell. Or, we 00387 * must just all agree on a single set of standard definitions of reference cells 00388 * for different cell types. 00389 * 00390 * </ul> 00391 */ 00392 class MeshBase 00393 : public ObjectWithClassVerbosity<MeshBase> 00394 , public ObjectWithInstanceID<MeshBase> 00395 { 00396 public: 00397 00398 /** \name Constructors/Destructor */ 00399 //@{ 00400 00401 /** \brief . */ 00402 MeshBase(int dim, const MPIComm& comm, 00403 const MeshEntityOrder& meshOrder); 00404 00405 /** \brief . */ 00406 virtual ~MeshBase(){;} 00407 00408 //@} 00409 00410 /** \name Topological Information */ 00411 //@{ 00412 00413 /** \brief Get the ordering convention used by this mesh */ 00414 const MeshEntityOrder& meshOrder() const {return order_;} 00415 00416 /** \brief Get the spatial dimension of the mesh 00417 * 00418 * <b>Postconditions:</b><ul> 00419 * <li><tt>0 < returnVal <= 3</tt> 00420 * </ul> 00421 */ 00422 int spatialDim() const {return dim_;} 00423 00424 /** \brief Get the number of local cells having dimension cellDim. 00425 * 00426 * \todo Change this to <tt>numLocalCells(cellDim)</tt>. 00427 */ 00428 virtual int numCells(int cellDim) const = 0 ; 00429 00430 /** \brief Return the number of facets of the given relative cell. 00431 * 00432 * \param cellDim 00433 * [in] The dimension of the relative cell 00434 * \param cellLID 00435 * [in] The LID of the relative cell 00436 * \param facetDim 00437 * [in] The dimension of the facets for the 00438 * relative cell 00439 * 00440 * <b>Postconditions:</b><ul> 00441 * <li><tt>returnVal >= 2</tt>: Every cell has two or more 00442 * facets of a given dimension. 00443 * </ul> 00444 */ 00445 virtual int numFacets(int cellDim, int cellLID, 00446 int facetDim) const = 0 ; 00447 00448 /** \brief Return the LID of a single facet cell (and optionally its 00449 * orientation) with respect to a relative parent cell. 00450 * 00451 * \param cellDim 00452 * [in] Dimension of the parent cell whose facets 00453 * are being obtained. 00454 * \param cellLID 00455 * [in] Local ID of the parent cell whose facets 00456 * are being obtained 00457 * \param facetDim 00458 * [in] Dimension of the desired facet. 00459 * \param facetIndex 00460 * [in] The relative index into the list of the 00461 * cell's facets. 00462 * \param facetOrientation 00463 * [out] The orientation of the facet w.r.t the 00464 * parent cell. 00465 * 00466 * \todo Change the facetOrientation argument to be a raw pointer that can 00467 * be NULL and therefore easily ignored. 00468 */ 00469 virtual int facetLID(int cellDim, int cellLID, 00470 int facetDim, int facetIndex, 00471 int& facetOrientation) const = 0 ; 00472 00473 /** \brief Return an array containing the LIDs of facets of dimension 00474 * facetDim for the given batch of relative parent cells. 00475 * 00476 * \param cellDim 00477 * [in] Dimension of the relative parent cells whose facets 00478 * are being obtained 00479 * \param cellLIDs 00480 * [in] Array of LIDs for the relative parent cells 00481 * \param facetDim 00482 * [in] Dimension of the desired facets 00483 * \param facetLIDs 00484 * [out] On exit this array gives the local facet IDs 00485 * for all of the cells given in cellLIDs in one flat array 00486 * with size = <tt>cellLIDs.size()*nf</tt>, where 00487 * <tt>nf=numFacets(cellDim,cellLIDs[j],facetDim)</tt>) where 00488 * <tt>j</tt> can be any index <tt>0 <= j < numCells(cellDim)</tt>. 00489 * Specifically, the local facet IDs for <tt>cellLIDs[k]</tt>, where 00490 * <tt>k=0...cellLIDs.size()-1</tt>, is given 00491 * by <tt>facetLIDs[k*nf+j]</tt>, where <tt>j=0...nf-1</tt>. 00492 * \param facetOrientations 00493 * [out] On exist, if <tt>facetDim > 0</tt>, this array gives 00494 * the integer orientation of the facet with respect to 00495 * its parent cell (see above definition of "Facet Orientation"). 00496 * Oon exist this array will be 00497 * resized to size = <tt>cellLIDs.size()*nf</tt>, where 00498 * <tt>nf=numFacets(cellDim,cellLID,facetDim)</tt>). 00499 * Specifically, the local facet orientation for the cell cellLIDs[k], where 00500 * <tt>k=0...cellLIDs.size()-1</tt>, is given 00501 * by <tt>facetOrientations[k*nf+j]</tt>, where <tt>j=0...nf-1</tt>. 00502 * If <tt>facetDim==0</tt> then this array is ignored. 00503 * 00504 * <b>Warning!</b> This function will only work for homogeneous cell types! 00505 * 00506 * \todo Change the facetOrientation argument to an <tt>Array<int>*</tt> 00507 * type so that it can be NULL and therefore ignored (which is a common use 00508 * case). 00509 */ 00510 virtual void getFacetLIDs(int cellDim, 00511 const Array<int>& cellLIDs, 00512 int facetDim, 00513 Array<int>& facetLIDs, 00514 Array<int>& facetOrientations) const = 0 ; 00515 00516 /** \brief Return an array containing the LIDs of the facets of 00517 * dimension facetDim of the single relative parent cell (optionally also 00518 * returning the facet orientations). 00519 * 00520 * \param cellDim 00521 * [in] Dimension of the parent cell whose facets 00522 * are being obtained. 00523 * \param cellLID 00524 * [in] Local ID of the parent cell whose facets 00525 * are being obtained 00526 * \param facetDim 00527 * [in] Dimension of the desired facets 00528 * \param facetLIDs 00529 * [out] On exit this array gives the local facet IDs 00530 * for the parent cell 00531 * with size = <tt>nf</tt>, where 00532 * <tt>nf=numFacets(cellDim,cellLID,facetDim)</tt>. 00533 * \param facetOrientations 00534 * [out] On exist, if <tt>facetDim > 0</tt>, this array gives 00535 * the integer orientation of the facet with respect to 00536 * its parent cell (see above definition of "Facet Orientation"). 00537 * On exist this array will be 00538 * resized to size = <tt>nf</tt>, where 00539 * <tt>nf=numFacets(cellDim,cellLID,facetDim)</tt>). 00540 * If <tt>facetDim==0</tt> this this array argument is ignored! 00541 * 00542 * The default implementation loops over calls to 00543 * <tt>facetLID()</tt>. Subclasses can provide a more efficient 00544 * implementation if desired. 00545 * 00546 * \todo Rename to something like getSingleCellFacetLIDs(...). 00547 * 00548 * \todo Change the facetOrientation argument to an <tt>Array<int>*</tt> 00549 * type so that it can be NULL and therefore ignored (which is a common use 00550 * case). 00551 */ 00552 void getFacetArray(int cellDim, int cellLID, int facetDim, 00553 Array<int>& facetLIDs, 00554 Array<int>& facetOrientations) const ; 00555 00556 /** \brief Return a view of an array of LIDs for a maximum-dimensional 00557 * cell's zero-dimensional facets (i.e. vertexes). 00558 * 00559 * \param maxCellLID 00560 * [in] Local ID of the maximum dimensional cell (i.e. element). 00561 * 00562 * <b>Preconditions:</b><ul> 00563 * <li>The cell <tt>maxCellLID</tt> must be a maximum-dimensional cell! 00564 * </ul> 00565 * 00566 * \returns A raw pointer into an array which stores the local process IDs 00567 * of the vertices. These IDs are the local process IDs, not the facet 00568 * indexes relative to the relative parent cell. The dimension of this 00569 * array is determined by the cell type given by 00570 * <tt>cellType(maxCellLID)</tt>. 00571 * 00572 * \todo Return this array as a Teuchos::ArrayView<int> object which will 00573 * have the dimension embedded in it and will have full range checking! 00574 * 00575 * \todo Rename to something like getMaxCellZeroFacetsLIDsView(...). 00576 */ 00577 virtual const int* elemZeroFacetView(int maxCellLID) const ; 00578 00579 /** \brief Return the number of maximal co-facets of the given cell. 00580 * 00581 * <b>Preconditions:</b><ul> 00582 * <li><tt>0 <= cellDim < spatialDim()</tt> 00583 * </ul> 00584 * 00585 * Note that if <tt>cellDim==spatialDim()-1</tt> and <tt>returnVal==1</tt> then 00586 * this cell must be a boundary cell (i.e. a boundary face in 3D, a boundary 00587 * edge in 1D, or a boundary node in 1D)! 00588 */ 00589 virtual int numMaxCofacets(int cellDim, int cellLID) const = 0; 00590 00591 /** \brief Return the LID of a maximal co-facet of a cell. 00592 * 00593 * \param cellDim 00594 * [in] Dimension of the cell whose co-facets are being obtained 00595 * \param cellLID 00596 * [in] Local ID of the cell whose co-facets are being obtained 00597 * \param cofacetIndex 00598 * [in] Index into the list of the cell's co-facets. 00599 * \param facetIndex 00600 * [out] Returns the local index into the facet array 00601 * for the relative cell (cellDim,cellLID) with respect 00602 * to the maximal co-facet. In other words 00603 * <tt>cellLID==facetLID(spatialDim(),returnVal,cellDim,facetIndex)</tt>. 00604 * 00605 * <b>Preconditions:</b><ul> 00606 * <li><tt>0 <= cellDim < spatialDim()</tt> 00607 * </ul> 00608 * 00609 * \returns The LID of the requested maximal co-facet. 00610 * 00611 * \todo Make the <tt>facetIndex</tt> an <tt>int*</tt> argument and give it a 00612 * default value of <tt>NUL</tt> so that it can be easily ignored! 00613 */ 00614 virtual int maxCofacetLID(int cellDim, int cellLID, 00615 int cofacetIndex, 00616 int& facetIndex) const = 0 ; 00617 /** 00618 * Get the LIDs of the maximal cofacets for a batch of codimension-one 00619 * cells. The default implementation simply loops over the cells in the 00620 * cellLID array, taking no advantage of any internal data structures. 00621 * 00622 * \param cellLIDs [in] array of LIDs of the cells whose cofacets are 00623 * being obtained 00624 * \param cofacets [out] 00625 * \param facetIndex [out] index of each calling cell 00626 * into the list of its maximal cofacet's facets 00627 */ 00628 virtual void getMaxCofacetLIDs(const Array<int>& cellLIDs, 00629 MaximalCofacetBatch& cofacets) const ; 00630 00631 /** \brief Return an array of the LIDs of all of the co-facets for a 00632 * given relative cell. 00633 * 00634 * \param cellDim 00635 * [in] Dimension of the relative cell whose co-facets are being obtained 00636 * \param cellLID 00637 * [in] Local index of the relative cell whose co-facets are being obtained 00638 * \param cofacetDim 00639 * [in] Dimension of the co-facets to get 00640 * \param cofacetLIDs 00641 * [out] Array containing the LIDs of the co-facets 00642 * for the given relative cell (cellDim,cellLID). 00643 * 00644 * \todo Change name to <tt>getCofacetArray()</tt> to be consistent with 00645 * <tt>getFacetArray()</tt>! 00646 */ 00647 virtual void getCofacets(int cellDim, int cellLID, 00648 int cofacetDim, Array<int>& cofacetLIDs) const = 0 ; 00649 00650 /** \brief Get the cell type of the given cell dimension. 00651 * 00652 * Note: This function by its very definition assumes that all cells of a 00653 * given dimension have the same cell type within a mesh! 00654 * 00655 * \todo This function must be changed in order to deal with mixed cell 00656 * types with the same cellDim! 00657 */ 00658 virtual CellType cellType(int cellDim) const = 0 ; 00659 00660 /** \brief Get the label of the given cell. 00661 * 00662 * \todo Change to getLabel(...)? 00663 */ 00664 virtual int label(int cellDim, int cellLID) const = 0 ; 00665 00666 /** \brief Get the labels for a batch of cells. 00667 * 00668 * \param cellDim 00669 * [in] Dimension of the parent cell whose facets 00670 * are being obtained 00671 * \param cellLIDs 00672 * [in] Array of cell LIDs 00673 * \param labels 00674 * [out] On exit, contains an array (<tt>size=cellLIDs.size()</tt>) 00675 * of the labels for each of the given cells. 00676 */ 00677 void getLabels(int cellDim, const Array<int>& cellLIDs, 00678 Array<int>& labels) const ; 00679 00680 /** \brief Set the label for the given cell. 00681 * 00682 * \todo Move this out of this base interface and into a mesh loading 00683 * interface? 00684 */ 00685 virtual void setLabel(int cellDim, int cellLID, int label) = 0 ; 00686 00687 /** Get the list of all labels defined for cells of the given dimension */ 00688 virtual Set<int> getAllLabelsForDimension(int cellDim) const ; 00689 00690 /** 00691 * Get the cells associated with a specified label. The array 00692 * cellLID will be filled with those cells of dimension cellDim 00693 * having the given label. 00694 */ 00695 virtual void getLIDsForLabel(int cellDim, int label, Array<int>& cellLIDs) const ; 00696 //@} 00697 00698 /** \name Geometric Information */ 00699 //@{ 00700 00701 /** \brief Return the position of a local vertex. 00702 * 00703 * \param vertexLID 00704 * [in] The LID of the vertex. 00705 * 00706 * <b>Preconditions:</b><ul> 00707 * <li><tt>0 <= vertexLID < this->numCells(0)</tt> 00708 * </ul> 00709 * 00710 * \todo Change the name of this function to 00711 * <tt>getVertexPosition(...)</tt>. 00712 */ 00713 virtual Point nodePosition(int vertexLID) const = 0 ; 00714 00715 /** \brief Return a const view into an raw array for the position of a local 00716 * vertex. 00717 * 00718 * \param vertexLID [in] The LID of the vertex 00719 * 00720 * <b>Preconditions:</b><ul> 00721 * <li><tt>0 <= vertexLID < this->numCells(0)</tt> 00722 * </ul> 00723 * 00724 * \returns a raw pointer into an array of doubles of length 00725 * <tt>spatialDim</tt>. 00726 * 00727 * \todo Return this array as a Teuchos::ArrayView<double> object which will 00728 * have the dimension embedded in it and will have full range checking! 00729 * 00730 * \todo Change this function name to <tt>getVertexPositionView()</tt>. 00731 */ 00732 virtual const double* nodePositionView(int vertexLID) const = 0 ; 00733 00734 /** \brief Return the centroid of a cell. 00735 * 00736 * The default implementation just averages the positions of the 00737 * zero-dimensional facets (i.e. vertexes). 00738 */ 00739 virtual Point centroid(int cellDim, int cellLID) const ; 00740 00741 /** \brief Compute (or get) the Jacobians for a batch of cells. 00742 * 00743 * \param cellDim 00744 * [in] Dimension of the cells whose Jacobians are to be computed 00745 * \param cellLID 00746 * [in] Local IDs of the cells for which Jacobians are to be computed 00747 * \param jBatch 00748 * [out] Batch of cell Jacobians. 00749 * 00750 * Warning! The default implementation returns an empty batch of cell 00751 * Jacobians! 00752 * 00753 * \todo Add a query function to tell if this feature is supported and then 00754 * add a precondition based on this query function! 00755 */ 00756 00757 virtual void getJacobians(int cellDim, const Array<int>& cellLID, 00758 CellJacobianBatch& jBatch) const { ; } 00759 00760 //bvbw virtual void getJacobians( 00761 // int cellDim, const Array<int>& cellLID, 00762 // CellJacobianBatch& jBatch 00763 // ) const 00764 // { jBatch.resize(0,0,0); } 00765 00766 /** \brief Compute the diameters of a batch of cells. 00767 * 00768 * \param cellDim 00769 * [in] Dimension of the cells whose diameters are to be computed 00770 * \param cellLIDs 00771 * [in] Local IDs of the cells for which diameters are to be computed 00772 * \param diameters 00773 * [out] Array (<tt>size = cellLIDs.size()</tt>) of cell diameters. 00774 * 00775 * Warning! The default implementation returns an empty array of cell 00776 * diameters! 00777 * ToDo: Change the default implementation to compute diameters based on 00778 * calls to the node position accessors. Going through the Mesh interface in 00779 * that way will be less efficient than a low-level implementation, but 00780 * would be a reasonable intermediate step for mesh developers. 00781 * - KL 5 Aug 2007. 00782 */ 00783 virtual void getCellDiameters( 00784 int cellDim, const Array<int>& cellLIDs, 00785 Array<double>& diameters 00786 ) const 00787 { diameters.resize(0); } 00788 00789 00790 /** 00791 * Get the outward normals for the batch of cells of dimension 00792 * spatialDim()-1. If any cell in the batch is not on the boundary, 00793 * an exception is thrown. 00794 * 00795 * \param cellLIDs [in] LIDs for the cells whose normals are to be 00796 * computed. 00797 * \param outwardNormals [out] Outward normal unit vectors for each 00798 * cell in the batch. 00799 */ 00800 virtual void outwardNormals( 00801 const Array<int>& cellLIDs, 00802 Array<Point>& outwardNormals 00803 ) const ; 00804 00805 /** 00806 * Get tangent vectors for a batch of edges 00807 * 00808 * \param cellLIDs [in] LIDs for the cells whose tangents are to be 00809 * computed. 00810 * \param tangentVectors [out] Unit tangents for each cell 00811 */ 00812 virtual void tangentsToEdges( 00813 const Array<int>& cellLIDs, 00814 Array<Point>& tangenVectors 00815 ) const ; 00816 00817 00818 /** \brief Map points from a reference cell to physical points for a batch 00819 * of cells. 00820 * 00821 * \param cellDim 00822 * [in] Dimension of the cells 00823 * \param cellLIDs 00824 * [in] Local IDs of a batch of cells 00825 * \param refPts 00826 * [in] Array of points on the single reference cell with respect 00827 * to the reference cell's coordinate system. Note that the 00828 * interpretation of these reference points is strictly determined 00829 * by the coordinate system of the cell type 00830 * <tt>cellType(cellDim)</tt> and must be clearly defined by this 00831 * interface. 00832 * \param physPts 00833 * [out] Array (<tt>size = cellLIDs.size()*refPts.size()</tt>) of 00834 * the physical points given in a flat array for the batch of 00835 * cells. Specifically, the physical points for each cell 00836 * <tt>cellLIDs[k]</tt>, for <tt>k=0...cellLIDs.size()-1</tt>, is 00837 * given by <tt>physPts[k*nrp+j]</tt>, for <tt>j=0...nrp-1</tt> 00838 * and <tt>nrp=refPts.size()</tt>. 00839 * 00840 * Warning! The default implementation returns an empty array of physical 00841 * points! 00842 * 00843 * \todo Add a query function to tell if this feature is supported and then 00844 * add a precondition based on this query function! 00845 */ 00846 virtual void pushForward( 00847 int cellDim, const Array<int>& cellLIDs, 00848 const Array<Point>& refPts, 00849 Array<Point>& physPts 00850 ) const 00851 { physPts.resize(0); } 00852 00853 //@} 00854 00855 /** \name Parallel Information */ 00856 //@{ 00857 00858 /** \brief Return the MPI communicator over which this mesh is distributed. */ 00859 const MPIComm& comm() const {return comm_;} 00860 00861 /** \brief Return the rank of the processor that owns the given cell. 00862 * 00863 * If <tt>returnVal==comm().getRank()</tt> then this cell is owned by this 00864 * process. 00865 */ 00866 virtual int ownerProcID(int cellDim, int cellLID) const = 0 ; 00867 00868 /** \brief Determine whether a given cell GID exists on this processor. */ 00869 virtual bool hasGID(int cellDim, int cellGID) const = 0 ; 00870 00871 /** \brief Find the LID of a cell given its GID. 00872 * 00873 * \param cellDim 00874 * [in] Dimension of the cell 00875 * \param cellGID 00876 * [in] Global ID of the cell 00877 * 00878 * <b>Preconditions:</b><ul> 00879 * <li>hasGID(cellDim,cellGID)==true</tt> 00880 * </ul> 00881 * 00882 * <b>Postconditions:</b><ul> 00883 * <li>0 <= returnVal < numCells(cellDim)</tt> 00884 * </ul> 00885 */ 00886 virtual int mapGIDToLID(int cellDim, int cellGID) const = 0 ; 00887 00888 /** \brief Find the global ID of a cell given its LID. 00889 * 00890 * \param cellDim 00891 * [in] Dimension of the cell 00892 * \param cellLID 00893 * [in] Local ID of the cell 00894 * 00895 * <b>Postconditions:</b><ul> 00896 * <li><tt>returnVal >= 0 </tt> 00897 * </ul> 00898 */ 00899 virtual int mapLIDToGID(int cellDim, int cellLID) const = 0 ; 00900 00901 /** \brief Return if cells of dimension cellDim have been assigned global 00902 * IDs or not. 00903 * 00904 * \param cellDim 00905 * [in] Dimension of the cell 00906 */ 00907 virtual bool hasIntermediateGIDs(int cellDim) const = 0 ; 00908 00909 /** \brief Assign global IDs to cells of dimension cellDim. 00910 * 00911 * \param cellDim 00912 * [in] Dimension of the cell 00913 * 00914 * <b>Postconditions:</b><ul> 00915 * <li><tt>hasIntermediateGIDs(cellDim)==true</tt> 00916 * </ul> 00917 */ 00918 virtual void assignIntermediateCellGIDs(int cellDim) = 0 ; 00919 00920 //@} 00921 00922 /** \name Reordering */ 00923 //@{ 00924 00925 /** \brief Set the reordering strategy to be used with this mesh. */ 00926 void setReorderer(const CellReorderer& reorderer) 00927 {reorderer_ = reorderer.createInstance(this);} 00928 00929 /** \brief Get the reordering strategy to be used with this mesh. */ 00930 const CellReordererImplemBase* reorderer() const 00931 {return reorderer_.get();} 00932 00933 //@} 00934 00935 00936 00937 /** \name Functions for Mesh with hanging nodes */ 00938 //@{ 00939 /** Function returns true if the mesh allows hanging nodes (by refinement), 00940 * false otherwise */ 00941 virtual bool allowsHangingHodes() const { return false; } 00942 00943 /** Function returns true if the specified element is a "hanging" element 00944 * false otherwise <br> 00945 * @param cellDim [in] should be between 0 , D-1 00946 * @param cellLID [in] the local ID of the element */ 00947 virtual bool isElementHangingNode(int cellDim , int cellLID) const { return false; } 00948 00949 /** Returns the index in the parent maxdim Cell of the refinement tree 00950 * @param maxCellLID [in] the LID of the cell */ 00951 virtual int indexInParent(int maxCellLID) const { return 0; } 00952 00953 /** How many children has a refined element. <br> 00954 * This function provides information of either we have bi or trisection */ 00955 virtual int maxChildren() const { return 0;} 00956 00957 /** Function returns the facets of the parent cell (needed for HN treatment) <br> 00958 * @param childCellLID [in] the LID of the maxdim cell, whos parents facets we want 00959 * @param dimFacets [in] the dimension of the facets which we want to have 00960 * @param facetsLIDs [out] the LID of the parents facets (all) in the defined order 00961 * @param parentCellLIDs [out] the maxdim parent cell LID */ 00962 virtual void returnParentFacets( int childCellLID , int dimFacets , 00963 Array<int> &facetsLIDs , int &parentCellLIDs ) const { } 00964 //@} 00965 00966 00967 00968 /** \name Store special weights in the mesh (for Adaptive Cell Integration) */ 00969 //@{ 00970 00971 /** returns the status of the special weights if they are valid <br> 00972 * These weights are usually computed for one setting of the curve (Adaptive Cell Integration)*/ 00973 virtual bool IsSpecialWeightValid() const {return validWeights_;} 00974 00975 /** specifies if the special weights are valid <br> 00976 * if this is false then usually the special weights have to be recomputed */ 00977 virtual void setSpecialWeightValid(bool val) const { validWeights_ = val;} 00978 00979 /** deletes all special weights so those have to be recreated*/ 00980 virtual void flushSpecialWeights() const; 00981 00982 /** verifies if the specified cell with the given dimension has special weights */ 00983 virtual bool hasSpecialWeight(int dim, int cellLID) const; 00984 00985 /** Sets the special weights */ 00986 virtual void setSpecialWeight(int dim, int cellLID, Array<double>& w) const; 00987 00988 /** Returns the special weights */ 00989 virtual void getSpecialWeight(int dim, int cellLID, Array<double>& w) const; 00990 //@} 00991 00992 00993 /** \name Store the intersection/quadrature points for the curve/surf integral <br> 00994 * for a curve or surf integral we need some quadrature points along the curve in one curve <br> 00995 * These */ 00996 //@{ 00997 00998 /** */ 00999 virtual bool IsCurvePointsValid() const {return curvePoints_Are_Valid_;} 01000 01001 /** */ 01002 virtual void setCurvePointsValid(bool val) const { curvePoints_Are_Valid_ = val;} 01003 01004 /** detletes all the points and its normals which have been stored */ 01005 virtual void flushCurvePoints() const; 01006 01007 /** verifies if the specified maxCell has already precalculated quadrature point for one curve */ 01008 virtual bool hasCurvePoints(int maxCellLID , int curveID) const; 01009 01010 /** Sets the points, curve derivatives and curve normals for one maxCell needed for curve/surf integral*/ 01011 virtual void setCurvePoints(int maxCellLID, int curveID , Array<Point>& points , Array<Point>& derivs , Array<Point>& normals) const; 01012 01013 /** Gets the points, curve derivatives and curve normals for one maxCell needed for curve/surf integral*/ 01014 virtual void getCurvePoints(int maxCellLID, int curveID , Array<Point>& points , Array<Point>& derivs , Array<Point>& normals) const; 01015 01016 private: 01017 virtual int mapCurveID_to_Index(int curveID) const; 01018 public: 01019 //@} 01020 01021 /** \name Output */ 01022 //@{ 01023 01024 /** \brief Set whether to stagger output in parallel. Set to true for readable 01025 * output in parallel debugging sessions. 01026 * 01027 * This should be normally, as it causes one synchronization point per 01028 * process. 01029 * 01030 * \todo Get rid of this once we update to use Teuchos::FancyOStream since 01031 * parallel outputting will be done in a readable way automatically! 01032 */ 01033 static bool& staggerOutput() {static bool rtn=false; return rtn;} 01034 01035 //@} 01036 01037 private: 01038 01039 int dim_; 01040 01041 MPIComm comm_; 01042 01043 MeshEntityOrder order_; 01044 01045 RCP<CellReordererImplemBase> reorderer_; 01046 01047 01048 01049 /** flag to indicate if the weights stored are valid */ 01050 mutable bool validWeights_; 01051 01052 /** Object to store the special weights for integration */ 01053 mutable Array < Sundance::Map< int , Array<double> > > specialWeights_; 01054 01055 01056 01057 /** true if the curve did not moved, false if those points are not reusable*/ 01058 mutable bool curvePoints_Are_Valid_; 01059 01060 /** how many curves are participating in curve integrals*/ 01061 mutable int nrCurvesForIntegral_; 01062 01063 /** store intersection informations to calculate the curve integral*/ 01064 mutable Array < Sundance::Map< int , Array<Point> > > curvePoints_; 01065 01066 /** store directional derivative informations to calculate the curve integral*/ 01067 mutable Array < Sundance::Map< int , Array<Point> > > curveDerivative_; 01068 01069 /** store normal directional used in the curve or in the surf integral <br> 01070 * in case of the surf integral it is the cross product from the integral */ 01071 mutable Array < Sundance::Map< int , Array<Point> > > curveNormal_; 01072 01073 /** map the curve ID to index in the previous arrays */ 01074 mutable Sundance::Map< int , int > curveID_to_ArrayIndex_; 01075 }; 01076 01077 01078 } // namespace Sundance 01079 01080 #endif