// Manually add a surface element of a given type to an existing mesh object DLL_HEADER void Ng_AddSurfaceElement (Ng_Mesh * mesh, Ng_Surface_Element_Type et, int * pi) { Mesh * m = (Mesh*)mesh; Element2d el (3); el.SetIndex (1); el.PNum(1) = pi[0]; el.PNum(2) = pi[1]; el.PNum(3) = pi[2]; m->AddSurfaceElement (el); }
/* Philippose Rajan - 11 June 2009 Function to calculate the surface normal at a given vertex of a surface element, with respect to that surface element. This function is used by the boundary layer generation function, in order to calculate the effective direction in which the prismatic layer should grow */ void GetSurfaceNormal(Mesh & mesh, Element2d & el, int Vertex, Vec3d & SurfaceNormal) { int Vertex_A; int Vertex_B; Vertex_A = Vertex + 1; if(Vertex_A > el.GetNP()) Vertex_A = 1; Vertex_B = Vertex - 1; if(Vertex_B <= 0) Vertex_B = el.GetNP(); Vec3d Vect_A,Vect_B; Vect_A = mesh.Point(el.PNum(Vertex_A)) - mesh.Point(el.PNum(Vertex)); Vect_B = mesh.Point(el.PNum(Vertex_B)) - mesh.Point(el.PNum(Vertex)); SurfaceNormal = Cross(Vect_A,Vect_B); SurfaceNormal.Normalize(); }
void WriteEdgeElementFormat (const Mesh & mesh, const string & filename) { cout << "write edge element format" << endl; const MeshTopology * top = &mesh.GetTopology(); int npoints = mesh.GetNP(); int nelements = mesh.GetNE(); int nsurfelem = mesh.GetNSE(); int nedges = top->GetNEdges(); int i, j; int inverttets = mparam.inverttets; int invertsurf = mparam.inverttrigs; Array<int> edges; ofstream outfile (filename.c_str()); outfile.precision(6); outfile.setf (ios::fixed, ios::floatfield); outfile.setf (ios::showpoint); // vertices with coordinates outfile << npoints << "\n"; for (i = 1; i <= npoints; i++) { const Point3d & p = mesh.Point(i); outfile.width(10); outfile << p.X() << " "; outfile.width(9); outfile << p.Y() << " "; outfile.width(9); outfile << p.Z() << "\n"; } // element - edge - list outfile << nelements << " " << nedges << "\n"; for (i = 1; i <= nelements; i++) { Element el = mesh.VolumeElement(i); if (inverttets) el.Invert(); outfile.width(4); outfile << el.GetIndex() << " "; outfile.width(8); outfile << el.GetNP(); for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile.width(8); outfile << el.PNum(j); } top->GetElementEdges(i,edges); outfile << endl << " "; outfile.width(8); outfile << edges.Size(); for (j=1; j <= edges.Size(); j++) { outfile << " "; outfile.width(8); outfile << edges[j-1]; } outfile << "\n"; // orientation: top->GetElementEdgeOrientations(i,edges); outfile << " "; for (j=1; j <= edges.Size(); j++) { outfile << " "; outfile.width(8); outfile << edges[j-1]; } outfile << "\n"; } // surface element - edge - list (with boundary conditions) outfile << nsurfelem << "\n"; for (i = 1; i <= nsurfelem; i++) { Element2d el = mesh.SurfaceElement(i); if (invertsurf) el.Invert(); outfile.width(4); outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; outfile.width(8); outfile << el.GetNP(); for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile.width(8); outfile << el.PNum(j); } top->GetSurfaceElementEdges(i,edges); outfile << endl << " "; outfile.width(8); outfile << edges.Size(); for (j=1; j <= edges.Size(); j++) { outfile << " "; outfile.width(8); outfile << edges[j-1]; } outfile << "\n"; } int v1, v2; // edge - vertex - list outfile << nedges << "\n"; for (i=1; i <= nedges; i++) { top->GetEdgeVertices(i,v1,v2); outfile.width(4); outfile << v1; outfile << " "; outfile.width(8); outfile << v2 << endl; } }
void WriteNeutralFormat (const Mesh & mesh, const string & filename) { cout << "write neutral, new" << endl; int np = mesh.GetNP(); int ne = mesh.GetNE(); int nse = mesh.GetNSE(); int nseg = mesh.GetNSeg(); int i, j; int inverttets = mparam.inverttets; int invertsurf = mparam.inverttrigs; ofstream outfile (filename.c_str()); outfile.precision(6); outfile.setf (ios::fixed, ios::floatfield); outfile.setf (ios::showpoint); outfile << np << "\n"; for (i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile.width(10); outfile << p.X() << " "; outfile.width(9); outfile << p.Y() << " "; if (mesh.GetDimension() == 3) { outfile.width(9); outfile << p.Z(); } outfile << "\n"; } if (mesh.GetDimension() == 3) { outfile << ne << "\n"; for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (inverttets) el.Invert(); outfile.width(4); outfile << el.GetIndex() << " "; for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile.width(8); outfile << el.PNum(j); } outfile << "\n"; } } outfile << nse << "\n"; for (i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (invertsurf) el.Invert(); outfile.width(4); outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile.width(8); outfile << el.PNum(j); } outfile << "\n"; } if (mesh.GetDimension() == 2) { outfile << nseg << "\n"; for (i = 1; i <= nseg; i++) { const Segment & seg = mesh.LineSegment(i); outfile.width(4); outfile << seg.si << " "; outfile << " "; outfile.width(8); outfile << seg[0]; outfile << " "; outfile.width(8); outfile << seg[1]; outfile << "\n"; } } }
/*! GMSH v2.xx mesh format export function * * This function extends the export capabilities of * Netgen to include the GMSH v2.xx File Format. * * Current features of this function include: * * 1. Exports Triangles, Quadrangles and Tetrahedra \n * 2. Supports upto second order elements of each type * */ void WriteGmsh2Format (const Mesh & mesh, const CSGeometry & geom, const string & filename) { ofstream outfile (filename.c_str()); outfile.precision(6); outfile.setf (ios::fixed, ios::floatfield); outfile.setf (ios::showpoint); int np = mesh.GetNP(); /// number of points in mesh int ne = mesh.GetNE(); /// number of 3D elements in mesh int nse = mesh.GetNSE(); /// number of surface elements (BC) int i, j, k, l; /* * 3D section : Volume elements (currently only tetrahedra) */ if ((ne > 0) && (mesh.VolumeElement(1).GetNP() <= 10) && (mesh.SurfaceElement(1).GetNP() <= 6)) { cout << "Write GMSH v2.xx Format \n"; cout << "The GMSH v2.xx export is currently available for elements upto 2nd Order\n" << endl; int inverttets = mparam.inverttets; int invertsurf = mparam.inverttrigs; /// Prepare GMSH 2.0 file (See GMSH 2.0 Documentation) outfile << "$MeshFormat\n"; outfile << (float)2.0 << " " << (int)0 << " " << (int)sizeof(double) << "\n"; outfile << "$EndMeshFormat\n"; /// Write nodes outfile << "$Nodes\n"; outfile << np << "\n"; for (i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << i << " "; /// node number outfile << p.X() << " "; outfile << p.Y() << " "; outfile << p.Z() << "\n"; } outfile << "$EndNodes\n"; /// write elements (both, surface elements and volume elements) outfile << "$Elements\n"; outfile << ne + nse << "\n"; //// number of elements + number of surfaces BC for (i = 1; i <= nse; i++) { int elType = 0; Element2d el = mesh.SurfaceElement(i); if(invertsurf) el.Invert(); if(el.GetNP() == 3) elType = GMSH_TRIG; //// GMSH Type for a 3 node triangle if(el.GetNP() == 6) elType = GMSH_TRIG6; //// GMSH Type for a 6 node triangle if(elType == 0) { cout << " Invalid surface element type for Gmsh 2.0 3D-Mesh Export Format !\n"; return; } outfile << i; outfile << " "; outfile << elType; outfile << " "; outfile << "2"; //// Number of tags (2 => Physical and elementary entities) outfile << " "; outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; /// that means that physical entity = elementary entity (arbitrary approach) outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile << el.PNum(triGmsh[j]); } outfile << "\n"; } for (i = 1; i <= ne; i++) { int elType = 0; Element el = mesh.VolumeElement(i); if (inverttets) el.Invert(); if(el.GetNP() == 4) elType = GMSH_TET; //// GMSH Element type for 4 node tetrahedron if(el.GetNP() == 10) elType = GMSH_TET10; //// GMSH Element type for 10 node tetrahedron if(elType == 0) { cout << " Invalid volume element type for Gmsh 2.0 3D-Mesh Export Format !\n"; return; } outfile << nse + i; //// element number (Remember to add on surface elements) outfile << " "; outfile << elType; outfile << " "; outfile << "2"; //// Number of tags (2 => Physical and elementary entities) outfile << " "; outfile << 100000 + el.GetIndex(); /// that means that physical entity = elementary entity (arbitrary approach) outfile << " "; outfile << 100000 + el.GetIndex(); /// volume number outfile << " "; for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile << el.PNum(tetGmsh[j]); } outfile << "\n"; } outfile << "$EndElements\n"; } /* * End of 3D section */ /* * 2D section : available for triangles and quadrangles * upto 2nd Order */ else if(ne == 0) /// means that there's no 3D element { cout << "\n Write Gmsh v2.xx Surface Mesh (triangle and/or quadrangles upto 2nd Order)" << endl; /// Prepare GMSH 2.0 file (See GMSH 2.0 Documentation) outfile << "$MeshFormat\n"; outfile << (float)2.0 << " " << (int)0 << " " << (int)sizeof(double) << "\n"; outfile << "$EndMeshFormat\n"; /// Write nodes outfile << "$Nodes\n"; outfile << np << "\n"; for (i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << i << " "; /// node number outfile << p.X() << " "; outfile << p.Y() << " "; outfile << p.Z() << "\n"; } outfile << "$EndNodes\n"; /// write triangles & quadrangles outfile << "$Elements\n"; outfile << nse << "\n"; for (k = 1; k <= nse; k++) { int elType = 0; const Element2d & el = mesh.SurfaceElement(k); if(el.GetNP() == 3) elType = GMSH_TRIG; //// GMSH Type for a 3 node triangle if(el.GetNP() == 6) elType = GMSH_TRIG6; //// GMSH Type for a 6 node triangle if(el.GetNP() == 4) elType = GMSH_QUAD; //// GMSH Type for a 4 node quadrangle if(el.GetNP() == 8) elType = GMSH_QUAD8; //// GMSH Type for an 8 node quadrangle if(elType == 0) { cout << " Invalid surface element type for Gmsh 2.0 2D-Mesh Export Format !\n"; return; } outfile << k; outfile << " "; outfile << elType; outfile << " "; outfile << "2"; outfile << " "; outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; /// that means that physical entity = elementary entity (arbitrary approach) outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; for (l = 1; l <= el.GetNP(); l++) { outfile << " "; if((elType == GMSH_TRIG) || (elType == GMSH_TRIG6)) { outfile << el.PNum(triGmsh[l]); } else if((elType == GMSH_QUAD) || (elType == GMSH_QUAD8)) { outfile << el.PNum(quadGmsh[l]); } } outfile << "\n"; } outfile << "$EndElements\n"; } /* * End of 2D section */ else { cout << " Invalid element type for Gmsh v2.xx Export Format !\n"; } } // End: WriteGmsh2Format
void WriteJCMFormat (const Mesh & mesh, const NetgenGeometry & geom, const string & filename) { if (mesh.GetDimension() != 3) { cout <<"\n Error: Dimension 3 only supported by this output format!"<<endl; return; } int bc_at_infinity = 0; int i, j, jj, ct(0), counter; double dx1, dx2, dx3, dy1, dy2, dy3, dz1, dz2, dz3, vol; // number of points int np = mesh.GetNP(); // Identic points Array<int,1> identmap1, identmap2, identmap3; mesh.GetIdentifications().GetMap(1, identmap1); mesh.GetIdentifications().GetMap(2, identmap2); mesh.GetIdentifications().GetMap(3, identmap3); // number of volume elements int ne = mesh.GetNE(); int ntets = 0; int nprisms = 0; for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (el.GetNP() == 4) { ntets++; // Check that no two points on a tetrahedron are identified with each other for (j = 1; j <= 4; j++) for (jj = 1; jj <=4; jj++) { if (identmap1.Elem(el.PNum(j)) == el.PNum(jj)) { cout << "\n Error: two points on a tetrahedron identified (1) with each other" << "\n REFINE MESH !" << endl; return; } if (identmap2.Elem(el.PNum(j)) == el.PNum(jj)) { cout << "\n Error: two points on a tetrahedron identified (2) with each other" << "\n REFINE MESH !" << endl; return; } if (identmap3.Elem(el.PNum(j)) == el.PNum(jj)) { cout << "\n Error: two points on a tetrahedron identified (3) with each other" << "\n REFINE MESH !" << endl; return; } } } else if (el.GetNP() == 6) nprisms++; } if ( ne != (ntets+nprisms)) { cout<< "\n Error in determining number of volume elements!\n" << "\n Prisms and tetrahedra only implemented in the JCMwave format!\n"<<endl; return; } if (nprisms > 0) cout << " Please note: Boundaries at infinity have to carry the bc-attribute '-bc=" << bc_at_infinity <<"'."<<endl; // number of surface elements int nse = mesh.GetNSE(); // number of boundary triangles int nbtri = 0; // number of boundary quadrilaterals int nbquad = 0; // array with 1 if point on any tetra, 0 else // this is needed in order to arrange the prism points in the right order Array<int,1> pointsOnTetras; pointsOnTetras.SetSize (mesh.GetNP()); pointsOnTetras = 0; for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (el.GetNP() == 4) { for (j = 1; j <= 4; j++) pointsOnTetras.Set(int (el.PNum(j)),1); } } // number of boundary triangles and boundary quadrilaterals for (i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (el.GetNP() == 3 && ( mesh.GetFaceDescriptor (el.GetIndex()).DomainIn()==0 || mesh.GetFaceDescriptor (el.GetIndex()).DomainOut()==0 ) ) nbtri++; else if (el.GetNP() == 4 && ( mesh.GetFaceDescriptor (el.GetIndex()).DomainIn()==0 || mesh.GetFaceDescriptor (el.GetIndex()).DomainOut()==0 ) ) nbquad++; } ofstream outfile (filename.c_str()); outfile.precision(6); outfile.setf (ios::fixed, ios::floatfield); outfile.setf (ios::showpoint); outfile << "/* <BLOBHead>\n"; outfile << "__BLOBTYPE__=Grid\n"; outfile << "__OWNER__=JCMwave\n"; outfile << "<I>SpaceDim=3\n"; outfile << "<I>ManifoldDim=3\n"; outfile << "<I>NRefinementSteps=0\n"; outfile << "<I>NPoints="<<np<<"\n"; outfile << "<I>NTetrahedra="<<ntets<<"\n"; outfile << "<I>NPrisms="<<nprisms<<"\n"; outfile << "<I>NBoundaryTriangles="<<nbtri<<"\n"; outfile << "<I>NBoundaryQuadrilaterals="<<nbquad<<"\n"; outfile << "*/\n"; outfile << "\n"; outfile << "# output from Netgen\n\n"; int nDomains=mesh.GetNDomains(); for (i=1; i<=nDomains; i++) { if (mesh.GetMaterial(i)) outfile << "#" << mesh.GetMaterial(i) << ": Material ID = " << i << "\n"; } outfile << "# Points\n"; cout << " Please note: The unit of length in the .geo file is assumed to be 'microns'."<<endl; for (i = 1; i <= np; i++) { const Point<3> & p = mesh.Point(i); outfile << i << "\n"; outfile << p(0) << "e-6\n"; outfile << p(1) << "e-6\n"; outfile << p(2) << "e-6\n\n"; } outfile << "\n"; outfile << "# Tetrahedra\n"; counter = 0; for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (el.GetNP() == 4) { counter++; dx1 = mesh.Point(el.PNum(2))(0) - mesh.Point(el.PNum(1))(0); dx2 = mesh.Point(el.PNum(3))(0) - mesh.Point(el.PNum(1))(0); dx3 = mesh.Point(el.PNum(4))(0) - mesh.Point(el.PNum(1))(0); dy1 = mesh.Point(el.PNum(2))(1) - mesh.Point(el.PNum(1))(1); dy2 = mesh.Point(el.PNum(3))(1) - mesh.Point(el.PNum(1))(1); dy3 = mesh.Point(el.PNum(4))(1) - mesh.Point(el.PNum(1))(1); dz1 = mesh.Point(el.PNum(2))(2) - mesh.Point(el.PNum(1))(2); dz2 = mesh.Point(el.PNum(3))(2) - mesh.Point(el.PNum(1))(2); dz3 = mesh.Point(el.PNum(4))(2) - mesh.Point(el.PNum(1))(2); vol = (dy1*dz2-dz1*dy2)*dx3 + (dz1*dx2-dx1*dz2)*dy3 + (dx1*dy2-dy1*dx2)*dz3; if ( vol > 0 ) for (j = 1; j <= 4; j++) outfile << el.PNum(j)<<"\n"; else { for (j = 2; j >= 1; j--) outfile << el.PNum(j)<<"\n"; for (j = 3; j <= 4; j++) outfile << el.PNum(j)<<"\n"; } outfile << el.GetIndex() << "\n\n"; } } if ( counter != ntets) { cout<< "\n Error in determining number of tetras!\n"<<endl; return; } outfile << "\n"; outfile << "# Prisms\n"; counter = 0; for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (el.GetNP() == 6) { counter++; dx1 = mesh.Point(el.PNum(2))(0) - mesh.Point(el.PNum(1))(0); dx2 = mesh.Point(el.PNum(3))(0) - mesh.Point(el.PNum(1))(0); dx3 = mesh.Point(el.PNum(4))(0) - mesh.Point(el.PNum(1))(0); dy1 = mesh.Point(el.PNum(2))(1) - mesh.Point(el.PNum(1))(1); dy2 = mesh.Point(el.PNum(3))(1) - mesh.Point(el.PNum(1))(1); dy3 = mesh.Point(el.PNum(4))(1) - mesh.Point(el.PNum(1))(1); dz1 = mesh.Point(el.PNum(2))(2) - mesh.Point(el.PNum(1))(2); dz2 = mesh.Point(el.PNum(3))(2) - mesh.Point(el.PNum(1))(2); dz3 = mesh.Point(el.PNum(4))(2) - mesh.Point(el.PNum(1))(2); vol = (dy1*dz2-dz1*dy2)*dx3 + (dz1*dx2-dx1*dz2)*dy3 + (dx1*dy2-dy1*dx2)*dz3; if (pointsOnTetras.Get(el.PNum(1)) && pointsOnTetras.Get(el.PNum(2)) && pointsOnTetras.Get(el.PNum(3))) { if (vol > 0) for (j = 1; j <= 6; j++) outfile << el.PNum(j)<<"\n"; else { for (j = 3; j >= 1; j--) outfile << el.PNum(j)<<"\n"; for (j = 6; j >= 4; j--) outfile << el.PNum(j)<<"\n"; } } else if ( pointsOnTetras.Get(el.PNum(4)) && pointsOnTetras.Get(el.PNum(5)) && pointsOnTetras.Get(el.PNum(6)) ) { if ( vol < 0 ) { for (j = 4; j <= 6; j++) outfile << el.PNum(j)<<"\n"; for (j = 1; j <= 3; j++) outfile << el.PNum(j)<<"\n"; } else { for (j = 6; j >= 4; j--) outfile << el.PNum(j)<<"\n"; for (j = 3; j >= 1; j--) outfile << el.PNum(j)<<"\n"; } } else { cout << "\n Error in determining prism point numbering!\n"<<endl; return; } outfile << el.GetIndex() << "\n\n"; } } if ( counter != nprisms) { cout<< "\n Error in determining number of prisms!\n"<<endl; return; } int npid1 = 0; int npid2 = 0; int npid3 = 0; for (i=1; i<=np; i++) { if (identmap1.Elem(i)) npid1++; if (identmap2.Elem(i)) npid2++; if (identmap3.Elem(i)) npid3++; } outfile << "\n"; outfile << "# Boundary triangles\n"; outfile << "# Number of identified points in 1-direction: " << npid1 << "\n"; outfile << "# Number of identified points in 2-direction: " << npid2 << "\n"; outfile << "# Number of identified points in 3-direction: " << npid3 << "\n"; for (i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (el.GetNP() == 3 && (mesh.GetFaceDescriptor (el.GetIndex()).DomainIn()==0 || mesh.GetFaceDescriptor (el.GetIndex()).DomainOut()==0)) { outfile <<"# T\n"; for (j = 1; j <= 3; j++) outfile << el.PNum(j)<<"\n"; if (mesh.GetFaceDescriptor (el.GetIndex()).BCProperty()==bc_at_infinity) outfile << 1000 << "\n"; else outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << "\n"; if (mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() == bc_at_infinity) outfile << "-2\n\n"; else if (identmap1.Elem(el.PNum(1)) &&identmap1.Elem(el.PNum(2)) &&identmap1.Elem(el.PNum(3))) { outfile << "-1\n"; for (j = 1; j <= 3; j++) outfile << identmap1.Elem(el.PNum(j))<<"\n"; outfile << "\n"; } else if (identmap2.Elem(el.PNum(1)) &&identmap2.Elem(el.PNum(2)) &&identmap2.Elem(el.PNum(3))) { outfile << "-1\n"; for (j = 1; j <= 3; j++) outfile << identmap2.Elem(el.PNum(j))<<"\n"; outfile << "\n"; } else if (identmap3.Elem(el.PNum(1)) &&identmap3.Elem(el.PNum(2)) &&identmap3.Elem(el.PNum(3))) { outfile << "-1\n"; for (j = 1; j <= 3; j++) outfile << identmap3.Elem(el.PNum(j))<<"\n"; outfile << "\n"; } else outfile << "1\n\n"; } } outfile << "\n"; outfile << "# Boundary quadrilaterals\n"; for (i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (el.GetNP() == 4 && (mesh.GetFaceDescriptor (el.GetIndex()).DomainIn()==0 || mesh.GetFaceDescriptor (el.GetIndex()).DomainOut()==0)) { if (pointsOnTetras.Get(el.PNum(1)) && pointsOnTetras.Get(el.PNum(2))) ct = 0; else if (pointsOnTetras.Get(el.PNum(2)) && pointsOnTetras.Get(el.PNum(3))) ct = 1; else if (pointsOnTetras.Get(el.PNum(3)) && pointsOnTetras.Get(el.PNum(4))) ct = 2; else if (pointsOnTetras.Get(el.PNum(4)) && pointsOnTetras.Get(el.PNum(1))) ct = 3; else cout << "\nWarning: Quadrilateral with inconsistent points found!"<<endl; for (j = 1; j <= 4; j++) { jj = j + ct; if ( jj >= 5 ) jj = jj - 4; outfile << el.PNum(jj)<<"\n"; } outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << "\n"; if (mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() == bc_at_infinity) { outfile << "-2\n\n"; cout << "\nWarning: Quadrilateral at infinity found (this should not occur)!"<<endl; } else if ( identmap1.Elem(el.PNum(1)) && identmap1.Elem(el.PNum(2)) && identmap1.Elem(el.PNum(3)) && identmap1.Elem(el.PNum(4)) ) { outfile << "-1\n"; for (j = 1; j <= 4; j++) { jj = j + ct; if ( jj >= 5 ) jj = jj - 4; outfile << identmap1.Elem(el.PNum(jj))<<"\n"; } outfile << "\n"; } else if ( identmap2.Elem(el.PNum(1)) && identmap2.Elem(el.PNum(2)) && identmap2.Elem(el.PNum(3)) && identmap2.Elem(el.PNum(4)) ) { outfile << "-1\n"; for (j = 1; j <= 4; j++) { jj = j + ct; if ( jj >= 5 ) jj = jj - 4; outfile << identmap2.Elem(el.PNum(jj))<<"\n"; } outfile << "\n"; } else if ( identmap3.Elem(el.PNum(1)) && identmap3.Elem(el.PNum(2)) && identmap3.Elem(el.PNum(3)) && identmap3.Elem(el.PNum(4)) ) { outfile << "-1\n"; for (j = 1; j <= 4; j++) { jj = j + ct; if ( jj >= 5 ) jj = jj - 4; outfile << identmap3.Elem(el.PNum(jj))<<"\n"; } outfile << "\n"; } else outfile << "1\n\n"; } } cout << " JCMwave grid file written." << endl; }
// extern double teterrpow; MESHING3_RESULT MeshVolume (MeshingParameters & mp, Mesh& mesh3d) { int oldne; int meshed; Array<INDEX_2> connectednodes; if (&mesh3d.LocalHFunction() == NULL) mesh3d.CalcLocalH(mp.grading); mesh3d.Compress(); // mesh3d.PrintMemInfo (cout); if (mp.checkoverlappingboundary) if (mesh3d.CheckOverlappingBoundary()) throw NgException ("Stop meshing since boundary mesh is overlapping"); int nonconsist = 0; for (int k = 1; k <= mesh3d.GetNDomains(); k++) { PrintMessage (3, "Check subdomain ", k, " / ", mesh3d.GetNDomains()); mesh3d.FindOpenElements(k); /* bool res = mesh3d.CheckOverlappingBoundary(); if (res) { PrintError ("Surface is overlapping !!"); nonconsist = 1; } */ bool res = (mesh3d.CheckConsistentBoundary() != 0); if (res) { PrintError ("Surface mesh not consistent"); nonconsist = 1; } } if (nonconsist) { PrintError ("Stop meshing since surface mesh not consistent"); throw NgException ("Stop meshing since surface mesh not consistent"); } double globmaxh = mp.maxh; for (int k = 1; k <= mesh3d.GetNDomains(); k++) { if (multithread.terminate) break; PrintMessage (2, ""); PrintMessage (1, "Meshing subdomain ", k, " of ", mesh3d.GetNDomains()); (*testout) << "Meshing subdomain " << k << endl; mp.maxh = min2 (globmaxh, mesh3d.MaxHDomain(k)); mesh3d.CalcSurfacesOfNode(); mesh3d.FindOpenElements(k); if (!mesh3d.GetNOpenElements()) continue; Box<3> domain_bbox( Box<3>::EMPTY_BOX ); for (SurfaceElementIndex sei = 0; sei < mesh3d.GetNSE(); sei++) { const Element2d & el = mesh3d[sei]; if (el.IsDeleted() ) continue; if (mesh3d.GetFaceDescriptor(el.GetIndex()).DomainIn() == k || mesh3d.GetFaceDescriptor(el.GetIndex()).DomainOut() == k) for (int j = 0; j < el.GetNP(); j++) domain_bbox.Add (mesh3d[el[j]]); } domain_bbox.Increase (0.01 * domain_bbox.Diam()); for (int qstep = 1; qstep <= 3; qstep++) { // cout << "openquads = " << mesh3d.HasOpenQuads() << endl; if (mesh3d.HasOpenQuads()) { string rulefile = ngdir; const char ** rulep = NULL; switch (qstep) { case 1: rulefile += "/rules/prisms2.rls"; rulep = prismrules2; break; case 2: // connect pyramid to triangle rulefile += "/rules/pyramids2.rls"; rulep = pyramidrules2; break; case 3: // connect to vis-a-vis point rulefile += "/rules/pyramids.rls"; rulep = pyramidrules; break; } // Meshing3 meshing(rulefile); Meshing3 meshing(rulep); MeshingParameters mpquad = mp; mpquad.giveuptol = 15; mpquad.baseelnp = 4; mpquad.starshapeclass = 1000; mpquad.check_impossible = qstep == 1; // for prisms only (air domain in trafo) for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++) meshing.AddPoint (mesh3d[pi], pi); mesh3d.GetIdentifications().GetPairs (0, connectednodes); for (int i = 1; i <= connectednodes.Size(); i++) meshing.AddConnectedPair (connectednodes.Get(i)); for (int i = 1; i <= mesh3d.GetNOpenElements(); i++) { Element2d hel = mesh3d.OpenElement(i); meshing.AddBoundaryElement (hel); } oldne = mesh3d.GetNE(); meshing.GenerateMesh (mesh3d, mpquad); for (int i = oldne + 1; i <= mesh3d.GetNE(); i++) mesh3d.VolumeElement(i).SetIndex (k); (*testout) << "mesh has " << mesh3d.GetNE() << " prism/pyramid elements" << endl; mesh3d.FindOpenElements(k); } } if (mesh3d.HasOpenQuads()) { PrintSysError ("mesh has still open quads"); throw NgException ("Stop meshing since too many attempts"); // return MESHING3_GIVEUP; } if (mp.delaunay && mesh3d.GetNOpenElements()) { Meshing3 meshing((const char**)NULL); mesh3d.FindOpenElements(k); /* for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++) meshing.AddPoint (mesh3d[pi], pi); */ for (PointIndex pi : mesh3d.Points().Range()) meshing.AddPoint (mesh3d[pi], pi); for (int i = 1; i <= mesh3d.GetNOpenElements(); i++) meshing.AddBoundaryElement (mesh3d.OpenElement(i)); oldne = mesh3d.GetNE(); meshing.Delaunay (mesh3d, k, mp); for (int i = oldne + 1; i <= mesh3d.GetNE(); i++) mesh3d.VolumeElement(i).SetIndex (k); PrintMessage (3, mesh3d.GetNP(), " points, ", mesh3d.GetNE(), " elements"); } int cntsteps = 0; if (mesh3d.GetNOpenElements()) do { if (multithread.terminate) break; mesh3d.FindOpenElements(k); PrintMessage (5, mesh3d.GetNOpenElements(), " open faces"); cntsteps++; if (cntsteps > mp.maxoutersteps) throw NgException ("Stop meshing since too many attempts"); string rulefile = ngdir + "/tetra.rls"; PrintMessage (1, "start tetmeshing"); // Meshing3 meshing(rulefile); Meshing3 meshing(tetrules); Array<int, PointIndex::BASE> glob2loc(mesh3d.GetNP()); glob2loc = -1; for (PointIndex pi = mesh3d.Points().Begin(); pi < mesh3d.Points().End(); pi++) if (domain_bbox.IsIn (mesh3d[pi])) glob2loc[pi] = meshing.AddPoint (mesh3d[pi], pi); for (int i = 1; i <= mesh3d.GetNOpenElements(); i++) { Element2d hel = mesh3d.OpenElement(i); for (int j = 0; j < hel.GetNP(); j++) hel[j] = glob2loc[hel[j]]; meshing.AddBoundaryElement (hel); // meshing.AddBoundaryElement (mesh3d.OpenElement(i)); } oldne = mesh3d.GetNE(); mp.giveuptol = 15 + 10 * cntsteps; mp.sloppy = 5; meshing.GenerateMesh (mesh3d, mp); for (ElementIndex ei = oldne; ei < mesh3d.GetNE(); ei++) mesh3d[ei].SetIndex (k); mesh3d.CalcSurfacesOfNode(); mesh3d.FindOpenElements(k); // teterrpow = 2; if (mesh3d.GetNOpenElements() != 0) { meshed = 0; PrintMessage (5, mesh3d.GetNOpenElements(), " open faces found"); MeshOptimize3d optmesh(mp); const char * optstr = "mcmstmcmstmcmstmcm"; for (size_t j = 1; j <= strlen(optstr); j++) { mesh3d.CalcSurfacesOfNode(); mesh3d.FreeOpenElementsEnvironment(2); mesh3d.CalcSurfacesOfNode(); switch (optstr[j-1]) { case 'c': optmesh.CombineImprove(mesh3d, OPT_REST); break; case 'd': optmesh.SplitImprove(mesh3d, OPT_REST); break; case 's': optmesh.SwapImprove(mesh3d, OPT_REST); break; case 't': optmesh.SwapImprove2(mesh3d, OPT_REST); break; case 'm': mesh3d.ImproveMesh(mp, OPT_REST); break; } } mesh3d.FindOpenElements(k); PrintMessage (3, "Call remove problem"); RemoveProblem (mesh3d, k); mesh3d.FindOpenElements(k); } else { meshed = 1; PrintMessage (1, "Success !"); } } while (!meshed); PrintMessage (1, mesh3d.GetNP(), " points, ", mesh3d.GetNE(), " elements"); } mp.maxh = globmaxh; MeshQuality3d (mesh3d); return MESHING3_OK; }
void WriteDolfinFormat (const Mesh & mesh, const string & filename) { cout << "start writing dolfin export" << endl; int np = mesh.GetNP(); int ne = mesh.GetNE(); int nse = mesh.GetNSE(); int nsd = mesh.GetDimension(); ofstream outfile(filename.c_str()); outfile.precision(8); outfile.setf(ios::fixed, ios::floatfield); outfile.setf(ios::showpoint); if (nsd == 3) { outfile << "<?xml version=\"1.0\" encoding=\"UTF-8\"?>\n"; outfile << "<dolfin xmlns:dolfin=\"http://www.fenics.org/wiki/DOLFIN\">\n"; outfile << " <mesh celltype=\"tetrahedron\" dim=\"3\">\n"; outfile << " <vertices size=\""<<np<<"\">\n"; for (int i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << " <vertex index=\""<<i-1<<"\" x=\""<<p.X()<<"\" y=\""<<p.Y()<<"\" z=\""<<p.Z()<<"\"/>\n"; } outfile << " </vertices>\n"; outfile << " <cells size=\""<<ne<<"\">\n"; for (int i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (mparam.inverttets) el.Invert(); ELEMENT_TYPE type = el.GetType(); if (type == TET || type == TET10) outfile << " <tetrahedron index=\""<<i-1<<"\" v0=\""<<el.PNum(1)-1<<"\" v1=\""<<el.PNum(2)-1<<"\" v2=\""<<el.PNum(3)-1<<"\" v3=\""<<el.PNum(4)-1<<"\"/>\n"; else if (type == HEX) outfile << " <hexahedron index=\""<<i-1<<"\" v0=\""<<el.PNum(1)-1<<"\" v1=\""<<el.PNum(2)-1<<"\" v2=\""<<el.PNum(3)-1<<"\" v3=\""<<el.PNum(4)-1 <<"\" v4=\""<<el.PNum(5)-1<<"\" v5=\""<<el.PNum(6)-1<<"\" v6=\""<<el.PNum(7)-1<<"\" v7=\""<<el.PNum(8)-1<<"\"/>\n"; else cout << "Warning unsupported element." << endl; } outfile << " </cells>\n"; outfile << " </mesh>\n"; outfile << "</dolfin>"; } else if (nsd == 2) { outfile << "<?xml version=\"1.0\" encoding=\"UTF-8\"?>\n"; outfile << "<dolfin xmlns:dolfin=\"http://www.fenics.org/wiki/DOLFIN\">\n"; outfile << " <mesh celltype=\"triangle\" dim=\"2\">\n"; outfile << " <vertices size=\""<<np<<"\">\n"; for (int i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << " <vertex index=\""<<i-1<<"\" x=\""<<p.X()<<"\" y=\""<<p.Y()<<"\"/>\n"; } outfile << " </vertices>\n"; outfile << " <cells size=\""<<nse<<"\">\n"; for (int i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (mparam.inverttrigs) el.Invert(); ELEMENT_TYPE type = el.GetType(); if (type == TRIG || type == TRIG6) outfile << " <triangle index=\""<<i-1<<"\" v0=\""<<el.PNum(1)-1<<"\" v1=\""<<el.PNum(2)-1<<"\" v2=\""<<el.PNum(3)-1<<"\"/>\n"; else if (type == QUAD || type == QUAD8) outfile << " <quadilateral index=\""<<i-1<<"\" v0=\""<<el.PNum(1)-1<<"\" v1=\""<<el.PNum(2)-1<<"\" v2=\""<<el.PNum(3)-1<<"\" v3=\""<<el.PNum(4)-1<<"\"/>\n"; else cout << "Warning unsupported element." << endl; } outfile << " </cells>\n"; outfile << " </mesh>\n"; outfile << "</dolfin>"; } outfile.close(); cout << "done writing dolfin export" << endl; }
void WriteGmshFormat (const Mesh & mesh, const CSGeometry & geom, const string & filename) { ofstream outfile (filename.c_str()); outfile.precision(6); outfile.setf (ios::fixed, ios::floatfield); outfile.setf (ios::showpoint); int np = mesh.GetNP(); /// number of point int ne = mesh.GetNE(); /// number of element int nse = mesh.GetNSE(); /// number of surface element (BC) int i, j, k, l; /* * 3D section : Linear volume elements (only tetrahedra) */ if (ne > 0 && mesh.VolumeElement(1).GetNP() == 4) { cout << "Write GMSH Format \n"; cout << "The GMSH format is available for linear tetrahedron elements only in 3D\n" << endl; int inverttets = mparam.inverttets; int invertsurf = mparam.inverttrigs; /// Write nodes outfile << "$NOD\n"; outfile << np << "\n"; for (i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << i << " "; /// node number outfile << p.X() << " "; outfile << p.Y() << " "; outfile << p.Z() << "\n"; } outfile << "$ENDNOD\n"; /// write elements outfile << "$ELM\n"; outfile << ne + nse << "\n"; //// number of elements + number of surfaces BC for (i = 1; i <= nse; i++) { Element2d el = mesh.SurfaceElement(i); if (invertsurf) el.Invert(); outfile << i; outfile << " "; outfile << "2"; outfile << " "; outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; /// that means that physical entity = elementary entity (arbitrary approach) outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; outfile << "3"; outfile << " "; for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile << el.PNum(j); } outfile << "\n"; } for (i = 1; i <= ne; i++) { Element el = mesh.VolumeElement(i); if (inverttets) el.Invert(); outfile << nse + i; /// element number outfile << " "; outfile << "4"; /// element type i.e. Tetraedron == 4 outfile << " "; outfile << 100000 + el.GetIndex(); /// that means that physical entity = elementary entity (arbitrary approach) outfile << " "; outfile << 100000 + el.GetIndex(); /// volume number outfile << " "; outfile << "4"; /// number of nodes i.e. 4 for a tetrahedron for (j = 1; j <= el.GetNP(); j++) { outfile << " "; outfile << el.PNum(j); } outfile << "\n"; } outfile << "$ENDELM\n"; } /* * End of 3D section */ /* * 2D section : available for triangles and quadrangles */ else if (ne == 0) /// means that there's no 3D element { cout << "\n Write Gmsh Surface Mesh (triangle and/or quadrangles)" << endl; /// Write nodes outfile << "$NOD\n"; outfile << np << "\n"; for (i = 1; i <= np; i++) { const Point3d & p = mesh.Point(i); outfile << i << " "; /// node number outfile << p.X() << " "; outfile << p.Y() << " "; outfile << p.Z() << "\n"; } outfile << "$ENDNOD\n"; /// write triangles & quadrangles outfile << "$ELM\n"; outfile << nse << "\n"; for (k = 1; k <= nse; k++) { const Element2d & el = mesh.SurfaceElement(k); outfile << k; outfile << " "; outfile << (el.GetNP()-1); // 2 for a triangle and 3 for a quadrangle outfile << " "; outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; /// that means that physical entity = elementary entity (arbitrary approach) outfile << mesh.GetFaceDescriptor (el.GetIndex()).BCProperty() << " "; outfile << (el.GetNP()); // number of node per surfacic element outfile << " "; for (l = 1; l <= el.GetNP(); l++) { outfile << " "; outfile << el.PNum(l); } outfile << "\n"; } outfile << "$ENDELM$ \n"; } /* * End of 2D section */ else { cout << " Invalide element type for Gmsh volume Format !\n"; } }
void Refinement :: MakeSecondOrder (Mesh & mesh) { int nseg, nse, ne; mesh.ComputeNVertices(); mesh.SetNP(mesh.GetNV()); INDEX_2_HASHTABLE<int> between(mesh.GetNP() + 5); bool thinlayers = 0; for (ElementIndex ei = 0; ei < mesh.GetNE(); ei++) if (mesh[ei].GetType() == PRISM || mesh[ei].GetType() == PRISM12) thinlayers = 1; nseg = mesh.GetNSeg(); for (SegmentIndex si = 0; si < nseg; si++) { Segment & el = mesh.LineSegment(si); INDEX_2 i2 = INDEX_2::Sort (el[0], el[1]); if (between.Used(i2)) el[2] = between.Get(i2); else { Point<3> pb; EdgePointGeomInfo ngi; PointBetween (mesh.Point (el[0]), mesh.Point (el[1]), 0.5, el.surfnr1, el.surfnr2, el.epgeominfo[0], el.epgeominfo[1], pb, ngi); el[2] = mesh.AddPoint (pb, mesh.Point(el[0]).GetLayer(), EDGEPOINT); between.Set (i2, el[2]); } } // refine surface elements nse = mesh.GetNSE(); for (SurfaceElementIndex sei = 0; sei < nse; sei++) { int j; const Element2d & el = mesh.SurfaceElement(sei); int onp(0); Element2d newel; newel.SetIndex (el.GetIndex()); static int betw_trig[3][3] = { { 1, 2, 3 }, { 0, 2, 4 }, { 0, 1, 5 } }; static int betw_quad6[2][3] = { { 0, 1, 4 }, { 3, 2, 5 } }; static int betw_quad8[4][3] = { { 0, 1, 4 }, { 3, 2, 5 }, { 0, 3, 6 }, { 1, 2, 7 } }; int (*betw)[3] = NULL; switch (el.GetType()) { case TRIG: case TRIG6: { betw = betw_trig; newel.SetType (TRIG6); onp = 3; break; } case QUAD: case QUAD6: case QUAD8: { if (thinlayers) { betw = betw_quad6; newel.SetType (QUAD6); } else { betw = betw_quad8; newel.SetType (QUAD8); } onp = 4; break; } default: PrintSysError ("Unhandled element in secondorder:", int(el.GetType())); } for (j = 0; j < onp; j++) newel[j] = el[j]; int nnp = newel.GetNP(); for (j = 0; j < nnp-onp; j++) { int pi1 = newel[betw[j][0]]; int pi2 = newel[betw[j][1]]; INDEX_2 i2 = INDEX_2::Sort (pi1, pi2); if (between.Used(i2)) newel[onp+j] = between.Get(i2); else { Point<3> pb; PointGeomInfo newgi; PointBetween (mesh.Point (pi1), mesh.Point (pi2), 0.5, mesh.GetFaceDescriptor(el.GetIndex ()).SurfNr(), el.GeomInfoPi (betw[j][0]+1), el.GeomInfoPi (betw[j][1]+1), pb, newgi); newel[onp+j] = mesh.AddPoint (pb, mesh.Point(pi1).GetLayer(), SURFACEPOINT); between.Set (i2, newel[onp+j]); } } mesh.SurfaceElement(sei) = newel; } // int i, j; // refine volume elements ne = mesh.GetNE(); for (int i = 1; i <= ne; i++) { const Element & el = mesh.VolumeElement(i); int onp = 0; Element newel; newel.SetIndex (el.GetIndex()); static int betw_tet[6][3] = { { 0, 1, 4 }, { 0, 2, 5 }, { 0, 3, 6 }, { 1, 2, 7 }, { 1, 3, 8 }, { 2, 3, 9 } }; static int betw_prism[6][3] = { { 0, 2, 6 }, { 0, 1, 7 }, { 1, 2, 8 }, { 3, 5, 9 }, { 3, 4, 10 }, { 4, 5, 11 }, }; int (*betw)[3] = NULL; switch (el.GetType()) { case TET: case TET10: { betw = betw_tet; newel.SetType (TET10); onp = 4; break; } case PRISM: case PRISM12: { betw = betw_prism; newel.SetType (PRISM12); onp = 6; break; } default: PrintSysError ("MakeSecondOrder, illegal vol type ", el.GetType()); } for (int j = 1; j <= onp; j++) newel.PNum(j) = el.PNum(j); int nnp = newel.GetNP(); for (int j = 0; j < nnp-onp; j++) { INDEX_2 i2(newel[betw[j][0]], newel[betw[j][1]]); i2.Sort(); if (between.Used(i2)) newel.PNum(onp+1+j) = between.Get(i2); else { newel.PNum(onp+1+j) = mesh.AddPoint (Center (mesh.Point(i2.I1()), mesh.Point(i2.I2())), mesh.Point(i2.I1()).GetLayer(), INNERPOINT); between.Set (i2, newel.PNum(onp+1+j)); } } mesh.VolumeElement (i) = newel; } // makes problems after linear mesh refinement, since // 2nd order identifications are not removed // update identification tables for (int i = 1; i <= mesh.GetIdentifications().GetMaxNr(); i++) { Array<int,PointIndex::BASE> identmap; mesh.GetIdentifications().GetMap (i, identmap); for (INDEX_2_HASHTABLE<int>::Iterator it = between.Begin(); it != between.End(); it++) { INDEX_2 i2; int newpi; between.GetData (it, i2, newpi); INDEX_2 oi2(identmap.Get(i2.I1()), identmap.Get(i2.I2())); oi2.Sort(); if (between.Used (oi2)) { int onewpi = between.Get(oi2); mesh.GetIdentifications().Add (newpi, onewpi, i); } } /* for (int j = 1; j <= between.GetNBags(); j++) for (int k = 1; k <= between.GetBagSize(j); k++) { INDEX_2 i2; int newpi; between.GetData (j, k, i2, newpi); INDEX_2 oi2(identmap.Get(i2.I1()), identmap.Get(i2.I2())); oi2.Sort(); if (between.Used (oi2)) { int onewpi = between.Get(oi2); mesh.GetIdentifications().Add (newpi, onewpi, i); } } */ } // mesh.mglevels++; int oldsize = mesh.mlbetweennodes.Size(); mesh.mlbetweennodes.SetSize(mesh.GetNP()); for (int i = oldsize; i < mesh.GetNP(); i++) mesh.mlbetweennodes[i] = INDEX_2(0,0); /* for (i = 1; i <= between.GetNBags(); i++) for (j = 1; j <= between.GetBagSize(i); j++) { INDEX_2 oldp; int newp; between.GetData (i, j, oldp, newp); mesh.mlbetweennodes.Elem(newp) = oldp; } */ for (INDEX_2_HASHTABLE<int>::Iterator it = between.Begin(); it != between.End(); it++) { mesh.mlbetweennodes[between.GetData (it)] = between.GetHash(it); } mesh.ComputeNVertices(); mesh.RebuildSurfaceElementLists(); // ValidateSecondOrder (mesh); }