SundanceBasicSimplicialMesh.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                              Sundance
00005 //                 Copyright (2005) Sandia Corporation
00006 // 
00007 // Copyright (year first published) Sandia Corporation.  Under the terms 
00008 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government 
00009 // retains certain rights in this software.
00010 // 
00011 // This library is free software; you can redistribute it and/or modify
00012 // it under the terms of the GNU Lesser General Public License as
00013 // published by the Free Software Foundation; either version 2.1 of the
00014 // License, or (at your option) any later version.
00015 //  
00016 // This library is distributed in the hope that it will be useful, but
00017 // WITHOUT ANY WARRANTY; without even the implied warranty of
00018 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00019 // Lesser General Public License for more details.
00020 //                                                                                 
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License along with this library; if not, write to the Free Software
00023 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00024 // USA                                                                                
00025 // Questions? Contact Kevin Long (krlong@sandia.gov), 
00026 // Sandia National Laboratories, Livermore, California, USA
00027 // 
00028 // ************************************************************************
00029 /* @HEADER@ */
00030 
00031 #ifndef SUNDANCE_BASICSIMPLICIALMESH_H
00032 #define SUNDANCE_BASICSIMPLICIALMESH_H
00033 
00034 
00035 
00036 #include "SundanceDefs.hpp"
00037 #include "SundanceMeshBase.hpp"
00038 #include "SundanceSet.hpp"
00039 #include "SundanceBasicVertexView.hpp"
00040 #include "SundanceArrayOfTuples.hpp"
00041 #include "SundanceIncrementallyCreatableMesh.hpp"
00042 #include "Teuchos_Array.hpp"
00043 #include "Teuchos_Hashtable.hpp"
00044 
00045 namespace Sundance
00046 {
00047 
00048 /**
00049  * A no-frills parallel simplicial mesh
00050  */
00051 class BasicSimplicialMesh : public IncrementallyCreatableMesh
00052 {
00053 public:
00054   /** */
00055   BasicSimplicialMesh(int dim, const MPIComm& comm, 
00056     const MeshEntityOrder& order);
00057 
00058   /** */
00059   virtual ~BasicSimplicialMesh(){;}
00060 
00061   /** 
00062    * Get the number of cells having dimension dim
00063    */
00064   virtual int numCells(int dim) const  ;
00065 
00066   /** 
00067    * Return the position of the i-th node
00068    */
00069   virtual Point nodePosition(int i) const {return points_[i];}
00070 
00071   /** 
00072    * Return a view of the i-th node's position
00073    */
00074   const double* nodePositionView(int i) const {return &(points_[i][0]);}
00075 
00076      
00077 
00078   /** 
00079    * Compute the jacobians of a batch of cells, returning the 
00080    * result via reference argument
00081    *
00082    * @param cellDim dimension of the cells whose Jacobians are to
00083    * be computed
00084    * @param cellLID local indices of the cells for which Jacobians
00085    * are to be computed
00086    * @param jBatch reference to the resulting Jacobian batch
00087    */
00088   virtual void getJacobians(int cellDim, const Array<int>& cellLID,
00089     CellJacobianBatch& jBatch) const  ;
00090 
00091   /** 
00092    * Compute the diameters of a batch of cells,
00093    * result via reference argument
00094    *
00095    * @param cellDim dimension of the cells whose diameters are to
00096    * be computed
00097    * @param cellLID local indices of the cells for which diameters
00098    * are to be computed
00099    * @param diameters reference to the array of cell diameters
00100    */
00101   virtual void getCellDiameters(int cellDim, const Array<int>& cellLID,
00102     Array<double>& diameters) const ;
00103 
00104 
00105   /**
00106    * Map reference quadrature points to physical points on the
00107    * given cells. 
00108    */
00109   virtual void pushForward(int cellDim, const Array<int>& cellLID,
00110     const Array<Point>& refQuadPts,
00111     Array<Point>& physQuadPts) const ;
00112 
00113       
00114 
00115   /** 
00116    * Return the rank of the processor that owns the given cell
00117    */
00118   virtual int ownerProcID(int cellDim, int cellLID) const  ;
00119 
00120       
00121 
00122   /** 
00123    * Return the number of facets of the given cell
00124    */
00125   virtual int numFacets(int cellDim, int cellLID, 
00126     int facetDim) const  ;
00127 
00128   /** 
00129    * Return the local ID of a facet cell
00130    * @param cellDim dimension of the cell whose facets are being obtained
00131    * @param cellLID local index of the cell whose
00132    * facets are being obtained
00133    * @param facetDim dimension of the desired facet
00134    * @param facetIndex index into the list of the cell's facets
00135    * @param facetOrientation orientation of the facet as seen from
00136    * the given cell (returned via reference)
00137    */
00138   virtual int facetLID(int cellDim, int cellLID,
00139     int facetDim, int facetIndex,
00140     int& facetOrientation) const  ;
00141 
00142 
00143   /** 
00144    * Return by reference argument an array containing
00145    * the LIDs of the facetDim-dimensional facets of the
00146    * given batch of cells 
00147    */
00148   virtual void getFacetLIDs(int cellDim, 
00149     const Array<int>& cellLID,
00150     int facetDim,
00151     Array<int>& facetLID,
00152     Array<int>& facetOrientations) const ;
00153 
00154   /** 
00155    * Return a view of an element's zero-dimensional facets
00156    */
00157   const int* elemZeroFacetView(int cellLID) const 
00158     {return &(elemVerts_.value(cellLID, 0));}
00159 
00160   /** 
00161    * Return the number of maximal cofacets of the given cell
00162    */
00163   virtual int numMaxCofacets(int cellDim, int cellLID) const  ;
00164 
00165   /** 
00166    * Return the local ID of a maximal cofacet cell
00167    * @param cellDim dimension of the cell whose cofacets are being obtained
00168    * @param cellLID local index of the cell whose
00169    * cofacets are being obtained
00170    * @param cofacetIndex which maximal cofacet to get
00171    * @param cofacetIndex index of the cell cellLID into the list of the 
00172    * maximal cell's facets
00173    */
00174   virtual int maxCofacetLID(int cellDim, int cellLID,
00175     int cofacetIndex,
00176     int& facetIndex) const  ;
00177   /** 
00178    * Get the LIDs of the maximal cofacets for a batch of codimension-one
00179    * cells.
00180    *
00181    * \param cellLIDs [in] array of LIDs of the cells whose cofacets are 
00182    * being obtained
00183    * \param cofacets [out] the batch of cofacets
00184    */
00185   virtual void getMaxCofacetLIDs(const Array<int>& cellLIDs,
00186     MaximalCofacetBatch& cofacets) const ;
00187         
00188       
00189   /** 
00190    * Find the cofacets of the given cell
00191    * @param cellDim dimension of the cell whose cofacets are being obtained
00192    * @param cellLID local index of the cell whose
00193    * cofacets are being obtained
00194    * @param cofacetDim dimension of the cofacets to get
00195    * @param cofacetLIDs LIDs of the cofacet
00196    */
00197   void getCofacets(int cellDim, int cellLID,
00198     int cofacetDim, Array<int>& cofacetLIDs) const ;
00199       
00200   /** 
00201    * Find the local ID of a cell given its global index
00202    */
00203   virtual int mapGIDToLID(int cellDim, int globalIndex) const  ;
00204 
00205   /** 
00206    * Determine whether a given cell GID exists on this processor
00207    */
00208   virtual bool hasGID(int cellDim, int globalIndex) const ;
00209 
00210     
00211 
00212   /** 
00213    * Find the global ID of a cell given its local index
00214    */
00215   virtual int mapLIDToGID(int cellDim, int localIndex) const  ;
00216 
00217   /**
00218    * Get the type of the given cell
00219    */
00220   virtual CellType cellType(int cellDim) const  ;
00221 
00222 
00223   /** Get the label of the given cell */
00224   virtual int label(int cellDim, int cellLID) const ;
00225 
00226   /** Get the labels for a batch of cells */
00227   virtual void getLabels(int cellDim, const Array<int>& cellLID, 
00228     Array<int>& labels) const ;
00229 
00230       
00231 
00232   /** Get the list of all labels defined for cells of the given dimension */
00233   virtual Set<int> getAllLabelsForDimension(int cellDim) const ;
00234 
00235   /** 
00236    * Get the cells associated with a specified label. The array 
00237    * cellLID will be filled with those cells of dimension cellDim
00238    * having the given label.
00239    */
00240   virtual void getLIDsForLabel(int cellDim, int label, Array<int>& cellLIDs) const ;
00241 
00242   /** \name Incremental creation methods */
00243   //@{
00244   /** 
00245    * Add new new vertex to the mesh.
00246    * \param globalIndex the GID of the new vertex
00247    * \param x the spatial position of the new vertex
00248    * \param ownerProcID the processor that "owns" this vertex 
00249    * \param label a label for this vertex (optionally used in setting 
00250    * loads, boundary conditions, etc)
00251    * \return the LID of the vertex.
00252    */
00253   virtual int addVertex(int globalIndex, const Point& x,
00254     int ownerProcID, int label);
00255 
00256   /** 
00257    * Add a new element to the mesh.
00258    * \param globalIndex the GID of the new element
00259    * \param vertexGIDs tuple of GIDs for the vertices defining this element. 
00260    * \param ownerProcID the processor that "owns" this element
00261    * \param label a label for this element (optionally used in setting loads, 
00262    * material properties, etc)
00263    * \return the LID of the element
00264    */
00265   virtual int addElement(int globalIndex, const Array<int>& vertexGIDs,
00266     int ownerProcID, int label);
00267 
00268   /** Set the label of the given cell */
00269   virtual void setLabel(int cellDim, int cellLID, int label)
00270     {
00271       labels_[cellDim][cellLID] = label;
00272     }
00273 
00274   /** Optional preallocation of space for an estimated number of vertices */
00275   virtual void estimateNumVertices(int numVertices);
00276 
00277   /** Optional preallocation of space for an estimated number of elements */
00278   virtual void estimateNumElements(int numElements);
00279 
00280   /** Coordinate intermediate cell definitions across processors  */
00281   virtual void assignIntermediateCellGIDs(int cellDim) ;
00282 
00283   /** */
00284   virtual bool hasIntermediateGIDs(int dim) const 
00285     {
00286       if (dim==1) return hasEdgeGIDs_;
00287       return hasFaceGIDs_;
00288     }
00289       
00290   //@}
00291 private:
00292 
00293   /** 
00294    * Add a new face, first checking to see if it already exists. 
00295    * This function is called from within addElement(),
00296    * not by the user, and is therefore private.
00297    *
00298    * \param vertLID{1,2,3} The LIDs for the three vertices of the face 
00299    * \param edgeLID{1,2,3} The LIDs for the three edges of the face 
00300    * \param rotation An integer specifying the permutation that transforms
00301    * the sorted ordering of vertGID{1,2,3} to the order given by the
00302    * element definition. Initial value is not used; value is
00303    * returned by reference argument.
00304    * \return LID of this face
00305    */
00306   int addFace(int vertLID1, int vertLID2, int vertLID3,
00307     int edgeLID1, int edgeLID2, int edgeLID3,
00308     int elemLID,
00309     int elemGID,
00310     int& rotation);
00311 
00312   /**
00313    * Add a new edge, first checking to see if it already exists. 
00314    * This function is called from within addElement(),
00315    * not by the user, and is therefore private.
00316    *
00317    * \param vertLID{1,2} 
00318    * \param elemLID LID of the element that is adding the edge
00319    * \param elemGID GID of the element that is adding the edge
00320    * \param myFacetNumber facet number of the edge within the element
00321    * \return LID of this edge
00322    */
00323   int addEdge(int vertLID1, int vertLID2, int elemLID, int elemGID, 
00324     int myFacetNumber);
00325 
00326   /** 
00327    * Check for the presence of the edge (vertLID1, vertLID2) in the mesh.
00328    * \return the edge LID if the edge exists, -1 otherwise
00329    */
00330   int checkForExistingEdge(int vertLID1, int vertLID2);
00331 
00332   /** 
00333    * Check for the presence of the face (vertLID1, int vertLID2, int vertLID3),
00334    * returning information on orientation of the face as defined by the calling
00335    * element. 
00336    */
00337   int checkForExistingFace(int vertLID1, int vertLID2, int vertLID3,
00338     int edgeLID1, int edgeLID2, int edgeLID3,
00339     int* sortedVertGIDs,
00340     int* reorderedVertLIDs,
00341     int* reorderedEdgeLIDs,
00342     int& rotation) ;
00343 
00344   /** 
00345    * Check whether the face defined by the given vertices exists
00346    * in the mesh. Returns -1 if the face does not exist.
00347    * Called during the synchronization of intermediate cell GIDs.
00348    * \param vertGID{1,2,3} the global indices of the vertices defining the face
00349    * \return the LID of the face
00350    */
00351   int lookupFace(int vertGID1, int vertGID2, int vertGID3) ;
00352 
00353   /** */
00354   void synchronizeGIDNumbering(int dim, int localCount) ;
00355       
00356   /** */
00357   void resolveEdgeOwnership(int cellDim);
00358 
00359   /** */
00360   std::string cellStr(int dim, const int* verts) const ;
00361 
00362   /** */
00363   std::string cellToStr(int dim, int cellLID) const ;
00364 
00365   /** */
00366   std::string printCells(int dim) const ;
00367 
00368   /** */
00369   void synchronizeNeighborLists();
00370       
00371       
00372 
00373   /** Number of cells of each dimension. */
00374   Array<int> numCells_;
00375     
00376   /** coordinates of vertices. The index into the array is the vertex LID .*/
00377   Array<Point> points_;
00378 
00379   /** pairs of local vertex indices for the edges, each pair ordered
00380    * from lower to higher <i>global</i> vertex index in order to define
00381    * an absolute edge orientation. Because global vertex indices are used, all
00382    * processors will agree on this orientation, regardless of the orientation
00383    * of the edge as seen by the element first defining it. 
00384    * The first index into this 2D array is the edge LID, the second the vertex 
00385    * number within the edge. 
00386    * */
00387   ArrayOfTuples<int> edgeVerts_;
00388 
00389   /** Tuples of local vertex indices for the faces, with each tuple ordered from
00390    * lowest to highest <i>global</i> index in order to define an absolute edge
00391    * orientation. Because global vertex indices are used, all
00392    * processors will agree on this orientation, regardless of the orientation
00393    * of the face as seen by the element first defining it. 
00394    * The first index into this 2D array is the face LID, the second the
00395    * vertex number within the face. */
00396   ArrayOfTuples<int> faceVertLIDs_;
00397 
00398   /** Tuples of global vertex indices for the faces, with each tuple ordered from
00399    * lowest to highest <i>global</i> index in order to define an absolute edge
00400    * orientation. Because global vertex indices are used, all
00401    * processors will agree on this orientation, regardless of the orientation
00402    * of the face as seen by the element first defining it. 
00403    * The first index into this 2D array is the face LID, the second the
00404    * vertex number within the face. 
00405    *
00406    * Notice that we duplicate the face vertex storage, storing both the
00407    * vertex LIDs and vertex GIDs for each face. This lets us do quick comparison
00408    * with the sorted GID array in order to identify pre-existing faces, while
00409    * also making it possible to retrieve face vertex LID information without
00410    * doing hashtable lookups.
00411    *
00412    */
00413   ArrayOfTuples<int> faceVertGIDs_;
00414     
00415   /** Tuples of local indices for the edges of all faces. The first index
00416    * into this 2D array is the face LID, the second the edge number. */
00417   ArrayOfTuples<int> faceEdges_;
00418     
00419   /** Tuples of edge signs for the faces. The edge sign indicates 
00420    * whether the orientation of the edge as given by moving around the face 
00421    * is parallel or antiparallel to the absolute orientation of the edge. */
00422   ArrayOfTuples<int> faceEdgeSigns_;
00423     
00424   /** tuples of local vertex indices for the elements. The first index into this
00425    * 2D array is the element LID, the second is the vertex number.  */
00426   ArrayOfTuples<int> elemVerts_;
00427 
00428   /** tuples of local edge indices for the elements. The first index into
00429    * this 2D array is the element LID, the second is the edge number. */
00430   ArrayOfTuples<int> elemEdges_;
00431 
00432   /** tuples of edge orientations for the elements, indicating whether 
00433    * the orientation of the edge as given by moving around the element 
00434    * is parallel or antiparallel to the absolute orientation of the edge. 
00435    * The first index into this 2D array is the element LID, the second 
00436    * the edge number. 
00437    * */
00438   ArrayOfTuples<int> elemEdgeSigns_;
00439 
00440   /** tuples of face LIDs for the elements. The first index is the
00441    * element LID, the second the face number. */
00442   ArrayOfTuples<int> elemFaces_;
00443 
00444   /** tuples of face rotations for the elements, defined relative to the 
00445    * absolute orientation of the face. */
00446   ArrayOfTuples<int> elemFaceRotations_;
00447 
00448   /** table for mapping vertex set -> face index */
00449   Hashtable<VertexView, int> vertexSetToFaceIndexMap_;
00450 
00451   /** array of face cofacets for the edges. The first index
00452    * is the edge LID, the second the cofacet number. */
00453   Array<Array<int> > edgeFaces_;
00454 
00455   /** array of element cofacets for the edges. The first index
00456    * is the edge LID, the second the cofacet number. */
00457   Array<Array<int> > edgeCofacets_;
00458 
00459   /** array of element cofacets for the faces. The first index is the
00460    * face LID, the second the cofacet number. */
00461   Array<Array<int> > faceCofacets_;
00462 
00463   /** array of edge cofacets for the vertices. The first index is the 
00464    * vertex LID, the second the edge cofacet number. */
00465   Array<Array<int> > vertEdges_;
00466 
00467   /** array of face cofacet LIDs for the vertices. The first index is the 
00468    * vertex LID, the second the cofacet number. */
00469   Array<Array<int> > vertFaces_;
00470 
00471   /** array of maximal cofacets for the vertices. The first index is the
00472    * vertex LID, the second the cafacet number. */
00473   Array<Array<int> > vertCofacets_;
00474 
00475   /** array of edge partners for the vertices. The partners are other
00476    * vertices sharing an edge with the specified vertex. */
00477   Array<Array<int> > vertEdgePartners_;
00478 
00479   /** map from local to global cell indices. The first index into this
00480    * 2D array is the cell dimension, the second the cell LID. */
00481   Array<Array<int> > LIDToGIDMap_;
00482 
00483   /** map from global to local cell indices. The array index is the 
00484    * cell dimension. The hashtable key is the cell GID, the value the
00485    * cell LID. */
00486   Array<Hashtable<int, int> > GIDToLIDMap_;
00487     
00488   /** Array of labels for the cells */
00489   Array<Array<int> > labels_;
00490 
00491   /** Array of owning processor IDs for the cells. Each cell is owned by
00492    * a single processor that is responsible for assigning global indices,
00493    * DOF numbers, and so on. */
00494   Array<Array<int> > ownerProcID_;
00495       
00496   /** 
00497    * Pointer to the pointer at the base of the face vertex GID array. This is
00498    * used to get small-array "views" of the face vertex GID array without
00499    * making copies, resulting in a significant performance improvement
00500    * in the vertex set hashtable lookups to identify pre-existing faces. 
00501    *
00502    * We use double rather than single indirection here because as elements
00503    * are added, the face vertex GID array will often be resized, thus changing
00504    * the base pointer. Each vertex set "view" keeps a pointer to the base pointer,
00505    * so that it always remains synchronized with the array of face vertices. 
00506    * 
00507    * IMPORTANT: any time faceVertGIDs_ is resized, faceVertGIDBase_[0] must
00508    * be reset to the base of the faceVertGIDs_ array so that the vertex
00509    * sets are pointing to the right place.
00510    */
00511   Array<int*> faceVertGIDBase_;
00512 
00513 
00514   /** flag indicating whether the edge GIDs have been synchronized */
00515   bool hasEdgeGIDs_;
00516 
00517   /** flag indicating whether the face GIDs have been synchronized */
00518   bool hasFaceGIDs_;
00519 
00520   /** Set of all neighbor processors sharing data with this one */
00521   Set<int> neighbors_;
00522 
00523   /** Whether the neighbor lists have been synchronized across 
00524    * processors */
00525   bool neighborsAreSynchronized_;
00526 
00527 
00528   /** 
00529    * Sort the 3 vertices of a face in order of increasing GID
00530    * \param vertGID{1,2,3} the three vertices
00531    * \param sortedVertGIDs array in which the sorted vertices are to be returned
00532    */
00533   void getSortedFaceVertices(int vertGID1, int vertGID2, int vertGID3,
00534     int* sortedVertGIDs) const ;
00535 
00536 
00537   /** 
00538    * Sort the the vertices of a face in order or increasing GID, reordering
00539    * the LIDs of the vertices and associated edges following the GID ordering.
00540    * Also, return a rotation flag describing the permutation between
00541    * the unsorted and sorted arrays. 
00542    */
00543   void getSortedFaceVertices(int a, int b, int c,
00544     int la, int lb, int lc,
00545     int eAB, int eBC, int eCA, 
00546     int* sortedVertGIDs,
00547     int* reorderedVertLIDs,
00548     int* reorderedEdgeLIDs,
00549     int& rotation) const ;
00550 
00551   /** Helper function to stuff three numbers into an array */
00552   inline void fillSortedArray(int a, int b, int c, int* s) const
00553     {
00554       s[0] = a;
00555       s[1] = b;
00556       s[2] = c;
00557     }
00558 
00559 
00560 };
00561 
00562 }
00563 
00564 
00565 
00566 
00567 #endif

Site Contact