1 /*
2 XLiFE++ is an extended library of finite elements written in C++
3     Copyright (C) 2014  Lunéville, Eric; Kielbasiewicz, Nicolas; Lafranche, Yvon; Nguyen, Manh-Ha; Chambeyron, Colin
4 
5     This program is free software: you can redistribute it and/or modify
6     it under the terms of the GNU General Public License as published by
7     the Free Software Foundation, either version 3 of the License, or
8     (at your option) any later version.
9     This program is distributed in the hope that it will be useful,
10     but WITHOUT ANY WARRANTY; without even the implied warranty of
11     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
12     GNU General Public License for more details.
13     You should have received a copy of the GNU General Public License
14     along with this program.  If not, see <http://www.gnu.org/licenses/>.
15  */
16 
17 /*!
18   \file Mesh.hpp
19   \authors D. Martin, Y. Lafranche, E. Lunéville, N. Kielbasiewicz
20   \since 26 may 2010
21   \date 17 oct 2014
22 
23   \brief Definition of the xlifepp::Mesh class
24 
25   xlifepp::Mesh represents the main class to deal with mesh data
26   either created by xlife++ meshing tools or read from mesh files
27   generated by external mesh tools (gmsh, TetMesh-GHS3D, ...)
28 
29   This xlifepp::Mesh class supports any kind of 1D, 2D, 3D meshes made of
30     - segments, triangles, quadrangles, tetraedron, hexaedron, prisms, pyramids
31     - of order 1 (plane elements) 2, 3, ... (curved elements)
32     - elements of same dimension can be mixed
33 
34   These elements are GeomElements object in xlife++ terminology
35   and should not be confused with Element object supporting finite element interpolation in addition
36   Note that the dimension of elements may be less than the space dimension
37 
38   A Mesh consists of :
39     - a vector of points : vertices of elements and additonnal nodes for element of order greater than 1
40     - a list of xlifepp::GeomElement : list of vertices, a xlifepp::RefElement describing geometrical interpolation, ...
41     - a list of sides (edge in 2D, faces in 3D) seen as xlifepp::GeomElement (xlifepp::MeshSideElement)
42     - a list of side of sides (edges in 3D) seen as xlifepp::GeomElement (xlifepp::MeshSideElement)
43     - a list of vertices indices
44     - a xlifepp::Geometry object collecting global geometric data such as bounding box, boundary parametrisations,...
45     - a list of xlifepp::GeomDomain describing physical meaning part of the mesh, in particular boundaries
46 
47   Contrary to a xlifepp::MeshElement which carries a lot of data, a MeshSideElement contains only its parents (GeomElement)
48   and the related side numbers
49 
50   The class manages a pointer (firstOrderMesh_p) to its underlying first order mesh. If the mesh is itself a first order mesh
51   (meaning that all elements are first order elements), firstOrderMesh_p is equal to the mesh pointer (this), else firstOrderMesh_p=0
52   until its underlying first order mesh is built. The underlying first order mesh is required by some tools visualization
53   and therefore it is generally built by mesh export functions.
54 */
55 
56 #ifndef MESH_HPP
57 #define MESH_HPP
58 
59 #include "config.h"
60 #include "GeomElement.hpp"
61 #include "GeomDomain.hpp"
62 #include "Geometry.hpp"
63 #include "geometries_utils.hpp"
64 #include "geometries1D.hpp"
65 #include "geometries2D.hpp"
66 #include "geometries3D.hpp"
67 #include "subdivision/subUtil/HexahedronMesh.hpp"
68 #include "subdivision/subUtil/TetrahedronMesh.hpp"
69 #include "subdivision/subUtil/QuadrangleMesh.hpp"
70 #include "subdivision/subUtil/TriangleMesh.hpp"
71 #include "utils.h"
72 
73 namespace xlifepp
74 {
75 
76 /*!
77   \class CrackData
78   store data related to a crack : domain name, domain id and dimension
79  */
80 
81 class CrackData
82 {
83   public:
84     string_t domainName;    //!< name of the domain supporting the crack
85     number_t id;            //!< id of the domain supporting the crack
86     number_t dim;           //!< dimension of the crack
87     MinimalBox minimalBox;  //!< minimal box of the geometry supporting the crack
88     number_t domIdToOpen;   //!< stricly positive is open crack
CrackData()89     CrackData()             //! default constructor
90     {
91       domainName="Crack";
92       id=1;
93       dim=1;
94       domIdToOpen=0;
95       minimalBox=MinimalBox();
96     }
97     //! constructor
CrackData(string_t name,number_t i,number_t d,const MinimalBox & mb,number_t domId=0)98     CrackData(string_t name, number_t i, number_t d, const MinimalBox& mb, number_t domId = 0)
99     : domainName(name), id(i), dim(d), minimalBox(mb), domIdToOpen(domId) {}
100     bool operator<(const CrackData& cd) const; //!< comparison operator
101     std::vector<real_t> computeNormal() const; //!< compute the normal vector for a 1D (in plane) or 2D (in space) cracked domain
102     void print(std::ostream&) const;                        //!< print mesh data
print(PrintStream & os) const103     void print(PrintStream& os) const {print(os.currentStream());}
104     friend std::ostream& operator<<(std::ostream&, const CrackData&);
105 };
106 
107 std::ostream& operator<<(std::ostream&, const CrackData&); //!< print operator
108 
109 //! finds a crack data in a list
110 short int findId(std::vector<CrackData>::const_iterator it_b, std::vector<CrackData>::const_iterator it_e, number_t id);
111 
112 /*!
113   \class Mesh
114   handles mesh data either created by internal mesh tools
115   or read from file generated by external mesh tools
116  */
117 class Mesh
118 {
119   public:
120     Geometry* geometry_p;                       //!< global geometric information (space dimension, variable names, bounding box, ...)
121     std::set<CrackData> cracks;                 //!< list of domains to be cracked (with gmsh generator)
122     std::vector<Point> nodes;                   //!< list of all mesh nodes (vertices and other nodes)
123     mutable number_t lastIndex_;                //!< last index of GeomElements (= max of id numbers over all GeomElements in the mesh)
124 
125   protected:
126     string_t name_;                             //!< a name for documentation purposes
127     string_t comment_;                          //!< comment documentation purposes ("load from file data.msh",...)
128     std::vector<GeomElement*> elements_;        //!< list of geometric elements (MeshElement)
129     std::vector<GeomDomain*> domains_;          //!< list of geometric domains (part of the mesh)
130     std::vector<number_t> vertices_;            //!< list of vertices indexed by their number in nodes vector
131     bool isMadeOfSimplices_;                    //!< only composed of segments (1D mesh), triangles (2D mesh) or tetrahedra (3D mesh)
132     dimen_t order_;                             //!< order of the mesh = max of element orders
133 
134     // built if required
135     std::vector<GeomElement*> sides_;           //!< list of sides (edges in 2D, face in 3D), built if required
136     std::vector<GeomElement*> sideOfSides_;     //!< list of sides of sides (edges in 3D), built if required
137     std::vector<std::vector<GeoNumPair> > vertexElements_;  //!< for each vertex v, list of elements having vertex v, built if required
138     mutable Mesh* firstOrderMesh_p;             //!< underlying first order mesh when mesh has elements of order greater than one
139 
140   public:
141     // constructors
142     Mesh(); //!< void constructor
143     //! Mesh from a file of MeshFormat type
144     Mesh(const string_t& filename, const string_t& na, IOFormat mft = _undefFormat, number_t nodesDim=0);
145     void copyAllButNodes(const Mesh&); //!< deep copy of all mesh attributes, nodes excepted
146     void copy(const Mesh&);            //!< deep copy of all mesh attributes
Mesh(const Mesh & m)147     Mesh(const Mesh& m) { copy(m); }   //!< copy constructor
148     Mesh& operator=(const Mesh&);      //!< assign operator
149     void clear();                      //!< clear all data and delete pointers
~Mesh()150     ~Mesh() { clear(); }               //!< destructor
151 
152     void build1DMesh(const Geometry& g, number_t order, MeshGenerator mg);
153     void buildMesh(const Geometry& g, ShapeType sh, number_t order, MeshGenerator mg, std::set<MeshOption> mos);
154     //! constructor from any 1D geometry
Mesh(const Geometry & g,number_t order=1,MeshGenerator mg=_defaultGenerator,const string_t & name="")155     Mesh(const Geometry& g, number_t order = 1, MeshGenerator mg = _defaultGenerator, const string_t& name = "")
156     : name_(name) { build1DMesh(g, order, mg); }
157     //! constructor from any 2D and 3D geometries
Mesh(const Geometry & g,ShapeType sh,number_t order=1,MeshGenerator mg=_defaultGenerator,const string_t & name="")158     Mesh(const Geometry& g, ShapeType sh, number_t order = 1, MeshGenerator mg = _defaultGenerator, const string_t& name = "")
159     : name_(name) {
160       std::set<MeshOption> mos; mos.insert(_defaultMeshOption);
161       buildMesh(g, sh, order, mg, mos);
162     }
163     //! constructor from any 2D and 3D geometries
Mesh(const Geometry & g,ShapeType sh,number_t order,MeshGenerator mg,MeshOption mo,const string_t & name="")164     Mesh(const Geometry& g, ShapeType sh, number_t order, MeshGenerator mg, MeshOption mo, const string_t& name = "")
165     : name_(name) {
166       std::set<MeshOption> mos; mos.insert(mo);
167       buildMesh(g, sh, order, mg, mos);
168     }
169     //! constructor from any 2D and 3D geometries
Mesh(const Geometry & g,ShapeType sh,number_t order,MeshGenerator mg,MeshOption mo1,MeshOption mo2,const string_t & name="")170     Mesh(const Geometry& g, ShapeType sh, number_t order, MeshGenerator mg, MeshOption mo1, MeshOption mo2, const string_t& name = "")
171     : name_(name) {
172       std::set<MeshOption> mos; mos.insert(mo1); mos.insert(mo2);
173       buildMesh(g, sh, order, mg, mos);
174     }
175     //! constructor from any 2D and 3D geometries
Mesh(const Geometry & g,ShapeType sh,number_t order,MeshGenerator mg,MeshOption mo1,MeshOption mo2,MeshOption mo3,const string_t & name="")176     Mesh(const Geometry& g, ShapeType sh, number_t order, MeshGenerator mg, MeshOption mo1, MeshOption mo2, MeshOption mo3,
177          const string_t& name = "")
178     : name_(name) {
179       std::set<MeshOption> mos; mos.insert(mo1); mos.insert(mo2); mos.insert(mo3);
180       buildMesh(g, sh, order, mg, mos);
181     }
182     //! constructor from a mesh to be subdivided
183     Mesh(const Mesh& msh, number_t nbsubdiv, number_t order = 1);
184     //! elementary mesh constructor with only one reference element for test purpose
185     Mesh(ShapeType shape);
186     //! constructor from extrusion
Mesh(const Mesh & ms,const Point & O,const Point & d,number_t nbl,number_t namingDomain=0,number_t namingSection=0,number_t namingSide=0,const string_t & meshName="")187     Mesh(const Mesh& ms, const Point& O, const Point& d , number_t nbl,
188          number_t namingDomain=0, number_t namingSection=0, number_t namingSide=0,
189          const string_t& meshName="")
190     : name_(meshName) {buildExtrusion(ms,O,d,nbl,namingDomain, namingSection, namingSide);}
191     //! build new mesh converting type of element of a given mesh, if possible
192     Mesh(const Mesh& mesh, ShapeType sh, const string_t name="");
193 
194     //! choose mesh generator from Geometry and element shape
195     MeshGenerator defaultMeshGenerator(const Geometry&, ShapeType = _noShape);
196 
197     //construct additionnal data on demand
createSideIndex(std::map<string_t,std::vector<GeoNumPair>> & sideIndex) const198     void createSideIndex(std::map<string_t,std::vector<GeoNumPair> >& sideIndex) const //! index all the sides of mesh
199     { sideIndex.clear(); xlifepp::createSideIndex(elements_, sideIndex);}
200     void createSideEltIndex(std::map<string_t, GeomElement*>& sideEltIndex) const;     //!< index all the side elements of mesh
201     void buildSides();                //!< build the sides elements
202     void buildSideOfSides();          //!< build the sides of sides elements
203     void buildVertexElements();       //!< build for each vertex v, the list of elements having v as vertex
204     void buildGeomData();             //!< compute element measures and orientation
205     void setShapeTypes();             //!< set the shape types for all mesh domains
206     void mergeDomainsWithSameName();  //!< merge all domains having the same name (internal tool)
207     void updateSideDomains();         //!< update parent sides of elements involved in all side domains
208 
209     //================================================
210     //                   accessors
211     //================================================
elements() const212     const std::vector<GeomElement*>& elements() const       //! access to the elements list
213       { return elements_; }
domains() const214     const std::vector<GeomDomain*>& domains()   const       //! access to the domains list
215       { return domains_; }
216     const GeomDomain& domain(const string_t&) const;        //!< access to a domain by its name (const)
217     const GeomDomain& domain(number_t) const;               //!< access to a domain by its number (const)
218     GeomDomain& domain(const string_t&);                    //!< access to a domain by its name (no const)
219     GeomDomain& domain(number_t);                           //!< access to a domain by its number (no const)
220     number_t domainNumber(const string_t&) const;           //!< access to a domain number by its name (const)
221     number_t domainNumber(const string_t&);                 //!< access to a domain number by its name (no const)
222     string_t domainName(number_t) const;                    //!< access to a domain name by its number (const)
223     string_t domainName(number_t);                          //!< access to a domain name by its number (no const)
sides() const224     const std::vector<GeomElement*>& sides()const           //! access to the sides list
225       { return sides_; }
sideOfSides() const226     const std::vector<GeomElement*>& sideOfSides()const     //! access to the sides of sides list
227       { return sideOfSides_; }
vertices() const228     const std::vector<number_t>& vertices()const            //! access to the vertices list
229       { return vertices_; }
vertexElements() const230     const std::vector<std::vector<GeoNumPair> >& vertexElements() const //! access to the vertexElements list
231       { return vertexElements_; }
name() const232     const string_t& name() const                            //! access to the mesh name
233       { return name_; }
name(const string_t & n)234     void name(const string_t& n)                                     //! set to the mesh name
235       { name_=n; }
isMadeOfSimplices() const236     bool isMadeOfSimplices()const                           //! true if the mesh is composed only with simplices
237       { return isMadeOfSimplices_; }
element(number_t i) const238     const GeomElement& element(number_t i) const            //! access to the i-th element (i>0; no control)
239       { return *elements_[i - 1]; }
domains(number_t i) const240     const GeomDomain& domains(number_t i) const             //! access to the i-th domain (i>0; no control)
241       { return *domains_[i - 1]; }
side(number_t i) const242     const GeomElement& side(number_t i) const               //! access to the i-th side (i>0, no control)
243       { return *sides_[i - 1]; }
sideOfSide(number_t i) const244     const GeomElement& sideOfSide(number_t i) const         //! access to the i-th side of side (i>0, no control)
245       { return *sideOfSides_[i - 1]; }
vertex(number_t i) const246     number_t vertex(number_t i) const                       //! access to the i-th vertex as a MeshSideElement (i>0, no control)
247       { return vertices_[i - 1]; }
spaceDim() const248     dimen_t spaceDim() const                                //! return the dimension of physical space (dimension of node)
249       { if(nodes.size()>0) return nodes[0].size(); else return 0; }
geometryDim() const250     dimen_t geometryDim() const                             //! return the dimension of geometry
251       { return geometry_p->dim(); }
meshDim() const252     dimen_t meshDim() const                                 //! return the dimension of mesh element (may be less than space dimension)
253       { return elements_[0]->elementDim(); }
nbOfElements() const254     number_t nbOfElements() const                           //!return the number of elements
255       { return elements_.size(); }
nbOfDomains() const256     number_t nbOfDomains() const                            //! return the number of domains
257       { return domains_.size(); }
nbOfSides() const258     number_t nbOfSides() const                              //! return the number of sides
259       { return sides_.size(); }
nbOfSideOfSides() const260     number_t nbOfSideOfSides() const                        //! return the number of side of sides
261       { return sideOfSides_.size(); }
nbOfVertices() const262     number_t nbOfVertices() const                           //! return the number of vertices
263       { return vertices_.size(); }
nbOfNodes() const264     number_t nbOfNodes() const                              //! return the number of nodes
265       { return nodes.size(); }
order() const266     dimen_t order() const                                   //! return the order of the mesh (max of element orders)
267       { return order_; }
firstOrderMesh()268     Mesh& firstOrderMesh()                                  //! return underlying first order mesh as reference
269       { return *firstOrderMesh_p; }
firstOrderMesh() const270     const Mesh& firstOrderMesh() const                      //! return underlying first order mesh as const reference
271       { return *firstOrderMesh_p; }
272 
273     //================================================
274     //                  mesh tools
275     //================================================
276     Mesh* createFirstOrderMesh() const;                       //!< create underlying first order mesh if necessary
277     void merge(const Mesh&, bool mergeSharedBoundary = true, const string_t& nd="main domain"); //!< "add" mesh to the current one
278 
renameDomain(number_t n,const string_t & newname)279     void renameDomain(number_t n, const string_t& newname)               //! rename a domain given by its number
280       { domain(n).rename(newname); }
renameDomain(const string_t & oldname,const string_t & newname)281     void renameDomain(const string_t& oldname, const string_t& newname)  //! rename a domain given by its name
282       { domain(oldname).rename(newname); }
283     void removeDomain(const string_t&);                       //!< remove a domain given by its name, use with cautious
284 
285     //================================================
286     //               print facilities
287     //================================================
288     void printInfo(std::ostream& os=std::cout) const;            //!< print general informations about what the mesh is made of
printInfo(PrintStream & os) const289     void printInfo(PrintStream& os) const
290      {printInfo(os.currentStream());}
291     void print(std::ostream&) const;                             //!< print mesh data
print(PrintStream & os) const292     void print(PrintStream& os) const {print(os.currentStream());}
293     void addSuffix(const string_t& s);                           //!< add a suffix to all names (mesh name and domain names)
294     friend std::ostream& operator<<(std::ostream&, const Mesh&);
295 
296     //================================================
297     //           input/output facilities
298     //================================================
299     void saveToFile(const string_t &, IOFormat mf=_undefFormat, bool withDomains= false) const; //!< save mesh to file
300     void saveToVtk(const string_t& fname, bool withDomains) const;       //!< save to vtk file format
301     void vtkExport(std::ostream& out) const;                             //!< export mesh as vtk file format to ostream
302     void vtkExport(const GeomDomain& dom, std::ostream& out) const;      //!< export mesh domain as vtk file format
303     void saveToMsh(const string_t& fname, bool withDomains=false) const; //!< save to msh file format
304     void mshExport(std::ostream& out) const;                             //!< export mesh as msh file format to ostream
305     void mshExport(const GeomDomain& dom, std::ostream& out) const;      //!< export mesh domain as msh file format
306     void saveToMel(const string_t& fname, bool withDomains) const;       //!< save to mel file format
307     void melExport(std::ostream& out) const;                             //!< export mesh as mel file format to ostream
308     void exportNodes(const string_t&) const;                             //!< export nodes to a file
309 
310     //================================================
311     //         transformations facilities
312     //================================================
313     //! apply a geometrical transformation on a Mesh
314     Mesh& transform(const Transformation& t);
315     //! apply a translation on a Mesh (vector version)
316     Mesh& translate(std::vector<real_t> u = std::vector<real_t>(3,0.));
317     //! apply a translation on a Mesh (3 reals version)
318     Mesh& translate(real_t ux, real_t uy = 0., real_t uz = 0.);
319     //! apply a rotation on a Mesh
320     Mesh& rotate2d(const Point& c = Point(0.,0.), real_t angle = 0.);
321     //! apply a rotation on a Mesh
322     Mesh& rotate3d(const Point& c = Point(0.,0.,0.), std::vector<real_t> u = std::vector<real_t>(3,0.), real_t angle = 0.);
323     //! apply a rotation on a Mesh
324     Mesh& rotate3d(real_t ux, real_t uy, real_t angle);
325     //! apply a rotation on a Mesh
326     Mesh& rotate3d(real_t ux, real_t uy, real_t uz, real_t angle);
327     //! apply a rotation on a Mesh
328     Mesh& rotate3d(const Point& c, real_t ux, real_t uy, real_t angle);
329     //! apply a rotation on a Mesh
330     Mesh& rotate3d(const Point& c, real_t ux, real_t uy, real_t uz, real_t angle);
331     //! apply a homothety on a Mesh
332     Mesh& homothetize(const Point& c = Point(0.,0.,0.), real_t factor = 1.);
333     //! apply a homothety on a Mesh
334     Mesh& homothetize(real_t factor);
335     //! apply a point reflection on a Mesh
336     Mesh& pointReflect(const Point& c = Point(0.,0.,0.));
337     //! apply a reflection2d on a Mesh
338     Mesh& reflect2d(const Point& c = Point(0.,0.), std::vector<real_t> u = std::vector<real_t>(2,0.));
339     //! apply a reflection2d on a Mesh
340     Mesh& reflect2d(const Point& c, real_t ux, real_t uy = 0.);
341     //! apply a reflection3d on a Mesh
342     Mesh& reflect3d(const Point& c = Point(0.,0.,0.), std::vector<real_t> u = std::vector<real_t>(3,0.));
343     //! apply a reflection3d on a Mesh
344     Mesh& reflect3d(const Point& c, real_t ux, real_t uy, real_t uz = 0.);
345 
346   private:
347     //================================================
348     // build functions using a structured algorithm
349     //================================================
350 
351     //! Pk regular mesh of a segment
352     void meshP1Segment(Segment&, number_t, const std::vector<string_t>&);
353     //! Pk regular mesh of a parallelogram
354     void meshP1Parallelogram(Parallelogram&, number_t, number_t, const std::vector<string_t>&);
355     //! Qk regular mesh of a parallelogram
356     void meshQ1Parallelogram(Parallelogram&, number_t, number_t, const std::vector<string_t>&);
357     //! P1 regular mesh of a parallelepiped
358     void meshP1Parallelepiped(Parallelepiped& , number_t, number_t, number_t, const std::vector<string_t>&);
359     //! Q1 regular mesh of a parallelepiped
360     void meshQ1Parallelepiped(Parallelepiped& , number_t, number_t, number_t, const std::vector<string_t>&);
361     //! Prisme regular mesh of a parallelepiped
362     void meshPr1Parallelepiped(Parallelepiped& , number_t, number_t, number_t, const std::vector<string_t>&);
363     //! Pyramid regular mesh of a parallelepiped
364     void meshPy1Parallelepiped(Parallelepiped& , number_t, number_t, number_t, const std::vector<string_t>&);
365 
366     //! P1 mesh of a polygonal line
367     void meshP1Line(SetOfPoints& sp);
368     //! Utility function used to build 1D mesh
369     void complete1Dmesh(const string_t& domname, const std::vector<string_t>& sideNames);
370 
371     //================================================
372     // build functions using a subdivision algorithm
373     //================================================
374 
375     //! Pk mesh of a (part of a) sphere (surface or volume)
376     void subdvMesh(Ball& sph, const ShapeType shape, const int nboctants,
377                    const number_t nbsubdiv=0, const number_t order=1,
378                    const number_t type=1, const string_t& TeXFilename = "");
379 
380     //! Pk mesh of a (part of a) cube (volume)
381     void subdvMesh(Cube& cube, const ShapeType shape, const int nboctants,
382                    const number_t nbsubdiv=0, const number_t order=1,
383                    const string_t& TeXFilename = "");
384 
385     //! Pk mesh of a truncated cone, a cone or a cylinder (volume)
386     void subdvMesh(RevTrunk& cyl, const ShapeType shape, const number_t nbSubDomains=1,
387                    const number_t nbsubdiv=0, const number_t order=1,
388                    const number_t type=1, const string_t& TeXFilename = "");
389 
390     //! Pk or Qk mesh of a portion of disk (surface)
391     void subdvMesh(Disk& carc, const ShapeType shape,
392                    const number_t nbsubdiv, const number_t order,
393                    const number_t type, const string_t& TeXFilename = "");
394 
395     //! Pk or Qk mesh starting from an initial mesh (surface or volume)
396     void subdvMesh(const std::vector<Point>& pts, const std::vector<std::vector<number_t> >& elems,
397                    const std::vector<std::vector<number_t> >& bounds, const ShapeType shape,
398                    const number_t nbsubdiv=0, const number_t order=1, const string_t& TeXFilename = "");
399 
400   protected:
401     //! Construction of a mesh of order k in the volume of a 3D domain from a mesh built with the help of a subdivision algorithm
402     template <class ST_> void build3Delements(subdivision::GeomFigureMesh<ST_> *TM_p, const ShapeType elShape, number_t nbsubdom = 1);
403     //! Construction of a mesh of order k in 2D or over a surface in R3 from a mesh built with the help of a subdivision algorithm
404     template <class ST_> void build2Delements(subdivision::GeomFigureMesh<ST_> *TM_p, const ShapeType elShape, number_t nbsubdom = 1);
405     //! Construction of an extruded mesh from section mesh
406     void buildExtrusion(const Mesh& ms, const Point& O, const Point& d , number_t nbl,
407                         number_t namingDomain=0, number_t namingSection=0, number_t namingSide=0);
408     //! Construction of an pyramid mesh from hexaedron mesh
409     void buildPyramidFromHexadron(const Mesh& hexMesh);
410     //! Compute the bounding box containing the mesh
411     const std::vector<RealPair> computeBB();
412 
413   private:
414     //! Utilitary function to copy information from subdivision mesh to XLiFE++
415     template <class ST_> void copyPtsEltsDoms(subdivision::GeomFigureMesh<ST_> *TM_p, const ShapeType elShape,
416                                               const dimen_t elementDim, const dimen_t spaceDim, number_t& el_k,
417                                               const std::vector<string_t>& subdomainNames, const string_t& domName);
418     //! Manages construction of domains of a mesh computed using a subdivision algorithm.
419     void manageDomains(subdivision::SubdivisionMesh *SM_p, number_t nbsubdom, dimen_t elementDim,
420                        std::vector<string_t>& boundaryNames, std::vector<string_t>& interfaceNames,
421                        std::vector<string_t>& subdomainNames, string_t& domName);
422 
423     //================================================
424     //           input/output facilities
425     //================================================
426     //! Read mesh from input file filename according to Gmsh mesh format version 2.2 at least
427     void loadGmsh(const string_t& filename, number_t nodesDim=0);
428     //! Read mesh from input file filename according to Gmsh geometry file format (calls and runs gmsh)
429     void loadGeo(const string_t& filename, number_t nodesDim=0);
430     //! Read mesh from input file filename according to Medit Format
431     void loadMedit(const string_t& filename, number_t nodesDim=0);
432     //! Read mesh from input file filename according to Melina mesh format
433     void loadMelina(const string_t& filename, number_t nodesDim=0);
434     //! Read mesh from input file filename according to VTK mesh format
435     void loadVtk(const string_t& filename, number_t nodesDim=0);
436     //! loads mesh from input stream according to POLYDATA dataset in vtk mesh format
437     dimen_t loadPolyDataVtk(std::istream& data, number_t nodesDim=0);
438     //! loads mesh from input stream according to RECTILINEAR_GRID dataset in vtk mesh format
439     dimen_t loadRectilinearVtk(std::istream& data, number_t nodesDim=0);
440     //! loads mesh from input stream according to STRUCTURED_GRID dataset in vtk mesh format
441     dimen_t loadStructuredGridVtk(std::istream& data, number_t nodesDim=0);
442     //! loads mesh from input stream according to UNSTRUCTURED_GRID dataset in vtk mesh format
443     dimen_t loadUnstructuredGridVtk(std::istream& data, number_t nodesDim=0);
444     //! Read mesh from input file filename according to VTP mesh format
445     void loadVtp(const string_t& filename, number_t nodesDim=0);
446     //! Read mesh from input file filename according to VTR mesh format
447     void loadVtr(const string_t& filename, number_t nodesDim=0);
448     //! Read mesh from input file filename according to VTS mesh format
449     void loadVts(const string_t& filename, number_t nodesDim=0);
450     //! Read mesh from input file filename according to VTU mesh format
451     void loadVtu(const string_t& filename, number_t nodesDim=0);
452     //! Read mesh from input file filename according to PLY mesh format
453     void loadPly(const string_t& filename, number_t nodesDim=0);
454 
455     //! writing a 1D geometry in a geo file
456     void saveToGeo(Geometry& g, number_t order, const string_t& filename);
457     //! writing a 2D/3D geometry in a geo file
458     void saveToGeo(Geometry& g, ShapeType sh, number_t order, std::set<MeshOption> mos, const string_t& filename);
459 
460   public:
461     //! adaptative domains for levelset methods
462     void adaptDomains(number_t mainDomNum, const std::vector<number_t>& elementsColorMap);
463 };
464 
465 std::ostream& operator<<(std::ostream&, const Mesh&); //!< print operator
466 
467 //! set the characteristic length
setDefaultCharacteristicLength(real_t h)468 inline void setDefaultCharacteristicLength(real_t h) { theDefaultCharacteristicLength=h; }
469 
470 //! merge two meshes in a new one
471 Mesh merge(const Mesh&, const Mesh&);
472 
473 //! returns the id msh number of a shape, according to the finite elements order
474 number_t mshType(ShapeType sht, number_t order);
475 
476 //@{
477 //! input/output function
saveToFile(const string_t & fname,const Mesh & m,IOFormat iof=_undefFormat,bool withDomains=false)478 inline void saveToFile(const string_t& fname, const Mesh& m, IOFormat iof=_undefFormat, bool withDomains=false)
479 { m.saveToFile(fname, iof, withDomains); }
480 void saveToFile(const string_t&, const GeomDomain&, IOFormat iof=_undefFormat);
481 //@}
482 
483 //! export a split mesh of a domain to msh format
484 void mshExport(const GeomDomain& dom, const std::vector<Point>& coords,
485                const splitvec_t& elementsInfo, std::ostream& out);
486 //! export a split mesh of a domain to Matlab - Octave format
487 void mtlbExport(const GeomDomain& dom, const std::vector<Point>& coords,
488                 const splitvec_t& elementsInfo, std::ostream& out);
489 //! export a split mesh of a domain to vtk format
490 void vtkExport(const GeomDomain& dom, const std::vector<Point>& coords,
491                const splitvec_t& elementsInfo, std::ostream& out);
492 
493 //================================================
494 //         transformations facilities
495 //================================================
496 //! apply a geometrical transformation on a Mesh (external)
497 Mesh transform(const Mesh& m, const Transformation& t);
498 //! apply a translation on a Mesh (vector version) (external)
499 Mesh translate(const Mesh& m, std::vector<real_t> u = std::vector<real_t>(3,0.));
500 //! apply a translation on a Mesh (3 reals version) (external)
501 Mesh translate(const Mesh& m, real_t ux, real_t uy = 0., real_t uz = 0.);
502 //! apply a rotation on a Mesh (external)
503 Mesh rotate2d(const Mesh& m, const Point& c = Point(0.,0.), real_t angle = 0.);
504 //! apply a rotation on a Mesh (external)
505 Mesh rotate3d(const Mesh& m, const Point& c = Point(0.,0.,0.), std::vector<real_t> u = std::vector<real_t>(3,0.), real_t angle = 0.);
506 //! apply a rotation on a Mesh (external)
507 Mesh rotate3d(const Mesh& m, real_t ux, real_t uy, real_t angle);
508 //! apply a rotation on a Mesh (external)
509 Mesh rotate3d(const Mesh& m, real_t ux, real_t uy, real_t uz, real_t angle);
510 //! apply a rotation on a Mesh (external)
511 Mesh rotate3d(const Mesh& m, const Point& c, real_t ux, real_t uy, real_t angle);
512 //! apply a rotation on a Mesh (external)
513 Mesh rotate3d(const Mesh& m, const Point& c, real_t ux, real_t uy, real_t uz, real_t angle);
514 //! apply a homothety on a Mesh (external)
515 Mesh homothetize(const Mesh& m, const Point& c = Point(0.,0.,0.), real_t factor = 1.);
516 //! apply a homothety on a Mesh (external)
517 Mesh homothetize(const Mesh& m, real_t factor);
518 //! apply a point reflection on a Mesh (external)
519 Mesh pointReflect(const Mesh& m, const Point& c = Point(0.,0.,0.));
520 //! apply a reflection2d on a Mesh (external)
521 Mesh reflect2d(const Mesh& m, const Point& c = Point(0.,0.), std::vector<real_t> u = std::vector<real_t>(2,0.));
522 //! apply a reflection2d on a Mesh (external)
523 Mesh reflect2d(const Mesh& m, const Point& c, real_t ux, real_t uy = 0.);
524 //! apply a reflection3d on a Mesh (external)
525 Mesh reflect3d(const Mesh& m, const Point& c = Point(0.,0.,0.), std::vector<real_t> u = std::vector<real_t>(3,0.));
526 //! apply a reflection3d on a Mesh (external)
527 Mesh reflect3d(const Mesh& m, const Point& c, real_t ux, real_t uy, real_t uz = 0.);
528 
529 } // end of namespace xlifepp
530 
531 #endif // MESH_HPP
532