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 * SundanceHNMesh3D.h 00043 * 00044 * Created on: May 30, 2010 00045 * Author: benk 00046 */ 00047 00048 #ifndef SUNDANCEHNMESH3D_H_ 00049 #define SUNDANCEHNMESH3D_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 3D 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 HNMesh3D : public MeshBase{ 00081 00082 public: 00083 /** 00084 * The Ctor for the dummy grid with hanging nodes */ 00085 HNMesh3D(int dim, 00086 const MPIComm& comm, 00087 const MeshEntityOrder& order); 00088 00089 /** 00090 * The Ctor for the HNMesh3D grid in 3D*/ 00091 void createMesh( 00092 double position_x, 00093 double position_y, 00094 double position_z, 00095 double offset_x, 00096 double offset_y, 00097 double offset_z, 00098 int resolution_x, 00099 int resolution_y, 00100 int resolution_z, 00101 const RefinementClass& refineClass, 00102 const MeshDomainDef& meshDomain ); 00103 00104 /** Dtor */ 00105 virtual ~HNMesh3D(){;} 00106 00107 /** 00108 * Get the number of cells having dimension dim 00109 */ 00110 virtual int numCells(int dim) const; 00111 00112 /** 00113 * Return the position of the i-th node 00114 */ 00115 virtual Point nodePosition(int i) const; 00116 00117 /** 00118 * Return a view of the i-th node's position 00119 */ 00120 const double* nodePositionView(int i) const; 00121 00122 00123 00124 /** 00125 * Compute the jacobians of a batch of cells, returning the 00126 * result via reference argument 00127 * 00128 * @param cellDim dimension of the cells whose Jacobians are to 00129 * be computed 00130 * @param cellLID local indices of the cells for which Jacobians 00131 * are to be computed 00132 * @param jBatch reference to the resulting Jacobian batch 00133 */ 00134 virtual void getJacobians(int cellDim, const Array<int>& cellLID, 00135 CellJacobianBatch& jBatch) const; 00136 00137 /** 00138 * Compute the diameters of a batch of cells, 00139 * result via reference argument 00140 * 00141 * @param cellDim dimension of the cells whose diameters are to 00142 * be computed 00143 * @param cellLID local indices of the cells for which diameters 00144 * are to be computed 00145 * @param diameters reference to the array of cell diameters 00146 */ 00147 virtual void getCellDiameters(int cellDim, const Array<int>& cellLID, 00148 Array<double>& cellDiameters) const; 00149 00150 00151 /** 00152 * Map reference quadrature points to physical points on the 00153 * given cells. 00154 */ 00155 virtual void pushForward(int cellDim, const Array<int>& cellLID, 00156 const Array<Point>& refQuadPts, 00157 Array<Point>& physQuadPts) const; 00158 00159 /** 00160 * Return the rank of the processor that owns the given cell 00161 */ 00162 virtual int ownerProcID(int cellDim, int cellLID) const; 00163 00164 /** 00165 * Return the number of facets of the given cell 00166 */ 00167 virtual int numFacets(int cellDim, int cellLID, 00168 int facetDim) const; 00169 00170 /** 00171 * Return the local ID of a facet cell 00172 * @param cellDim dimension of the cell whose facets are being obtained 00173 * @param cellLID local index of the cell whose 00174 * facets are being obtained 00175 * @param facetDim dimension of the desired facet 00176 * @param facetIndex index into the list of the cell's facets 00177 * @param facetOrientation orientation of the facet as seen from 00178 * the given cell (returned via reference) 00179 */ 00180 virtual int facetLID(int cellDim, int cellLID, 00181 int facetDim, int facetIndex, 00182 int& facetOrientation) const; 00183 00184 /** 00185 * Return by reference argument an array containing&(elemVerts_.value(cellLID, 0)) 00186 * the LIDs of the facetDim-dimensional facets of the 00187 * given batch of cells 00188 */ 00189 virtual void getFacetLIDs(int cellDim, 00190 const Array<int>& cellLID, 00191 int facetDim, 00192 Array<int>& facetLID, 00193 Array<int>& facetSign) const; 00194 00195 /** 00196 * Return a view of an element's zero-dimensional facets, 00197 * @return an array of integers with the indexes of the points which for it 00198 */ 00199 const int* elemZeroFacetView(int cellLID) const; 00200 00201 /** 00202 * Return the number of maximal cofacets of the given cell 00203 */ 00204 virtual int numMaxCofacets(int cellDim, int cellLID) const; 00205 00206 /** 00207 * Return the local ID of a maximal cofacet cell 00208 * @param cellDim dimension of the cell whose cofacets are being obtained 00209 * @param cellLID local index of the cell whose 00210 * cofacets are being obtained 00211 * @param cofacetIndex which maximal cofacet to get 00212 * @param facetIndex index of the cell cellLID into the list of the 00213 * maximal cell's facets 00214 */ 00215 virtual int maxCofacetLID(int cellDim, int cellLID, 00216 int cofacetIndex, 00217 int& facetIndex) const; 00218 /** 00219 * Get the LIDs of the maximal cofacets for a batch of codimension-one 00220 * cells. 00221 * 00222 * \param cellLIDs [in] array of LIDs of the cells whose cofacets are 00223 * being obtained 00224 * \param cofacets [out] the batch of cofacets 00225 */ 00226 virtual void getMaxCofacetLIDs(const Array<int>& cellLIDs, 00227 MaximalCofacetBatch& cofacets) const; 00228 00229 00230 /** 00231 * Find the cofacets of the given cell 00232 * @param cellDim dimension of the cell whose cofacets are being obtained 00233 * @param cellLID local index of the cell whose 00234 * cofacets are being obtained 00235 * @param cofacetDim dimension of the cofacets to get 00236 * @param cofacetLIDs LIDs of the cofacet 00237 */ 00238 void getCofacets(int cellDim, int cellLID, 00239 int cofacetDim, Array<int>& cofacetLIDs) const; 00240 00241 /** 00242 * Find the local ID of a cell given its global index 00243 */ 00244 virtual int mapGIDToLID(int cellDim, int globalIndex) const; 00245 00246 /** 00247 * Determine whether a given cell GID exists on this processor 00248 */ 00249 virtual bool hasGID(int cellDim, int globalIndex) const; 00250 00251 00252 /** 00253 * Find the global ID of a cell given its local index 00254 */ 00255 virtual int mapLIDToGID(int cellDim, int localIndex) const; 00256 00257 /** 00258 * Get the type of the given cell 00259 */ 00260 virtual CellType cellType(int cellDim) const; 00261 00262 00263 /** Get the label of the given cell */ 00264 virtual int label(int cellDim, int cellLID) const; 00265 00266 /** Get the labels for a batch of cells */ 00267 virtual void getLabels(int cellDim, const Array<int>& cellLID, 00268 Array<int>& labels) const; 00269 00270 00271 00272 /** Get the list of all labels defined for cells of the given dimension */ 00273 virtual Set<int> getAllLabelsForDimension(int cellDim) const; 00274 00275 /** 00276 * Get the cells associated with a specified label. The array 00277 * cellLID will be filled with those cells of dimension cellDim 00278 * having the given label. 00279 */ 00280 virtual void getLIDsForLabel(int cellDim, int label, Array<int>& cellLIDs) const; 00281 00282 00283 /** Set the label of the given cell */ 00284 virtual void setLabel(int cellDim, int cellLID, int label); 00285 00286 /** Coordinate intermediate cell definitions across processors */ 00287 virtual void assignIntermediateCellGIDs(int cellDim); 00288 00289 /** */ 00290 virtual bool hasIntermediateGIDs(int dim) const; 00291 00292 00293 /** Function returns true if the mesh allows hanging nodes (by refinement), 00294 * false otherwise */ 00295 virtual bool allowsHangingHodes() const { return true; } 00296 00297 /** Function returns true if the specified element is a "hanging" element 00298 * false otherwise <br> 00299 * @param dim [in] should be between 0 , D-1 00300 * @param cellLID [in] the local ID of the element */ 00301 virtual bool isElementHangingNode(int cellDim , int cellLID) const ; 00302 00303 /** Returns the index in the parent maxdim Cell of the refinement tree 00304 * @param maxCellLID [in] the LID of the cell */ 00305 virtual int indexInParent(int maxCellLID) const ; 00306 00307 /** How many children has a refined element. <br> 00308 * This function provides information of either we have bi or trisection */ 00309 virtual int maxChildren() const {return 27;} 00310 00311 /** Function returns the facets of the maxdim parent cell (needed for HN treatment) <br> 00312 * @param childCellLID [in] the LID of the maxdim cell, whos parents facets we want 00313 * @param dimFacets [in] the dimension of the facets which we want to have 00314 * @param facetsLIDs [out] the LID of the parents facets (all) in the defined order 00315 * @param parentCellLIDs [out] the maxdim parent cell LID */ 00316 virtual void returnParentFacets( int childCellLID , int dimFacets , 00317 Array<int> &facetsLIDs , int &parentCellLIDs ) const; 00318 00319 private: 00320 00321 /** For HN , returns parent facets, if the facet is not leaf, then return -1 at that place */ 00322 int facetLID_tree(int cellDim, int cellLID, 00323 int facetDim, int facetIndex) const; 00324 00325 /** adds one vertex to the mesh */ 00326 void addVertex(int vertexLID , int ownerProc , bool isHanging , 00327 double coordx , double coordy , double coordz , const Array<int> &maxCoF); 00328 00329 /** adds one edge to the mesh */ 00330 void addEdge(int edgeLID , int ownerProc , bool isHanging , int edgeOrientation , 00331 const Array<int> &vertexLIDs , const Array<int> &maxCoF); 00332 00333 /** adds one edge to the mesh */ 00334 void addFace(int faceLID , int ownerProc , bool isHanging , int faceOrientation , 00335 const Array<int> &vertexLIDs , const Array<int> &edgeLIDs , 00336 const Array<int> &maxCoF); 00337 00338 /** adds one cell(3D) to the mesh <br> 00339 * cell must be always leaf*/ 00340 void addCell(int cellLID , int ownerProc , 00341 int indexInParent , int parentCellLID , int level , 00342 const Array<int> &faceLIDs , const Array<int> &edgeLIDs , 00343 const Array<int> &vertexLIDs) ; 00344 00345 /** creates the mesh on the coarsest level as it is specified */ 00346 void createCoarseMesh(); 00347 00348 /** Does one single refinement iteration. <br> 00349 * Iterates trough all cells which are owned by the processor and refines if necessary */ 00350 bool oneRefinementIteration(); 00351 00352 /** refine the given cell by cellID, (method assumes that this cell can be refined)*/ 00353 void refineCell(int cellLID); 00354 00355 /** Create Leaf numbering */ 00356 void createLeafNumbering(); 00357 00358 /** Create Leaf numbering, but with a better load balancing for parallel case */ 00359 void createLeafNumbering_sophisticated(); 00360 00361 /** estimate the load of one cell*/ 00362 int estimateCellLoad(int ID); 00363 00364 /** marks the cells recursivly in the tree (and facets) owned by one processor*/ 00365 void markCellsAndFacets(int cellID , int procID); 00366 00367 /** this updates the array with the local index of vertex, edge and face in the static array <br> 00368 * this method is only used in the coarse mesh creation */ 00369 void updateLocalCoarseNumbering(int ix , int iy , int iz , int Nx , int Ny); 00370 00371 /** Used for refinemement 00372 * @param cellDim the requested elements dimension (if exists) 00373 * @param usedge , whould we look for that element in edge or in face 00374 * @param parentID , ID of the parent edge or face 00375 * @param parentOffset , the offset in the parent */ 00376 int getHangingElement(int cellDim, bool useEdge, int parentID , int parentOffset); 00377 00378 /** Used for refinemement 00379 * @param cellDim the requested elements dimension (if exists) 00380 * @param cellID , ID of the new hanging cell 00381 * @param usedge , whould we look for that element in edge or in face 00382 * @param parentID , ID of the parent edge or face 00383 * @param parentOffset , the offset in the parent */ 00384 void addHangingElement(int cellDim, int cellID ,bool useEdge, int parentID , int parentOffset); 00385 00386 /** 00387 * @param cellDim 00388 * @param cellID */ 00389 int numMaxCofacets_ID(int cellDim, int cellID); 00390 00391 /** Number of processors */ 00392 int nrProc_; 00393 int myRank_; 00394 /** The communicator */ 00395 const MPIComm& _comm; 00396 /** */ 00397 double _pos_x; 00398 /** */ 00399 double _pos_y; 00400 /** */ 00401 double _pos_z; 00402 /** */ 00403 double _ofs_x; 00404 /** */ 00405 double _ofs_y; 00406 /** */ 00407 double _ofs_z; 00408 /** */ 00409 int _res_x; 00410 /** */ 00411 int _res_y; 00412 /** */ 00413 int _res_z; 00414 /** */ 00415 mutable RefinementClass refineClass_; 00416 /** */ 00417 mutable MeshDomainDef meshDomain_; 00418 00419 //------ Point storage ---- 00420 /** all the global points index is [ID] */ 00421 Array<Point> points_; 00422 /** [4] the nr of ID per dim*/ 00423 Array<int> nrElem_; 00424 /** [4] the nr of owned elements per dim*/ 00425 Array<int> nrElemOwned_; 00426 00427 //----- Facets ----- -1 -> MeshBoundary 00428 00429 /** [cellID][8]*/ 00430 Array< Array<int> > cellsPoints_; 00431 /** [cellID][12]*/ 00432 Array< Array<int> > cellsEdges_; 00433 /** [cellID][6]*/ 00434 Array< Array<int> > cellsFaces_; 00435 /** [faceID][4]*/ 00436 Array< Array<int> > faceEdges_; 00437 /** [faceID][4]*/ 00438 Array< Array<int> > facePoints_; 00439 /** [edgeID][2]*/ 00440 Array< Array<int> > edgePoints_; 00441 00442 /** stores the edge orientation {0,1,2} 00443 * [edgeID]*/ 00444 Array<short int> edgeOrientation_; 00445 /** stores the face orientation {0,1,2} 00446 * [faceID]*/ 00447 Array<short int> faceOrientation_; 00448 00449 // ----- MaxCofacets ----; -1 -> MeshBoundary , -2 -> ProcBoundary 00450 00451 /** [faceID][2] */ 00452 Array< Array<int> > faceMaxCoF_; 00453 /** [edgeID][4] */ 00454 Array< Array<int> > edgeMaxCoF_; 00455 /** [pointID][8]*/ 00456 Array< Array<int> > pointMaxCoF_; 00457 00458 //------ Element (processor) ownership ----- 00459 00460 /** contains the ownership of the local elements [dim][ID] */ 00461 Array< Array<short int> > elementOwner_; 00462 00463 //---- hierarchical storage ----- 00464 00465 /** [cellID] , the child index in the parent */ 00466 Array<short int> indexInParent_; 00467 /** [cellID] , the LID of the parent cell */ 00468 Array<int> parentCellLID_; 00469 /** [cellID] , actual level of the cell*/ 00470 Array<short int> cellLevel_; 00471 /** [cellID] , if the element is leaf */ 00472 Array<bool> isCellLeaf_; 00473 /** [cellID] , if the cell is complete outside the user defined mesh domain*/ 00474 Array<bool> isCellOut_; 00475 /** [cellID] , children of the cell*/ 00476 Array< Array<int> > cellsChildren_; 00477 00478 // ---- "hanging" info storage --- 00479 /** [pointID] , true if the node is hanging , false otherwise */ 00480 Array<bool> isPointHanging_; 00481 /** [edgeID] , true if the edge is hanging , false otherwise*/ 00482 Array<bool> isEdgeHanging_; 00483 /** [faceID] , true if the face is hanging , false otherwise*/ 00484 Array<bool> isFaceHanging_; 00485 00486 // ---- hanging element and refinement (temporary) storage --- 00487 00488 /** [edgeID] - > { h P0 LID , h P1 LID , h E0 LID , h E1 LID , h E2 LID } <br> 00489 * These elements will only be put on not hanging when we access it 3 times <br> 00490 * together 5 elements*/ 00491 Hashtable< int, Array<int> > edgeHangingElmStore_; 00492 00493 /** the counter for each edge which stores hanging node information, when the counter reaches 2 (3 read access)<br> 00494 * then the hanging elements can be marked as not hanging <br> 00495 * At the beginning this should be equal with the numMaxCofacets */ 00496 Array<short int> hangingAccessCount_; 00497 00498 /** [faceID] - > { h P0 LID , h P1 LID , h P2 LID , h P3 LID , <br> 00499 * h E0 LID , h E1 LID , h E2 LID , h E3 LID , h E4 LID , h E5 LID , (4+)<br> 00500 * h E6 LID , h E7 LID , h E8 LID , h E9 LID , h E10 LID, h E11 LID , (16+) <br> 00501 * h F0 LID , h F1 LID , h F2 LID , h F3 LID , h F4 LID , h F5 LID , <br> 00502 * h F6 LID , h F7 LID , h F8 LID } <br> 00503 * together 16+9 = 25 elements*/ 00504 Hashtable< int, Array<int> > faceHangingElmStore_; 00505 00506 /** Neighbor Cell can mark the cell to provoke refinement */ 00507 Array<short int> refineCell_; 00508 00509 // ---- leaf mapping , GID and LID --- (points need this because) 00510 00511 /** [leaf_vertexLID] , the value must be a positive number */ 00512 Array<int> vertexLeafToLIDMapping_; 00513 /** [leaf_edgeLID] , the value must be a positive number */ 00514 Array<int> edgeLeafToLIDMapping_; 00515 /** [leaf_faceLID] , the value must be a positive number */ 00516 Array<int> faceLeafToLIDMapping_; 00517 /** [leaf_cellLID] , the value must be a positive number */ 00518 Array<int> cellLeafToLIDMapping_; 00519 /** [vertexLID] if vertex is inside the domain then > 0 , -1 otherwise */ 00520 Array<int> vertexLIDToLeafMapping_; 00521 /** [edgeLID] if edge is leaf(or inside the domain) then > 0 , -1 otherwise */ 00522 Array<int> edgeLIDToLeafMapping_; 00523 /** [faceLID] if face is leaf(or inside the domain) then > 0 , -1 otherwise */ 00524 Array<int> faceLIDToLeafMapping_; 00525 /** [cellLID] if cell is leaf(or inside the domain) then > 0 , -1 otherwise */ 00526 Array<int> cellLIDToLeafMapping_; 00527 00528 /** leaf LID numbering */ 00529 int nrVertexLeafLID_; 00530 int nrEdgeLeafLID_; 00531 int nrFaceLeafLID_; 00532 int nrCellLeafLID_; 00533 00534 /** [leaf_vertexGID] , the value must be a positive number */ 00535 Array<int> vertexLeafToGIDMapping_; 00536 /** [leaf_edgeGID] , the value must be a positive number */ 00537 Array<int> edgeLeafToGIDMapping_; 00538 /** [leaf_faceGID] , the value must be a positive number */ 00539 Array<int> faceLeafToGIDMapping_; 00540 /** [leaf_cellGID] , the value must be a positive number */ 00541 Array<int> cellLeafToGIDMapping_; 00542 /** [vertexGID] if vertex is inside the domain then > 0 , -1 otherwise */ 00543 Array<int> vertexGIDToLeafMapping_; 00544 /** [edgeGID] if edge is leaf(or inside the domain) then > 0 , -1 otherwise */ 00545 Array<int> edgeGIDToLeafMapping_; 00546 /** [faceGID] if edge is leaf(or inside the domain) then > 0 , -1 otherwise */ 00547 Array<int> faceGIDToLeafMapping_; 00548 /** [cellGID] if cell is leaf(or inside the domain) then > 0 , -1 otherwise */ 00549 Array<int> cellGIDToLeafMapping_; 00550 00551 /** leaf GID numbering */ 00552 int nrVertexLeafGID_; 00553 int nrEdgeLeafGID_; 00554 int nrFaceLeafGID_; 00555 int nrCellLeafGID_; 00556 00557 // ------------- static data ------------------- 00558 00559 /** the offset in the X coordinate on the reference cell*/ 00560 static int offs_Points_x_[8]; 00561 00562 /** the offset in the Y coordinate on the reference cell*/ 00563 static int offs_Points_y_[8]; 00564 00565 /** the offset in the Z coordinate on the reference cell*/ 00566 static int offs_Points_z_[8]; 00567 00568 /** stores the facet information on the reference Cell*/ 00569 static int edge_Points_localIndex[12][2]; 00570 00571 /** the edge orientation (local orientation)*/ 00572 static int edge_Orientation[12]; 00573 static int edge_MaxCofacetIndex[3][4]; 00574 static int edge_MaxCof[12]; 00575 00576 /** face edge-facet information */ 00577 static int face_Edges_localIndex[6][4]; 00578 00579 /** face point-facet information */ 00580 static int face_Points_localIndex[6][4]; 00581 00582 /** face orientation (local orientation)*/ 00583 static int face_Orientation[6]; 00584 static int face_MaxCofacetIndex[3][2]; 00585 static int face_MaxCof[6]; 00586 00587 /** used for coarse mesh creation for the vertex indexing */ 00588 static int vInd[8]; 00589 00590 /** used for coarse mesh creation for the edge indexing */ 00591 static int eInd[12]; 00592 00593 /** used for coarse mesh creation for the face indexing */ 00594 static int fInd[6]; 00595 00596 // the X and the Y coordinates of the newly 00597 static double vertex_X[64]; 00598 00599 static double vertex_Y[64]; 00600 00601 static double vertex_Z[64]; 00602 00603 // face index is above 20 00604 static int vertexToParentEdge[64]; 00605 // 00606 static int vertexInParentIndex[64]; 00607 // 00608 static int edgeToParentEdge[144]; 00609 // 00610 static int edgeInParentIndex[144]; 00611 // 00612 static int faceToParentFace[108]; 00613 // 00614 static int faceInParentIndex[108]; 00615 00616 }; 00617 } 00618 00619 00620 #endif /* SUNDANCEHNMESH3D_H_ */