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 * SundanceHNMesh2D.h 00043 * 00044 * Created on: April 30, 2010 00045 * Author: benk 00046 */ 00047 00048 #ifndef SUNDANCEHNMESH2D_H_ 00049 #define SUNDANCEHNMESH2D_H_ 00050 00051 #include "SundanceDefs.hpp" 00052 #include "SundanceMeshBase.hpp" 00053 #include "SundanceSet.hpp" 00054 #include "SundancePoint.hpp" 00055 #include "SundanceCellType.hpp" 00056 00057 #include "Teuchos_Array.hpp" 00058 #include "Teuchos_Hashtable.hpp" 00059 00060 #include "SundanceRefinementBase.hpp" 00061 #include "SundanceRefinementClass.hpp" 00062 00063 #include "SundanceDomainDefinition.hpp" 00064 00065 00066 namespace Sundance 00067 { 00068 using namespace Sundance; 00069 using namespace Teuchos; 00070 00071 /** 00072 * Class for 2D hierarchical structured quad Mesh <br> 00073 * The main internal idea of the mesh ist that there are different numberings of the elements 00074 * ID -> each element has an ID (also those which are completely outside the mesh domain, and those 00075 * which are not leaf) <br> 00076 * In the code sometimes internaly LID is used instead of ID due to legacy !! :-P <br> 00077 * GID -> all leaf elements which are visible to Sundance (only the leaf elements should be visible) <br> 00078 * LID -> all the elements which have GID and belong and either belong to the local processor or are in the ghost cells <br>*/ 00079 00080 class HNMesh2D : public MeshBase{ 00081 00082 public: 00083 /** 00084 * The Ctor for the dummy grid with hanging nodes */ 00085 HNMesh2D(int dim, 00086 const MPIComm& comm, 00087 const MeshEntityOrder& order); 00088 00089 /** 00090 * The Ctor for the HNMesh2D grid in 2D*/ 00091 void createMesh( 00092 double position_x, 00093 double position_y, 00094 double offset_x, 00095 double offset_y, 00096 int resolution_x, 00097 int resolution_y, 00098 const RefinementClass& refineClass, 00099 const MeshDomainDef& meshDomain ); 00100 00101 /** Dtor */ 00102 virtual ~HNMesh2D(){;} 00103 00104 /** 00105 * Get the number of cells having dimension dim 00106 */ 00107 virtual int numCells(int dim) const; 00108 00109 /** 00110 * Return the position of the i-th node 00111 */ 00112 virtual Point nodePosition(int i) const; 00113 00114 /** 00115 * Return a view of the i-th node's position 00116 */ 00117 const double* nodePositionView(int i) const; 00118 00119 00120 00121 /** 00122 * Compute the jacobians of a batch of cells, returning the 00123 * result via reference argument 00124 * 00125 * @param cellDim dimension of the cells whose Jacobians are to 00126 * be computed 00127 * @param cellLID local indices of the cells for which Jacobians 00128 * are to be computed 00129 * @param jBatch reference to the resulting Jacobian batch 00130 */ 00131 virtual void getJacobians(int cellDim, const Array<int>& cellLID, 00132 CellJacobianBatch& jBatch) const; 00133 00134 /** 00135 * Compute the diameters of a batch of cells, 00136 * result via reference argument 00137 * 00138 * @param cellDim dimension of the cells whose diameters are to 00139 * be computed 00140 * @param cellLID local indices of the cells for which diameters 00141 * are to be computed 00142 * @param diameters reference to the array of cell diameters 00143 */ 00144 virtual void getCellDiameters(int cellDim, const Array<int>& cellLID, 00145 Array<double>& cellDiameters) const; 00146 00147 00148 /** 00149 * Map reference quadrature points to physical points on the 00150 * given cells. 00151 */ 00152 virtual void pushForward(int cellDim, const Array<int>& cellLID, 00153 const Array<Point>& refQuadPts, 00154 Array<Point>& physQuadPts) const; 00155 00156 /** 00157 * Return the rank of the processor that owns the given cell 00158 */ 00159 virtual int ownerProcID(int cellDim, int cellLID) const; 00160 00161 /** 00162 * Return the number of facets of the given cell 00163 */ 00164 virtual int numFacets(int cellDim, int cellLID, 00165 int facetDim) const; 00166 00167 /** 00168 * Return the local ID of a facet cell 00169 * @param cellDim dimension of the cell whose facets are being obtained 00170 * @param cellLID local index of the cell whose 00171 * facets are being obtained 00172 * @param facetDim dimension of the desired facet 00173 * @param facetIndex index into the list of the cell's facets 00174 * @param facetOrientation orientation of the facet as seen from 00175 * the given cell (returned via reference) 00176 */ 00177 virtual int facetLID(int cellDim, int cellLID, 00178 int facetDim, int facetIndex, 00179 int& facetOrientation) const; 00180 00181 /** 00182 * Return by reference argument an array containing&(elemVerts_.value(cellLID, 0)) 00183 * the LIDs of the facetDim-dimensional facets of the 00184 * given batch of cells 00185 */ 00186 virtual void getFacetLIDs(int cellDim, 00187 const Array<int>& cellLID, 00188 int facetDim, 00189 Array<int>& facetLID, 00190 Array<int>& facetSign) const; 00191 00192 /** 00193 * Return a view of an element's zero-dimensional facets, 00194 * @return an array of integers with the indexes of the points which for it 00195 */ 00196 const int* elemZeroFacetView(int cellLID) const; 00197 00198 /** 00199 * Return the number of maximal cofacets of the given cell 00200 */ 00201 virtual int numMaxCofacets(int cellDim, int cellLID) const; 00202 00203 /** 00204 * Return the local ID of a maximal cofacet cell 00205 * @param cellDim dimension of the cell whose cofacets are being obtained 00206 * @param cellLID local index of the cell whose 00207 * cofacets are being obtained 00208 * @param cofacetIndex which maximal cofacet to get 00209 * @param facetIndex index of the cell cellLID into the list of the 00210 * maximal cell's facets 00211 */ 00212 virtual int maxCofacetLID(int cellDim, int cellLID, 00213 int cofacetIndex, 00214 int& facetIndex) const; 00215 /** 00216 * Get the LIDs of the maximal cofacets for a batch of codimension-one 00217 * cells. 00218 * 00219 * \param cellLIDs [in] array of LIDs of the cells whose cofacets are 00220 * being obtained 00221 * \param cofacets [out] the batch of cofacets 00222 */ 00223 virtual void getMaxCofacetLIDs(const Array<int>& cellLIDs, 00224 MaximalCofacetBatch& cofacets) const; 00225 00226 00227 /** 00228 * Find the cofacets of the given cell 00229 * @param cellDim dimension of the cell whose cofacets are being obtained 00230 * @param cellLID local index of the cell whose 00231 * cofacets are being obtained 00232 * @param cofacetDim dimension of the cofacets to get 00233 * @param cofacetLIDs LIDs of the cofacet 00234 */ 00235 void getCofacets(int cellDim, int cellLID, 00236 int cofacetDim, Array<int>& cofacetLIDs) const; 00237 00238 /** 00239 * Find the local ID of a cell given its global index 00240 */ 00241 virtual int mapGIDToLID(int cellDim, int globalIndex) const; 00242 00243 /** 00244 * Determine whether a given cell GID exists on this processor 00245 */ 00246 virtual bool hasGID(int cellDim, int globalIndex) const; 00247 00248 00249 /** 00250 * Find the global ID of a cell given its local index 00251 */ 00252 virtual int mapLIDToGID(int cellDim, int localIndex) const; 00253 00254 /** 00255 * Get the type of the given cell 00256 */ 00257 virtual CellType cellType(int cellDim) const; 00258 00259 00260 /** Get the label of the given cell */ 00261 virtual int label(int cellDim, int cellLID) const; 00262 00263 /** Get the labels for a batch of cells */ 00264 virtual void getLabels(int cellDim, const Array<int>& cellLID, 00265 Array<int>& labels) const; 00266 00267 00268 00269 /** Get the list of all labels defined for cells of the given dimension */ 00270 virtual Set<int> getAllLabelsForDimension(int cellDim) const; 00271 00272 /** 00273 * Get the cells associated with a specified label. The array 00274 * cellLID will be filled with those cells of dimension cellDim 00275 * having the given label. 00276 */ 00277 virtual void getLIDsForLabel(int cellDim, int label, Array<int>& cellLIDs) const; 00278 00279 00280 /** Set the label of the given cell */ 00281 virtual void setLabel(int cellDim, int cellLID, int label); 00282 00283 /** Coordinate intermediate cell definitions across processors */ 00284 virtual void assignIntermediateCellGIDs(int cellDim); 00285 00286 /** */ 00287 virtual bool hasIntermediateGIDs(int dim) const; 00288 00289 00290 /** Function returns true if the mesh allows hanging nodes (by refinement), 00291 * false otherwise */ 00292 virtual bool allowsHangingHodes() const { return true; } 00293 00294 /** Function returns true if the specified element is a "hanging" element 00295 * false otherwise <br> 00296 * @param dim [in] should be between 0 , D-1 00297 * @param cellLID [in] the local ID of the element */ 00298 virtual bool isElementHangingNode(int cellDim , int cellLID) const ; 00299 00300 /** Returns the index in the parent maxdim Cell of the refinement tree 00301 * @param maxCellLID [in] the LID of the cell */ 00302 virtual int indexInParent(int maxCellLID) const ; 00303 00304 /** How many children has a refined element. <br> 00305 * This function provides information of either we have bi or trisection */ 00306 virtual int maxChildren() const {return 9;} 00307 00308 /** Function returns the facets of the maxdim parent cell (needed for HN treatment) <br> 00309 * @param childCellLID [in] the LID of the maxdim cell, whos parents facets we want 00310 * @param dimFacets [in] the dimension of the facets which we want to have 00311 * @param facetsLIDs [out] the LID of the parents facets (all) in the defined order 00312 * @param parentCellLIDs [out] the maxdim parent cell LID */ 00313 virtual void returnParentFacets( int childCellLID , int dimFacets , 00314 Array<int> &facetsLIDs , int &parentCellLIDs ) const; 00315 00316 private: 00317 00318 /** For HN , returns parent facets, if the facet is not leaf, then return -1 at that place */ 00319 int facetLID_tree(int cellDim, int cellLID, 00320 int facetDim, int facetIndex) const; 00321 00322 /** adds one vertex to the mesh */ 00323 void addVertex(int vertexLID , int ownerProc , bool isHanging , 00324 double coordx , double coordy , const Array<int> &maxCoF); 00325 00326 /** adds one edge to the mesh */ 00327 void addEdge(int edgeLID , int ownerProc , bool isHanging , int edgeVertex , 00328 bool isProcBoundary , bool isMeshBoundary , 00329 const Array<int> &vertexLIDs , const Array<int> &maxCoF); 00330 00331 /** adds one cell(2D) to the mesh <br> 00332 * cell must be always leaf*/ 00333 void addCell(int cellLID , int ownerProc , 00334 int indexInParent , int parentCellLID , int level , 00335 const Array<int> &edgeLIDs , const Array<int> &vertexLIDs); 00336 00337 /** creates the mesh on the coarsest level as it is specified */ 00338 void createCoarseMesh(); 00339 00340 /** Does one single refinement iteration. <br> 00341 * Iterates trough all cells which are owned by the processor and refines if necessary */ 00342 bool oneRefinementIteration(); 00343 00344 /** refine the given cell by cellID, (method assumes that this cell can be refined)*/ 00345 void refineCell(int cellLID); 00346 00347 /** Create Leaf numbering */ 00348 void createLeafNumbering(); 00349 00350 /** Create Leaf numbering, but with a better load balancing for parallel case */ 00351 void createLeafNumbering_sophisticated(); 00352 00353 /** estimate the load of one cell*/ 00354 int estimateCellLoad(int ID); 00355 00356 /** marks the cells recursivly in the tree (and facets) owned by one processor*/ 00357 void markCellsAndFacets(int cellID , int procID); 00358 00359 /** The dimension of the grid*/ 00360 int _dimension; 00361 /** Number of processors */ 00362 int nrProc_; 00363 int myRank_; 00364 /** The communicator */ 00365 const MPIComm& _comm; 00366 /** */ 00367 double _pos_x; 00368 /** */ 00369 double _pos_y; 00370 /** */ 00371 double _ofs_x; 00372 /** */ 00373 double _ofs_y; 00374 /** */ 00375 int _res_x; 00376 /** */ 00377 int _res_y; 00378 /** */ 00379 mutable RefinementClass refineClass_; 00380 /** */ 00381 mutable MeshDomainDef meshDomain_; 00382 00383 //------ Point storage ---- 00384 /** all the global points index is [ID] */ 00385 Array<Point> points_; 00386 /** [3] the nr of ID per dim*/ 00387 Array<int> nrElem_; 00388 /** [3] the nr of owned elements per dim*/ 00389 Array<int> nrElemOwned_; 00390 00391 //----- Facets ----- ; -1 -> MeshBoundary , -2 -> ProcBoundary 00392 00393 /** [cellID][4]*/ 00394 Array< Array<int> > cellsPoints_; 00395 /** [cellID][4]*/ 00396 Array< Array<int> > cellsEdges_; 00397 /** [edgeID][2]*/ 00398 Array< Array<int> > edgePoints_; 00399 /** Information from the edge needs to be stored in one vertex <br> 00400 * We use the traditional mapping , information is put into the 0th Vertex of the Edge <br>*/ 00401 Array<int> edgeVertex_; 00402 00403 // ----- MaxCofacets ----; -1 -> MeshBoundary , -2 -> ProcBoundary 00404 00405 /** [edgeID][2] */ 00406 Array< Array<int> > edgeMaxCoF_; 00407 /** [pointID][4]*/ 00408 Array< Array<int> > pointMaxCoF_; 00409 00410 // ---- Different mesh boundaries ------ 00411 Array<bool> edgeIsProcBonudary_; 00412 Array<bool> edgeIsMeshDomainBonudary_; 00413 00414 //------ Element (processor) ownership ----- 00415 /** contains the ownership of the local elements [dim][ID] */ 00416 Array< Array< short signed int > > elementOwner_; 00417 00418 //---- hierarchical storage ----- 00419 00420 /** [cellID] , the child index in the parent */ 00421 Array<short int> indexInParent_; 00422 /** [cellID] , the LID of the parent cell */ 00423 Array<int> parentCellLID_; 00424 /** [cellID] , actual level of the cell*/ 00425 Array<short int> cellLevel_; 00426 /** [cellID] , if the element is leaf */ 00427 Array<bool> isCellLeaf_; 00428 /** [cellID] , if the cell is complete outside the user defined mesh domain*/ 00429 Array<bool> isCellOut_; 00430 /** [cellID] , children of the cell*/ 00431 Array< Array<int> > cellsChildren_; 00432 // ---- "hanging" info storage --- 00433 /** [pointID] , true if the node is hanging , false otherwise */ 00434 Array<bool> isPointHanging_; 00435 /** [edgeID] , true if the edge is hanging , false otherwise*/ 00436 Array<bool> isEdgeHanging_; 00437 00438 // ---- hanging element and refinement (temporary) storage --- 00439 00440 /** [vertexID] - > { h P1 LID , h P2 LID , h E1 LID , h E2 LID , h E3 LID } */ 00441 Array< Hashtable< int, Array<int> > > hangElmStore_; 00442 /** Neighbor Cell can mark the cell to provoke refinement */ 00443 Array<short int> refineCell_; 00444 00445 // ---- leaf mapping , GID and LID --- (points do not need this, all points are also leaf points) 00446 00447 /** [leaf_vertexLID] , the value must be a positive number */ 00448 Array<int> vertexLeafToLIDMapping_; 00449 /** [leaf_edgeLID] , the value must be a positive number */ 00450 Array<int> edgeLeafToLIDMapping_; 00451 /** [leaf_cellLID] , the value must be a positive number */ 00452 Array<int> cellLeafToLIDMapping_; 00453 /** [vertexLID] if vertex is inside the domain then > 0 , -1 otherwise */ 00454 Array<int> vertexLIDToLeafMapping_; 00455 /** [edgeLID] if edge is leaf(or inside the domain) then > 0 , -1 otherwise */ 00456 Array<int> edgeLIDToLeafMapping_; 00457 /** [cellLID] if cell is leaf(or inside the domain) then > 0 , -1 otherwise */ 00458 Array<int> cellLIDToLeafMapping_; 00459 00460 /** leaf LID numbering*/ 00461 int nrVertexLeafLID_; 00462 int nrCellLeafLID_; 00463 int nrEdgeLeafLID_; 00464 00465 /** [leaf_vertexGID] , the value must be a positive number */ 00466 Array<int> vertexLeafToGIDMapping_; 00467 /** [leaf_edgeGID] , the value must be a positive number */ 00468 Array<int> edgeLeafToGIDMapping_; 00469 /** [leaf_cellGID] , the value must be a positive number */ 00470 Array<int> cellLeafToGIDMapping_; 00471 /** [vertexGID] if vertex is inside the domain then > 0 , -1 otherwise */ 00472 Array<int> vertexGIDToLeafMapping_; 00473 /** [edgeGID] if edge is leaf(or inside the domain) then > 0 , -1 otherwise */ 00474 Array<int> edgeGIDToLeafMapping_; 00475 /** [cellGID] if cell is leaf(or inside the domain) then > 0 , -1 otherwise */ 00476 Array<int> cellGIDToLeafMapping_; 00477 00478 /** leaf GID numbering*/ 00479 int nrVertexLeafGID_; 00480 int nrCellLeafGID_; 00481 int nrEdgeLeafGID_; 00482 00483 // ------------- static data ------------------- 00484 00485 /** the offset in the X coordinate on the reference cell*/ 00486 static int offs_Points_x_[4]; 00487 00488 /** the offset in the Y coordinate on the reference cell*/ 00489 static int offs_Points_y_[4]; 00490 00491 /** stores the facet information on the reference Cell*/ 00492 static int edge_Points_localIndex[4][2]; 00493 00494 }; 00495 } 00496 00497 00498 #endif /* SUNDANCEHNMESH2D_H_ */