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\section{Defining geometries}
18
19Whatever the geometry is, the common characteristics are a dimension, a shape, a name, a bounding box and a more convenient box for geometry inclusion than the previous one, and easier to define than the convex hull, the so-called "minimal box". This is the reason why we define 3 classes : {\class BoundingBox}, {\class MinimalBox} and {\class Geometry}.
20
21\subsection{{\classtitle BoundingBox} and {\classtitle MinimalBox} classes}
22
23The bounding box is the smallest box whose faces are parallel to axis containing the geometry. As a segment in 1D, a rectangle in 2D or a parallelepiped in 3D, it can be defined from 2, 3 or 4 points, as illustrated in the following figure :
24
25\begin{figure}[H]
26\begin{center}
27\inputTikZ{box}
28\end{center}
29\caption{Definition of a box from points}
30\label{fig.boxdef}
31\end{figure}
32
33For example, if we define a segment whose boundaries are the points $A(0,0,-1)$ and $B(1,2,3)$, the bounding box is $[0;1]\times[0;2]\times[-1;3]$. The {\class BoundingBox} class is then defined as a vector of pairs of real numbers :
34
35\begin{lstlisting}
36typedef std::pair<Real, Real> RealPair;
37class BoundingBox
38{
39  private :
40    std::vector<RealPair> bounds_; //!< Bounding values of the box [xmin,xmax]x[ymin,ymax]x[zmin,zmax]
41\end{lstlisting}
42
43It offers a various list of constructors, from set of real pairs or set of {\class Point} :
44
45\begin{lstlisting}[deletekeywords={[3] vp}]
46BoundingBox() {}
47BoundingBox(const std::vector<RealPair>&);
48BoundingBox(Real, Real, Real, Real, Real, Real);
49BoundingBox(Real, Real, Real, Real);
50BoundingBox(Real, Real);
51BoundingBox(const Point&, const Point&, const Point&, const Point&);
52BoundingBox(const Point&, const Point&, const Point&);
53BoundingBox(const Point&, const Point&);
54BoundingBox(const std::vector<Point>& vp);
55\end{lstlisting}
56
57There is also :
58
59\begin{itemize}
60\item a various list of accessors :
61\begin{lstlisting}
62std::vector<RealPair> bounds() const; // return bounds in all direction
63RealPair bounds(Dimen i) const;   // return bounds in direction i (1 to dim)
64RealPair xbounds() const;         // return bounds in first direction
65RealPair ybounds() const;         // return bounds in second direction
66RealPair zbounds() const;         // return bounds in third direction
67Point minPoint() const;           // return the min point (left,bottom,front)
68Point maxPoint() const;           // return the max point (right,top,back)
69\end{lstlisting}
70\item print facilities :
71\begin{lstlisting}
72void print(std::ostream&) const; // print BoundingBox
73String asString() const;         // format as string : [a,b]x[c,d]x[e,f]
74
75std::ostream& operator<<(std::ostream&, const BoundingBox&); //output box
76std::ostream& operator<<(std::ostream&, const RealPair&); //output pair of reals
77\end{lstlisting}
78
79\item some utility functions, such as merging 2 bounding boxes or giving its dimension :
80\begin{lstlisting}
81Dimen dim() const; // dimension of the bounding box
82BoundingBox& operator +=(const BoundingBox&); // merge a bounding box
83\end{lstlisting}
84\end{itemize}
85
86As the segment example shows, the dimension of the geometry and the dimension of the bounding box are independent. So if we want to check if a geometry is inside another geometry, this is not the right box to use, if we want to be independent of the true shape of the geometry. That is why the {\class MinimalBox} class is defined :
87
88\begin{lstlisting}
89class MinimalBox
90{
91  protected :
92    std::vector<Point> bounds_;
93\end{lstlisting}
94You can notice that it is {\class Point} that are stored, and only the 2, 3 or 4 points necessary to define the box, as shown in \autoref{fig.boxdef}
95
96The {\class MinimalBox} class offers quite the same constructors as the {\class BoundingBox} class :
97\begin{lstlisting}
98MinimalBox() {}
99MinimalBox(const std::vector<Point>&);
100MinimalBox(const std::vector<RealPair>&);
101MinimalBox(Real, Real, Real, Real, Real, Real);
102MinimalBox(Real, Real, Real, Real);
103MinimalBox(Real, Real);
104MinimalBox(const Point&, const Point&, const Point&, const Point&);
105MinimalBox(const Point&, const Point&, const Point&);
106MinimalBox(const Point&, const Point&);
107\end{lstlisting}
108
109It also offers :
110\begin{itemize}
111\item some accessors :
112\begin{lstlisting}
113std::vector<Point> bounds() const;
114Point boundPt(Dimen i) const; // return bounds in direction i (1 to dim)
115Point minPoint() const;
116Point maxPoint() const;
117\end{lstlisting}
118\item print facilities :
119\begin{lstlisting}
120void print(std::ostream&) const;
121String asString() const;
122
123
124std::ostream& operator<<(std::ostream&, const MinimalBox&); // output MinimalBox
125\end{lstlisting}
126\item some utility functions, the same as in {\class BoundingBox}, plus a function to return all vertices of the box, not only the points stored in {\var bounds\_}.
127\begin{lstlisting}
128Dimen dim() const;
129MinimalBox& operator +=(const MinimalBox&);
130std::vector<Point> vertices() const;
131\end{lstlisting}
132\end{itemize}
133
134\begin{remark}
135A minimal box of a segment will always be a segment.
136A minimal box of a 2D plane geometry will always be a rectangle.
137\end{remark}
138
139\medskip
140
141Now that both boxes classes are clarified, we can now look at the main class : {\class Geometry}.
142
143\subsection{The {\classtitle Geometry} class}
144
145A geometry can be of 3 natures : canonical geometries, such as {\class Rectangle} or {\class Ball}, so-called "loop" geometries defined by its set of geometrical boundaries, and so-called "composite" geometries defined by union of elementary geometries and holes. The definition of the {\class Geometry} class takes into account these different natures :
146
147\begin{lstlisting}[deletekeywords={[3] name, map}]
148enum ShapeType
149{
150  _noshape=0,
151  _point ,
152  _segment, segment=_segment,
153  _triangle, triangle=_triangle,
154  _quadrangle, quadrangle=_quadrangle,
155  _tetrahedron, tetrahedron=_tetrahedron,
156  _hexahedron, hexahedron=_hexahedron,
157  _prism, prism=_prism,
158  _pyramid, pyramid=_pyramid,
159  _ellArc, _circArc,
160  _polygon, _parallelogram, _rectangle, _square, _ellipse, _disk, _ellipsoidPart, _spherePart,
161  _setofelems, _trunkPart, _cylinderPart, _conePart,
162  _polyhedron, _parallelepiped, _cuboid, _cube, _ellipsoid, _ball, _trunk,
163  _revTrunk, _cylinder, _revCylinder, _cone, _revCone,
164  _composite, _loop, _extrusion
165};
166
167class Geometry
168{
169  public:
170    BoundingBox boundingBox;
171    MinimalBox minimalBox;
172    mutable bool force; // tag to force inclusion when geometry engines fails
173    mutable bool crackable; // tag to define if the geometry should be cracked or not
174  protected:
175    string_t domName_;   // geometry name
176    Dimen dim_;       // physical space dimension
177    ShapeType shape_; // geometrical shape
178    std::vector<string_t> sideNames_; //names of sides (boundaries)
179    std::vector<string_t> theNamesOfVariables_; // variable names (x,y,z),...
180    string_t teXFilename_;                      // name of the file containing the tex code drawing the geometry
181  private:
182    mutable std::map<number_t, Geometry*> components_; // for composite/loop geometry
183    std::map<number_t, std::vector<number_t> > geometries_; // for composite geometry
184    std::map<number_t, std::vector<number_t> > loops_; // for loop geometry
185    mutable Transformation* extrusion_; // for extrusion geometry
186    number_t layers_; // for extrusion geometry
187\end{lstlisting}
188\vspace{0.2cm}
189
190{\var name} and {\var sideNames\_} are devoted to store domain names or side domain names.
191
192\begin{figure}[H]
193\includeTikZPict[width=18cm]{geometries.png}
194\caption{{\lib Geometry} child classes}
195\label{diag_geometries}
196\end{figure}
197
198The {\class Geometry} is the base class for a large set of canonical geometries, so the management of "composite" and "loop" geometries implies 2 main vectors of {\class Geometry} pointers and a set of virtual functions to cast a {\class Geometry} into one of its child :
199
200\begin{lstlisting}[deletekeywords={[3] name}]
201//! access to child Xxxxxxx object (const)
202virtual const Xxxxxxx* xxxxxxx() const { error("bad_geometry", name, words("shape",shape_), words("shape",_yyyyyyy)); return 0; }
203//! access to child Xxxxxxx object
204virtual Xxxxxxx* xxxxxxx() { error("bad_geometry", name, words("shape",shape_), words("shape",_yyyyyyy)); return 0; }
205\end{lstlisting}
206
207{\tt xxxxxxx} must be one of the following keyword : {\tt segment, ellArc, circArc, polygon, triangle, quadrangle, parallelogram, rectangle, square, ellipse, disk, polyhedron, tetrahedron, hexahedron, parallelepiped, cuboid, cube, ellipsoid, ball, revTrunk, trunk, cylinder, revCylinder, prism, cone, revCone, pyramid, setOfElems }
208
209{\tt yyyyyyy} must be the corresponding {\class ShapeType} item.
210
211\medskip
212
213The {\class Geometry} class offers constructors from {\class BoundingBox}, dimension or both :
214
215\begin{lstlisting}[deletekeywords={[3] name, nx, ny, nz}]
216Geometry(): force(false), crackable(false), domName_(""), dim_(0), shape_(_noshape), extrusion_(0), layers_(0) { sideNames_.resize(0); }  // void constructor
217Geometry(const Geometry& g);                        // copy constructor
218// basic constructor from bounding box
219Geometry(const BoundingBox&, const string_t& na = "?", ShapeType sht=_noshape,
220     const string_t& nx = "x", const string_t& ny = "y", const string_t& nz = "z");
221Geometry(const BoundingBox&, dimen_t, const string_t& na = "?",
222     ShapeType sht=_noshape, const string_t& nx = "x", const string_t& ny = "y",
223     const string_t& nz = "z");
224Geometry(dimen_t, const string_t& na = "?", ShapeType sht=_noshape,
225     const string_t& nx = "x", const string_t& ny = "y", const string_t& nz = "z");
226virtual ~Geometry()  {}
227\end{lstlisting}
228
229It also offers :
230\begin{itemize}
231\item some accessors :
232\begin{lstlisting}[deletekeywords={[3] map}]
233dimen_t dim() const;
234ShapeType shape() const;
235void shape(ShapeType sh);
236string_t domName() const;
237void domName(const string_t& nm);
238string_t teXFilename() const;
239void teXFilename(const string_t& fn);
240const std::vector<string_t> sideNames() const;
241std::vector<string_t> sideNames();
242std::map<number_t, Geometry*> components();
243const std::map<number_t, Geometry*> components() const;
244const std::map<number_t, std::vector<number_t> > geometries() const;
245std::map<number_t, std::vector<number_t> > loops() const;
246number_t layers() const;
247Transformation* extrusion();
248const Transformation* extrusion() const;
249\end{lstlisting}
250\item some print facilities :
251\begin{lstlisting}
252void print(std::ostream&) const;           //!< print Geometry
253friend std::ostream& operator<<(std::ostream&, const Geometry&); //!< output Geometry
254\end{lstlisting}
255\item some functions to compute or update boxes :
256\begin{lstlisting}
257virtual void computeMB(); //! compute the minimal box
258void updateBB(const std::vector<RealPair>& bb) {boundingBox = BoundingBox(bb);}
259\end{lstlisting}
260As the minimal box depends on the true shape of the geometry, it is a virtual function reimplemented in child classes.
261\item some virtual functions to get data from canonical geometries, such as the vector of points defining the geometry, or the parameters for the subdivision mesh algorithm :
262\begin{lstlisting}
263virtual std::vector<const Point*> boundNodes() const;
264virtual std::vector<Point*> nodes();
265virtual std::vector<const Point*> nodes() const;
266virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
267virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
268\end{lstlisting}
269\item member functions to determine if a geometry is coplanar to another geometry or inside another geometry, and external functions or operators to define "composite", "loop" or "extrusions" geometries.
270\begin{lstlisting}[deletekeywords={[3] name}]
271bool isCoplanar(const Geometry& g) const;
272bool isInside(const Geometry&) const; //check if geometry is inside an other one
273
274friend Geometry operator+(const Geometry&, const Geometry&); //union
275Geometry& operator+=(const Geometry&);
276friend Geometry operator-(const Geometry&, const Geometry&); //difference (hole)
277Geometry& operator-=(const Geometry&);
278
279// definition of a geometry 2D/3D from its boundary 1D/2D
280friend Geometry surfaceFrom(const Geometry&, String name = "");
281friend Geometry volumeFrom(const Geometry&, String name = "");
282friend Geometry extrude(const Geometry& g, const Transformation& t, number_t layers, string_t domName, std::vector<string_t> sidenames);
283friend Geometry extrude(const Geometry& g, const Transformation& t, number_t layers, std::vector<string_t> sidenames);
284\end{lstlisting}
285
286Operators + and - work as follows :
287\begin{itemize}
288\item if both arguments are canonical geometries (neither "composite" nor "loop"), it puts pointers to them in the {\var components\_} attribute and if one geometry is inside the other, it puts it in the {\var geometries} attribute.
289\item If at least one of the argument is "composite" or "loop", it inserts the canonical geometry correctly, or merge the 4 maps.
290\end{itemize}
291
292{\cmd surfaceFrom} and {\cmd volumeFrom} allow to define surfaces and respectively volumes from their boundaries, so-called "loop" geometries. The mechanism is the same : taking a "composite" geometry as argument, it puts every component in the {\var borders\_} attribute and defines the loop in the {\var loops\_} attribute, so that we can define "composite" geometries of "loop" geometries, but it is for the moment forbidden.
293\end{itemize}
294
295\subsection{{\classtitle Curve}, {\classtitle Surface} and {\classtitle Volume} classes}
296
297
298The {\class Curve} class is the base class for every canonical 1D geometry. It has no attributes and only general constructors :
299
300\begin{lstlisting}[deletekeywords={[3] name}]
301class Curve : public Geometry
302{
303  public:
304    Curve();
305  protected:
306    void buildParam(const Parameter& p);
307    void buildDefaultParam(ParameterKey key);
308    std::set<ParameterKey> getParamsKeys();
309  public:
310    Curve(const Curve& c) : Geometry(c) {} //!< copy constructor
311    virtual Geometry* clone() const { return new Curve(*this); } //!< virtual copy constructor for Geometry
312    virtual ~Curve() {} //!< virtual destructor
313};
314\end{lstlisting}
315
316The {\class Surface} class is the base class for every canonical 2D geometry. Actually, it has no attributes and only constructors :
317
318\begin{lstlisting}[deletekeywords={[3] name}]
319class Surface : public Geometry
320{
321  public:
322    //! default constructor for 2D geometries
323    Surface();
324
325  protected:
326    void buildParam(const Parameter& p);
327    void buildDefaultParam(ParameterKey key);
328    std::set<ParameterKey> getParamsKeys();
329
330  public:
331    Surface(const Surface& s) : Geometry(s) {} //!< copy constructor
332    virtual Geometry* clone() const { return new Surface(*this); } //!< virtual copy constructor for Geometry
333    virtual Surface* cloneS() const { return new Surface(*this); } //!< virtual copy constructor for Surface
334    virtual ~Surface() {} //!< virtual destructor
335
336    virtual std::vector<Point> p() const
337    { error("not_yet_implemented","std::vector<Point> Surface::p() const"); return std::vector<Point>(); }
338    virtual std::vector<number_t> n() const
339    { error("not_yet_implemented","std::vector<Point> Surface::n() const"); return std::vector<number_t>(); }
340    virtual Point p(number_t i) const
341    { error("not_yet_implemented","Point Surface::p(Number i) const"); return Point(); }
342    virtual number_t n(number_t i) const
343    { error("not_yet_implemented","Number Surface::n(Number i) const"); return 0; }
344    virtual number_t& n(number_t i)
345    { error("not_yet_implemented","Number& Surface::n(Number i)"); return *new number_t(0); }
346    virtual std::vector<real_t> h() const
347    { error("not_yet_implemented","std::vector<Real> Surface::h() const"); return std::vector<real_t>(); }
348    virtual real_t h(number_t i) const
349    { error("not_yet_implemented","Number Surface::h(Number i) const"); return 0.; }
350    virtual real_t& h(number_t i)
351    { error("not_yet_implemented","Number& Surface::h(Number i)"); return *new real_t(0); }
352
353    bool isPolygon()
354    {
355      if (shape_ == _polygon || shape_ == _triangle || shape_ == _quadrangle || shape_ == _parallelogram || shape_ == _rectangle ||
356          shape_ == _square) { return true; }
357      return false;
358    }
359
360    //! format as string
361    virtual string_t asString() const { error("shape_not_handled", words("shape",shape_)); return string_t(); }
362};
363\end{lstlisting}
364
365The {\class Volume} class is the base class for every canonical 3D geometry. Actually, it has no attributes and only constructors :
366
367\begin{lstlisting}[deletekeywords={[3] name}]
368class Volume : public Geometry
369{
370  public:
371    //! default constructor for 3D geometries
372    Volume();
373
374  protected:
375    void buildParam(const Parameter& p);
376    void buildDefaultParam(ParameterKey key);
377    std::set<ParameterKey> getParamsKeys();
378
379  public:
380    Volume(const Volume& v) : Geometry(v) {} //!< copy constructor
381    virtual Geometry* clone() const { return new Volume(*this); } //!< virtual copy constructor for Geometry
382    virtual ~Volume() {} //!< virtual destructor
383
384    //! format as string
385    virtual string_t asString() const { error("shape_not_handled", words("shape",shape_)); return string_t(); }
386};
387\end{lstlisting}
388
389\subsection{Definition of canonical geometries}
390
391In this subsection, we will explain how to define canonical geometries. All constructors for these classes are writtent in the same way : with {\class Parameter} arguments, so that the user can use a key-value system.
392
393\subsubsection{The key-value system}
394
395Some keys are easy to determine : the center for rectangles, squares, ellipses, disks, cuboids, cubes, ellipsoids and balls, the radius for disks, balls, and cylinders, \ldots. To do so, global keys are defined (type {\class Parameter}). Here is the list of available keys with authorized data types and classes using them:
396
397\medskip
398
399\begin{tabular}{|p{3cm}|p{7cm}|p{6cm}|}
400\hline
401key & authorized types & concerned classes \\
402\hline
403\hline
404\param \_angle1 & single integer or real value & class {\class Disk} \\
405\hline
406\param \_angle2 & single integer or real value & class {\class Disk} \\
407\hline
408\param \_apex & single integer or real value, or {\class Point} & classes {\class Cone}, {\class Pyramid}, {\class RevCone} \\
409\hline
410\param \_apogee & single integer or real value, or {\class Point} & class {\class EllArc} \\
411\hline
412\hline
413\param \_basis & classes {\class Polygon}, {\class Triangle}, {\class Quadrangle}, {\class Parallelogram}, {\class Rectangle}, {\class Square}, {\class Ellipse}, {\class Disk} & classes {\class Trunk}, {\class Cylinder}, {\class Prism}, {\class Cone}, {\class Pyramid} \\
414\hline
415\hline
416\param \_center & single integer or real value, or {\class Point} & classes {\class EllArc}, {\class CircArc}, {\class Rectangle}, {\class Square}, {\class Ellipse}, {\class Disk}, {\class Cuboid}, {\class Cube}, {\class Ellipsoid}, {\class Ball}, {\class RevCone} \\
417\hline
418\param \_center1 & single integer or real value, or {\class Point} & classes {\class Trunk}, {\class Cylinder}, {\class Cone}, {\class RevTrunk}, {\class RevCylinder} \\
419\hline
420\param \_center2 & single integer or real value, or {\class Point} & classes {\class Trunk}, {\class Cylinder}, {\class RevTrunk}, {\class RevCylinder} \\
421\hline
422\hline
423\param \_direction & std::vector of real values or {\class Point} & classes {\class Cylinder}, {\class Prism} \\
424\hline
425\param \_domain\_name & single string & every class \\
426\hline
427\hline
428\param \_end\_shape & enum {\class GeometricEndShape} & class {\class RevCone} \\
429\hline
430\param \_end1\_shape & enum {\class GeometricEndShape} & classes{\class RevTrunk}, {\class RecCylinder} \\
431\hline
432\param \_end2\_shape & enum {\class GeometricEndShape} & classes {\class RevTrunk}, {\class RecCylinder} \\
433\hline
434\param \_end\_distance & single unsigned  integer or real positive value & class {\class RevCone} \\
435\hline
436\param \_end1\_distance & single unsigned  integer or real positive value & classes {\class RevTrunk}, {\class RecCylinder} \\
437\hline
438\param \_end2\_distance & single unsigned  integer or real positive value & classes {\class RevTrunk}, {\class RecCylinder} \\
439\hline
440\hline
441\param \_faces & vector of {\class Polygon} & class {\class Polyhedron} \\
442\hline
443\hline
444\param \_hsteps & single real value, std::vector of real values or {\class Reals} & every class but {\class Polyhedron}  \\
445\hline
446\hline
447\param \_length & single unsigned  integer or real positive value & classes {\class Square}, {\class Cube} \\
448\hline
449\hline
450\param \_nboctants & single unsigned integer value & classes {\class Cube}, {\class Ellipsoid}, {\class Ball} \\
451\hline
452\param \_nbsubdomains & single unsigned integer value & classes {\class RevTrunk}, {\class RevCylinder}, {\class RevCone} \\
453\hline
454\param \_nnodes & single unsigned integer value, std::vector of integer values  or {\class Numbers} & every class but {\class Polyhedron} \\
455\hline
456\hline
457\param \_origin & single integer or real value, or {\class Point} & classes {\class Rectangle}, {\class Square}, {\class Cuboid}, {\class Cube}, {\class Trunk}, {\class Cylinder}, {\class Prism} \\
458\hline
459\end{tabular}
460
461\begin{tabular}{|p{3cm}|p{7cm}|p{6cm}|}
462\hline
463key & authorized types & concerned classes \\
464\hline
465\hline
466\param \_radius & single unsigned  integer or real positive value & classes {\class Disk}, {\class Ball}, {\class RevCylinder}, {\class RevCone} \\
467\hline
468\param \_radius1 & single unsigned  integer or real positive value & class {\class RevTrunk} \\
469\hline
470\param \_radius2 & single unsigned  integer or real positive value & class {\class RevTrunk} \\
471\hline
472\hline
473\param \_scale & single unsigned  integer or real positive value & class {\class Trunk} \\
474\hline
475\param \_side\_names & single string or std::vector of strings & every class but {\class Polyhedron} \\
476\hline
477\hline
478\param \_type & single unsigned integer value & classes {\class Ellipse}, {\class Disk}, {\class Ellipsoid}, {\class Ball}, {\class RevTrunk}, {\class RevCylinder}, {\class RevCone} \\
479\hline
480\hline
481\param \_varnames & single string or std::vector of strings & every class \\
482\hline
483\param \_vertices & vector of {\class Point} & class {\class Polygon} \\
484\hline
485\param \_v1 & single integer or real value, or {\class Point} & every class, but {\class Polygon}, {\class Polyhedron}, {\class RevTrunk}, {\class RevCylinder}, {\class RevCone} \\
486\hline
487\param \_v2 & single integer or real value, or {\class Point} & every class, but {\class Polygon}, {\class Polyhedron}, {\class RevTrunk}, {\class RevCylinder}, {\class RevCone} \\
488\hline
489\param \_v3 & single integer or real value, or {\class Point} & classes {\class Triangle}, {\class Quadrangle}, {\class Tetrahedron}, {\class Hexahedron}, {\class Prism}, {\class Pyramid} \\
490\hline
491\param \_v4 & single integer or real value, or {\class Point} &  classes {\class Quadrangle}, {\class Parallelogram}, {\class Rectangle}, {\class Square}, {\class Tetrahedron}, {\class Hexahedron}, {\class Parallelepiped}, {\class Cuboid}, {\class Cube}, {\class Pyramid} \\
492\hline
493\param \_v5 & single integer or real value, or {\class Point} & classes {\class Hexahedron}, {\class Parallelepiped}, {\class Cuboid}, {\class Cube} \\
494\hline
495\param \_v6 & single integer or real value, or {\class Point} & class {\class Hexahedron} \\
496\hline
497\param \_v7 & single integer or real value, or {\class Point} & class {\class Hexahedron} \\
498\hline
499\param \_v8 & single integer or real value, or {\class Point} & class {\class Hexahedron} \\
500\hline
501\hline
502\param \_xlength & single unsigned  integer or real positive value & classes {\class Rectangle}, {\class Ellipse}, {\class Cuboid}, {\class Ellipsoid} \\
503\hline
504\param \_xmin & single integer or real value & classes {\class Rectangle}, {\class Cuboid} \\
505\hline
506\param \_xmax & single integer or real value & classes {\class Rectangle}, {\class Cuboid} \\
507\hline
508\hline
509\param \_ylength & single unsigned  integer or real positive value & classes {\class Rectangle}, {\class Ellipse}, {\class Cuboid}, {\class Ellipsoid}\\
510\hline
511\param \_ymin & single integer or real value & classes {\class Rectangle}, {\class Cuboid} \\
512\hline
513\param \_ymax & single integer or real value & classes {\class Rectangle}, {\class Cuboid} \\
514\hline
515\hline
516\param \_zlength & single unsigned  integer or real positive value & class {\class Cuboid} \\
517\hline
518\param \_zmin & single integer or real value & class {\class Cuboid} \\
519\hline
520\param \_zmax & single integer or real value & class {\class Cuboid} \\
521\hline
522\end{tabular}
523
524A constructor of a geometry will take several {\class Parameter} object and transfer the building work to the 3 main following routines:
525
526\begin{lstlisting}
527// main function to build a geometry
528void build(const std::vector<Parameter>& ps);
529// build the list of available keys
530std::set<ParameterKey> getParamsKeys();
531// deal with one key and do the appropriate work
532void buildParam(const Parameter& gp);
533// deal with unused keys having default values and do the appropriate work
534void buildDefaultParam(ParameterKey key);
535\end{lstlisting}
536
537For some geometries, there are some additional routines called by {\cmd build}.
538
539\subsubsection{The {\class Segment} class}
540
541A segment is just a straight line between 2 points, with a number of nodes.
542
543\begin{lstlisting}
544class Segment : public Curve
545{
546  private:
547    Point p1_, p2_;
548    Number n_;
549    std::vector<real_t> h_;
550\end{lstlisting}
551
552\begin{figure}[H]
553\begin{center}
554\inputTikZ{geo-segment}
555\end{center}
556\end{figure}
557
558It offers constructors, taking from 2 to 5 {\class Parameter}. Authorized keys and what they do are:
559
560\begin{center}
561\begin{tabular}{|c|c|c|c|}
562\hline
563attribute & key to use to define it & optional & default value \\
564\hline
565\hline
566p1\_, p2\_ & ({\param \_v1}, {\param \_v2}) & no & no \\
567n\_ & {\param \_nnodes} & yes & 2 \\
568 & or ({\param \_xmin}, {\param \_xmax}) & & \\
569n\_ & {\param \_nnodes} & yes & 2 \\
570h\_ & {\param \_hsteps} & yes & no \\
571domName\_ & {\param \_domain\_name} & yes & "" \\
572sideNames\_ & {\param \_side\_names} & yes & "" \\
573\hline
574\end{tabular}
575\end{center}
576
577It also offers :
578
579\begin{itemize}
580\item some accessors
581\item implementation of inherited virtual functions concerning segments :
582\begin{lstlisting}
583virtual void computeMB();
584virtual std::vector<Point*> nodes();
585virtual std::vector<const Point*> nodes() const;
586virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
587\end{lstlisting}
588\end{itemize}
589
590\subsubsection{The {\class EllArc} class}
591
592To define an elliptic arc, you need 4 points : the center of the ellipse, a point on the main axis (apogee) and the bounds of the arc. When omitted, the apogee point is defined as the first bound of the arc. An elliptic arc must be smalller than a half-ellipse, to be defined correctly. This is a \gmsh restriction.
593
594\begin{lstlisting}
595class EllArc : public Curve
596{
597  protected:
598    Point c_, a_, p1_, p2_;
599    Number n_;
600    std::vector<real_t> h_;
601\end{lstlisting}
602
603\begin{figure}[H]
604\begin{center}
605\inputTikZ{geo-ellipsearc}
606\end{center}
607\end{figure}
608
609It offers constructors, taking from 3 to 7 {\class Parameter}. Authorized keys and what they do are:
610
611\begin{center}
612\begin{tabular}{|c|c|c|c|}
613\hline
614attribute & key to use to define it & optional & default value \\
615\hline
616\hline
617c\_, a\_, p1\_, p2\_ & ({\param \_center}, {\param \_apogee}, {\param \_v1}, {\param \_v2}) & no & no \\
618 & or ({\param \_center}, {\param \_v1}, {\param \_v2}) & & \\
619n\_ & {\param \_nnodes} & yes & 2 \\
620h\_ & {\param \_hsteps} & yes & no \\
621domName\_ & {\param \_domain\_name} & yes & "" \\
622sideNames\_ & {\param \_side\_names} & yes & "" \\
623\hline
624\end{tabular}
625\end{center}
626
627It also offers :
628
629\begin{itemize}
630\item some accessors
631\item implementation of inherited virtual functions concerning elliptic arcs :
632\begin{lstlisting}
633virtual void computeMB();
634virtual std::vector<Point*> nodes();
635virtual std::vector<const Point*> nodes() const;
636virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
637\end{lstlisting}
638\end{itemize}
639
640\subsubsection{The {\class CircArc} class}
641
642To define a circular arc, you need 3 points : the center of the circle and the bounds of the arc. A circular arc must be smalller than a half-circle, to be defined correctly. This is a \gmsh restriction.
643
644\begin{lstlisting}
645class CircArc : public Curve
646{
647  protected:
648    Point c_, p1_, p2_;
649    number_t n_;
650    std::vector<real_t> h_;
651\end{lstlisting}
652
653\begin{figure}[H]
654\begin{center}
655\inputTikZ{geo-circlearc}
656\end{center}
657\end{figure}
658
659It offers constructors, taking from 3 to 6 {\class Parameter}. Authorized keys and what they do are:
660
661\begin{center}
662\begin{tabular}{|c|c|c|c|}
663\hline
664attribute & key to use to define it & optional & default value \\
665\hline
666\hline
667c\_, p1\_, p2\_ & {\param \_center}, {\param \_v1} and {\param \_v2} & no & no \\
668n\_ & {\param \_nnodes} & yes & 2 \\
669h\_ & {\param \_hsteps} & yes & no \\
670domName\_ & {\param \_domain\_name} & yes & "" \\
671sideNames\_ & {\param \_side\_names} & yes & "" \\
672\hline
673\end{tabular}
674\end{center}
675
676It also offers :
677
678\begin{itemize}
679\item some accessors
680\item implementation of inherited virtual functions concerning circular arcs :
681\begin{lstlisting}
682virtual void computeMB();
683virtual std::vector<Point*> nodes();
684virtual std::vector<const Point*> nodes() const;
685virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
686\end{lstlisting}
687\end{itemize}
688
689\subsubsection{The {\class Polygon} class}
690
691To define a polygon, you give the ordered list of vertices.
692
693\begin{lstlisting}
694class Polygon : public Surface
695{
696  protected:
697    std::vector<Point> p_;
698    std::vector<real_t> h_;
699    std::vector<number_t> n_;
700\end{lstlisting}
701
702\begin{figure}[H]
703\begin{center}
704\inputTikZ{geo-polygon}
705\end{center}
706\end{figure}
707
708It offers constructors, taking from 1 to 4 {\class Parameter}. Authorized keys and what they do are:
709
710\begin{center}
711\begin{tabular}{|c|c|c|c|}
712\hline
713attribute & key to use to define it & optional & default value \\
714\hline
715\hline
716p\_ & {\param \_vertices} & no & no\\
717n\_ & {\param \_nnodes} & yes & 2 \\
718h\_ & {\param \_hsteps} & yes & no \\
719domName\_ & {\param \_domain\_name} & yes & "" \\
720sideNames\_ & {\param \_side\_names} & yes & "" \\
721\hline
722\end{tabular}
723\end{center}
724
725It also offers :
726
727\begin{itemize}
728\item some accessors
729\item implementation of inherited virtual functions concerning polygon :
730\begin{lstlisting}
731virtual void computeMB();
732virtual std::vector<const Point*> boundNodes() const;
733virtual std::vector<Point*> nodes();
734virtual std::vector<const Point*> nodes() const;
735virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
736virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
737\end{lstlisting}
738\end{itemize}
739
740\subsubsection{The {\class Triangle} class}
741
742To define a triangle, you give the 3 vertices.
743
744\begin{figure}[H]
745\begin{center}
746\inputTikZ{geo-triangle}
747\end{center}
748\end{figure}
749
750It offers constructors, taking from 3 to 6 {\class Parameter}. Authorized keys and what they do are:
751
752\begin{center}
753\begin{tabular}{|c|c|c|c|}
754\hline
755attribute & key to use to define it & optional & default value \\
756\hline
757\hline
758p\_ & {\param \_v1}, {\param \_v2} and {\param \_v3} & no & no\\
759n\_ & {\param \_nnodes} & yes & 2 \\
760h\_ & {\param \_hsteps} & yes & no \\
761domName\_ & {\param \_domain\_name} & yes & "" \\
762sideNames\_ & {\param \_side\_names} & yes & "" \\
763\hline
764\end{tabular}
765\end{center}
766
767It also offers :
768
769\begin{itemize}
770\item some accessors
771\item implementation of inherited virtual functions concerning triangles :
772\begin{lstlisting}
773virtual void computeMB();
774virtual std::vector<const Point*> boundNodes() const;
775virtual std::vector<Point*> nodes();
776virtual std::vector<const Point*> nodes() const;
777virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
778virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
779\end{lstlisting}
780\end{itemize}
781
782\subsubsection{The {\class Quadrangle} class}
783
784To define a quadrangle, you give the 4 vertices.
785
786\begin{figure}[H]
787\begin{center}
788\inputTikZ{geo-quadrangle}
789\end{center}
790\end{figure}
791
792It offers constructors, taking from 4 to 7 {\class Parameter}. Authorized keys and what they do are:
793
794\begin{center}
795\begin{tabular}{|c|c|c|c|}
796\hline
797attribute & key to use to define it & optional & default value \\
798\hline
799\hline
800p\_ & {\param \_v1}, {\param \_v2}, {\param \_v3} and {\param \_v4} & no & no\\
801n\_ & {\param \_nnodes} & yes & 2 \\
802h\_ & {\param \_hsteps} & yes & no \\
803domName\_ & {\param \_domain\_name} & yes & "" \\
804sideNames\_ & {\param \_side\_names} & yes & "" \\
805\hline
806\end{tabular}
807\end{center}
808
809It also offers :
810
811\begin{itemize}
812\item some accessors
813\item implementation of inherited virtual functions concerning quadrangles :
814\begin{lstlisting}
815virtual void computeMB();
816virtual std::vector<const Point*> boundNodes() const;
817virtual std::vector<Point*> nodes();
818virtual std::vector<const Point*> nodes() const;
819virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
820virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
821\end{lstlisting}
822\end{itemize}
823
824\subsubsection{The {\class Parallelogram} class}
825
826To define a parallelogram, you give 3 vertices. $p_3$ is indeed unnecessary. It can be calculated from the 3 other vertices.
827
828\begin{figure}[H]
829\begin{center}
830\inputTikZ{geo-parallelogram}
831\end{center}
832\end{figure}
833
834It offers constructors, taking from 3 to 6 {\class Parameter}. Authorized keys and what they do are:
835
836\begin{center}
837\begin{tabular}{|c|c|c|c|}
838\hline
839attribute & key to use to define it & optional & default value \\
840\hline
841\hline
842p\_ & {\param \_v1}, {\param \_v2} and {\param \_v4} & no & no\\
843n\_ & {\param \_nnodes} & yes & 2 \\
844h\_ & {\param \_hsteps} & yes & no \\
845domName\_ & {\param \_domain\_name} & yes & "" \\
846sideNames\_ & {\param \_side\_names} & yes & "" \\
847\hline
848\end{tabular}
849\end{center}
850
851It also offers :
852
853\begin{itemize}
854\item some accessors
855\item implementation of inherited virtual functions concerning parallelograms :
856\begin{lstlisting}
857virtual void computeMB();
858virtual std::vector<const Point*> boundNodes() const;
859virtual std::vector<Point*> nodes();
860virtual std::vector<const Point*> nodes() const;
861virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
862virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
863\end{lstlisting}
864\end{itemize}
865
866\subsubsection{The {\class Rectangle} class}
867
868To define a rectangle, you give 3 vertices, as for {\class Parallelogram}. You can also give the center (or the origin) and the lengths, or you can gives bounds. To store temporarily these values, the {\class Rectangle} class has internal attributes.
869
870\begin{lstlisting}
871class Rectangle : public Parallelogram
872{
873  protected:
874    Point center_, origin_;
875    bool isCenter_, isOrigin_;
876    real_t xlength_, ylength_;
877  private:
878    real_t xmin_, xmax_, ymin_, ymax_;
879    bool isBounds_;
880\end{lstlisting}
881
882\begin{figure}[H]
883\begin{center}
884\inputTikZ{geo-rectangle}
885\end{center}
886\end{figure}
887
888It offers constructors, taking from 3 to 7 {\class Parameter}. Authorized keys and what they do are:
889
890\begin{center}
891\begin{tabular}{|c|c|c|c|}
892\hline
893attribute & key to use to define it & optional & default value \\
894\hline
895\hline
896p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v4}) & no & no \\
897 & or ({\param \_center}, {\param \_xlength}, {\param \_ylength}) & & \\
898 & or ({\param \_origin}, {\param \_xlength}, {\param \_ylength}) & & \\
899 & or ({\param \_xmin}, {\param \_xmax}, {\param \_ymin}, {\param \_ymax}) & & \\
900n\_ & {\param \_nnodes} & yes & 2 \\
901h\_ & {\param \_hsteps} & yes & no \\
902domName\_ & {\param \_domain\_name} & yes & "" \\
903sideNames\_ & {\param \_side\_names} & yes & "" \\
904\hline
905\end{tabular}
906\end{center}
907
908It also offers :
909
910\begin{itemize}
911\item some accessors
912\item implementation of inherited virtual functions concerning  rectangles :
913\begin{lstlisting}
914virtual void computeMB();
915virtual std::vector<const Point*> boundNodes() const;
916virtual std::vector<Point*> nodes();
917virtual std::vector<const Point*> nodes() const;
918virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
919virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
920\end{lstlisting}
921\end{itemize}
922
923\subsubsection{The {\class Square} class}
924
925To define a square, you give 3 vertices. You can also give the center (or the origin) and the length.
926
927\begin{figure}[H]
928\begin{center}
929\inputTikZ{geo-square}
930\end{center}
931\end{figure}
932
933It offers constructors, taking from 2 to 6 {\class Parameter}. Authorized keys and what they do are:
934
935\begin{center}
936\begin{tabular}{|c|c|c|c|}
937\hline
938attribute & key to use to define it & optional & default value \\
939\hline
940\hline
941p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v4}) & no & no \\
942 & or ({\param \_center}, {\param \_xlength}, {\param \_ylength}) & & \\
943 & or ({\param \_origin}, {\param \_xlength}, {\param \_ylength}) & & \\
944n\_ & {\param \_nnodes} & yes & 2 \\
945h\_ & {\param \_hsteps} & yes & no \\
946domName\_ & {\param \_domain\_name} & yes & "" \\
947sideNames\_ & {\param \_side\_names} & yes & "" \\
948\hline
949\end{tabular}
950\end{center}
951
952It also offers :
953
954\begin{itemize}
955\item some accessors
956\item implementation of inherited virtual functions concerning squares :
957\begin{lstlisting}
958virtual void computeMB();
959virtual std::vector<const Point*> boundNodes() const;
960virtual std::vector<Point*> nodes();
961virtual std::vector<const Point*> nodes() const;
962virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
963virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
964\end{lstlisting}
965\end{itemize}
966
967\subsubsection{The {\class Ellipse} class}
968
969To define an elliptic surface, you have to define a center and 2 apogees ($c$, $p1$ and $p2$ in the following figure):
970
971\begin{figure}[H]
972\begin{center}
973\inputTikZ{geo-ellipse}
974\end{center}
975\end{figure}
976
977You can also define an ellipse from its center and the axis lengths, when axes are x-axis and y-axis.
978The main attributes are center\_, p1\_, p2\_, p3\_, p4\_, those that will be used in meshing.
979
980\begin{lstlisting}
981class Ellipse : public Surface
982{
983  protected:
984    Point center_, p1_, p2_, p3_, p4_;
985    real_t xlength_, ylength_;
986    bool isAxis_;
987    number_t n1_, n2_, n3_, n4_;
988    std::vector<real_t> h_;
989    real_t thetamin_;
990    real_t thetamax_;
991    dimen_t type_;
992\end{lstlisting}
993
994It offers constructors, taking from 3 to 9 {\class Parameter}. Authorized keys and what they do are:
995
996\begin{center}
997\begin{tabular}{|c|c|c|c|}
998\hline
999attribute & key to use to define it & optional & default value \\
1000\hline
1001\hline
1002center\_, p1\_, p2\_, p3\_, p4\_ & ({\param \_center}, {\param \_v1} and {\param \_v2}) & no & no \\
1003 & or ({\param \_center}, {\param \_xlength}, {\param \_ylength}) & & \\
1004n1\_, n2\_, n3\_, n4\_ & {\param \_nnodes} & yes & 2 \\
1005h\_ & {\param \_hsteps} & yes & no \\
1006domName\_ & {\param \_domain\_name} & yes & "" \\
1007sideNames\_ & {\param \_side\_names} & yes & "" \\
1008thetamin\_ & {\param \_angle1} & yes & 0 \\
1009thetamax\_ & {\param \_angle2} & yes & 360 \\
1010type\_ & {\param \_type} & yes & 1 \\
1011\hline
1012\end{tabular}
1013\end{center}
1014
1015It also offers :
1016
1017\begin{itemize}
1018\item some accessors
1019\item implementation of inherited virtual functions concerning elliptic surfaces and disks :
1020\begin{lstlisting}
1021virtual void computeMB();
1022virtual std::vector<const Point*> boundNodes() const;
1023virtual std::vector<Point*> nodes();
1024virtual std::vector<const Point*> nodes() const;
1025virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1026virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1027\end{lstlisting}
1028\end{itemize}
1029
1030\subsubsection{The {\class Disk} class}
1031
1032To define a disk, you have to do the same way as {\class Ellipse}.
1033
1034\begin{figure}[H]
1035\begin{center}
1036\inputTikZ{geo-disk}
1037\end{center}
1038\end{figure}
1039
1040It offers constructors, taking from 2 to 9 {\class Parameter}. Authorized keys and what they do are:
1041
1042\begin{center}
1043\begin{tabular}{|c|c|c|c|}
1044\hline
1045attribute & key to use to define it & optional & default value \\
1046\hline
1047\hline
1048center\_, p1\_, p2\_, p3\_, p4\_ & ({\param \_center}, {\param \_v1} and {\param \_v2}) & no & no \\
1049 & or ({\param \_center}, {\param \_radius}) & & \\
1050n1\_, n2\_, n3\_, n4\_ & {\param \_nnodes} & yes & 2 \\
1051h\_ & {\param \_hsteps} & yes & no \\
1052domName\_ & {\param \_domain\_name} & yes & "" \\
1053sideNames\_ & {\param \_side\_names} & yes & "" \\
1054thetamin\_ & {\param \_angle1} & yes & 0 \\
1055thetamax\_ & {\param \_angle2} & yes & 360 \\
1056type\_ & {\param \_type} & yes & 1 \\
1057\hline
1058\end{tabular}
1059\end{center}
1060
1061It also offers :
1062
1063\begin{itemize}
1064\item some accessors
1065\item implementation of inherited virtual functions concerning elliptic surfaces and disks :
1066\begin{lstlisting}
1067virtual void computeMB();
1068virtual std::vector<const Point*> boundNodes() const;
1069virtual std::vector<Point*> nodes();
1070virtual std::vector<const Point*> nodes() const;
1071virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1072virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1073\end{lstlisting}
1074\end{itemize}
1075
1076\subsubsection{The {\class Polyhedron} class}
1077
1078A polyhedron is defined bi its polygonal faces.
1079
1080\begin{lstlisting}
1081class Polyhedron : public Volume
1082{
1083  private:
1084    std::vector<Polygon*> faces_;
1085  protected:
1086    std::vector<Point> p_;
1087    std::vector<number_t> n_;
1088    std::vector<real_t> h_;
1089\end{lstlisting}
1090
1091p\_, n\_ and h\_ are defined here but used in child classes.
1092
1093\begin{figure}[H]
1094\begin{center}
1095\inputTikZ{geo-polyhedron}
1096\end{center}
1097\end{figure}
1098
1099It offers constructors, taking from 1 to 2 {\class Parameter}. Authorized keys and what they do are:
1100
1101\begin{center}
1102\begin{tabular}{|c|c|c|c|}
1103\hline
1104attribute & key to use to define it & optional & default value \\
1105\hline
1106\hline
1107faces\_ & {\param \_faces} & no & no \\
1108domName\_ & {\param \_domain\_name} & yes & "" \\
1109\hline
1110\end{tabular}
1111\end{center}
1112
1113It also offers :
1114
1115\begin{itemize}
1116\item some accessors
1117\item implementation of inherited virtual functions concerning hexahedron, parallelepipeds and cubes :
1118\begin{lstlisting}
1119virtual void computeMB();
1120virtual std::vector<const Point*> boundNodes() const;
1121virtual std::vector<Point*> nodes();
1122virtual std::vector<const Point*> nodes() const;
1123virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1124virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1125\end{lstlisting}
1126\end{itemize}
1127
1128\subsubsection{The {\class Tetrahedron} class}
1129
1130To define a tetrahedron, you just have to give the 4 vertices.
1131
1132\begin{figure}[H]
1133\begin{center}
1134\inputTikZ{geo-tetrahedron}
1135\end{center}
1136\end{figure}
1137
1138It offers constructors, taking from 4 to 7 {\class Parameter}. Authorized keys and what they do are:
1139
1140\begin{center}
1141\begin{tabular}{|c|c|c|c|}
1142\hline
1143attribute & key to use to define it & optional & default value \\
1144\hline
1145\hline
1146p\_ & {\param \_v1}, {\param \_v2}, {\param \_v3} and {\param \_v4} & no & no\\
1147n\_ & {\param \_nnodes} & yes & 2 \\
1148h\_ & {\param \_hsteps} & yes & no \\
1149domName\_ & {\param \_domain\_name} & yes & "" \\
1150sideNames\_ & {\param \_side\_names} & yes & "" \\
1151\hline
1152\end{tabular}
1153\end{center}
1154
1155It also offers :
1156
1157\begin{itemize}
1158\item some accessors
1159\item implementation of inherited virtual functions concerning tetrahedron :
1160\begin{lstlisting}
1161virtual void computeMB();
1162virtual std::vector<const Point*> boundNodes() const;
1163virtual std::vector<Point*> nodes();
1164virtual std::vector<const Point*> nodes() const;
1165virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1166virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1167\end{lstlisting}
1168\end{itemize}
1169
1170\subsubsection{The {\class Hexahedron} class}
1171
1172To define a hexahedron, you just have to give the 8 vertices.
1173
1174\begin{figure}[H]
1175\begin{center}
1176\inputTikZ{geo-hexahedron}
1177\end{center}
1178\end{figure}
1179
1180It offers constructors, taking from 8 to 11 {\class Parameter}. Authorized keys and what they do are:
1181
1182\begin{center}
1183\begin{tabular}{|c|c|c|c|}
1184\hline
1185attribute & key to use to define it & optional & default value \\
1186\hline
1187\hline
1188p\_ & {\param \_v1}, {\param \_v2}, {\param \_v3}, {\param \_v4}, {\param \_v5}, {\param \_v6}, {\param \_v7} and {\param \_v8} & no & no\\
1189n\_ & {\param \_nnodes} & yes & 2 \\
1190h\_ & {\param \_hsteps} & yes & no \\
1191domName\_ & {\param \_domain\_name} & yes & "" \\
1192sideNames\_ & {\param \_side\_names} & yes & "" \\
1193\hline
1194\end{tabular}
1195\end{center}
1196
1197It also offers :
1198
1199\begin{itemize}
1200\item some accessors
1201\item implementation of inherited virtual functions concerning hexahedra  :
1202\begin{lstlisting}
1203virtual void computeMB();
1204virtual std::vector<const Point*> boundNodes() const;
1205virtual std::vector<Point*> nodes();
1206virtual std::vector<const Point*> nodes() const;
1207virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1208virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1209\end{lstlisting}
1210\end{itemize}
1211
1212\subsubsection{The {\class Parallelepiped} class}
1213
1214To define a parallelepiped, you just have to give 4 vertices ($p_1$, $p_2$, $p_4$ and $p_5$ in the following figure):
1215
1216\begin{figure}[H]
1217\begin{center}
1218\inputTikZ{geo-parallelepiped}
1219\end{center}
1220\end{figure}
1221
1222It offers constructors, taking from 4 to 7 {\class Parameter}. Authorized keys and what they do are:
1223
1224\begin{center}
1225\begin{tabular}{|c|c|c|c|}
1226\hline
1227attribute & key to use to define it & optional & default value \\
1228\hline
1229\hline
1230p\_ & {\param \_v1}, {\param \_v2}, {\param \_v4} and {\param \_v5} & no & no\\
1231n\_ & {\param \_nnodes} & yes & 2 \\
1232h\_ & {\param \_hsteps} & yes & no \\
1233domName\_ & {\param \_domain\_name} & yes & "" \\
1234sideNames\_ & {\param \_side\_names} & yes & "" \\
1235\hline
1236\end{tabular}
1237\end{center}
1238
1239It also offers :
1240
1241\begin{itemize}
1242\item some accessors
1243\item implementation of inherited virtual functions concerning parallelepipeds :
1244\begin{lstlisting}
1245virtual void computeMB();
1246virtual std::vector<const Point*> boundNodes() const;
1247virtual std::vector<Point*> nodes();
1248virtual std::vector<const Point*> nodes() const;
1249virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1250virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1251\end{lstlisting}
1252\end{itemize}
1253
1254\subsubsection{The {\class Cuboid} class}
1255
1256To define a cuboid, you just have to give 4 vertices, as for {\class Parallelepiped}. And as for {\class Rectangle}, you can define it by its center or its origin and its lengths, or define it by its bounds
1257
1258\begin{figure}[H]
1259\begin{center}
1260\inputTikZ{geo-cuboid}
1261\end{center}
1262\end{figure}
1263
1264It offers constructors, taking from 4 to 9 {\class Parameter}. Authorized keys and what they do are:
1265
1266\begin{center}
1267\begin{tabular}{|c|c|c|c|}
1268\hline
1269attribute & key to use to define it & optional & default value \\
1270\hline
1271\hline
1272p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v4}, {\param \_v5}) & no & no \\
1273 & or ({\param \_center}, {\param \_xlength}, {\param \_ylength}, {\param \_zlength}) & & \\
1274 & or ({\param \_origin}, {\param \_xlength}, {\param \_ylength}, {\param \_zlength}) & & \\
1275 & or ({\param \_xmin}, {\param \_xmax}, {\param \_ymin}, {\param \_ymax}, {\param \_zmin}, {\param \_zmax}) & & \\
1276n\_ & {\param \_nnodes} & yes & 2 \\
1277h\_ & {\param \_hsteps} & yes & no \\
1278domName\_ & {\param \_domain\_name} & yes & "" \\
1279sideNames\_ & {\param \_side\_names} & yes & "" \\
1280\hline
1281\end{tabular}
1282\end{center}
1283
1284It also offers :
1285
1286\begin{itemize}
1287\item some accessors
1288\item implementation of inherited virtual functions concerning cuboids :
1289\begin{lstlisting}
1290virtual void computeMB();
1291virtual std::vector<const Point*> boundNodes() const;
1292virtual std::vector<Point*> nodes();
1293virtual std::vector<const Point*> nodes() const;
1294virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1295virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1296\end{lstlisting}
1297\end{itemize}
1298
1299\subsubsection{The {\class Cube} class}
1300
1301To define a cube, you just have to give 4 vertices. You can also give the center or the origin and the length.
1302
1303\begin{lstlisting}
1304class Cube : public Cuboid
1305{
1306  private:
1307    dimen_t nboctants_;
1308\end{lstlisting}
1309
1310\begin{figure}[H]
1311\begin{center}
1312\inputTikZ{geo-cube}
1313\end{center}
1314\end{figure}
1315
1316It offers constructors, taking from 2 to 7 {\class Parameter}. Authorized keys and what they do are:
1317
1318\begin{center}
1319\begin{tabular}{|c|c|c|c|}
1320\hline
1321attribute & key to use to define it & optional & default value \\
1322\hline
1323\hline
1324p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v4}, {\param \_v5}) & no & no \\
1325 & or ({\param \_center}, {\param \_length}) & & \\
1326 & or ({\param \_origin}, {\param \_length}) & & \\
1327n\_ & {\param \_nnodes} & yes & 2 \\
1328h\_ & {\param \_hsteps} & yes & no \\
1329domName\_ & {\param \_domain\_name} & yes & "" \\
1330sideNames\_ & {\param \_side\_names} & yes & "" \\
1331nboctants\_ & {\param \_nboctants} & yes & 8 \\
1332\hline
1333\end{tabular}
1334\end{center}
1335
1336It also offers :
1337
1338\begin{itemize}
1339\item some accessors
1340\item implementation of inherited virtual functions concerning cubes :
1341\begin{lstlisting}
1342virtual void computeMB();
1343virtual std::vector<const Point*> boundNodes() const;
1344virtual std::vector<Point*> nodes();
1345virtual std::vector<const Point*> nodes() const;
1346virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1347virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1348\end{lstlisting}
1349\end{itemize}
1350
1351\subsubsection{The {\class Ellipsoid} class}
1352
1353To define an ellipsoid, you do the same way as for an ellipse, namely defining a center and 3 apogees ($c$, $p1$, $p2$ and $p_6$ in the following figure):
1354
1355\begin{figure}[H]
1356\begin{center}
1357\inputTikZ{geo-ellipsoid}
1358\end{center}
1359\end{figure}
1360
1361\begin{lstlisting}
1362class Ellipsoid : public Volume
1363{
1364  protected:
1365    Point center_, p1_, p2_, p3_, p4_, p5_, p6_;
1366    real_t xlength_, ylength_, zlength_;
1367    bool isAxis_;
1368
1369    number_t n1_, n2_, n3_, n4_, n5_, n6_, n7_, n8_, n9_, n10_, n11_, n12_;
1370    std::vector<real_t> h_;
1371
1372    dimen_t nboctants_;
1373    dimen_t type_;
1374\end{lstlisting}
1375
1376It offers constructors, taking from 4 to 9 {\class Parameter}. Authorized keys and what they do are:
1377
1378\begin{center}
1379\begin{tabular}{|c|c|c|c|}
1380\hline
1381attribute & key to use to define it & optional & default value \\
1382\hline
1383\hline
1384p\_ & ({\param \_center}, {\param \_v1}, {\param \_v2}, {\param \_v6}) & no & no \\
1385 & or ({\param \_center}, {\param \_xlength}, {\param \_ylength}, {\param \_zlength}) & & \\
1386 & or ({\param \_origin}, {\param \_xlength}, {\param \_ylength}, {\param \_zlength}) & & \\
1387n\_ & {\param \_nnodes} & yes & 2 \\
1388h\_ & {\param \_hsteps} & yes & no \\
1389domName\_ & {\param \_domain\_name} & yes & "" \\
1390sideNames\_ & {\param \_side\_names} & yes & "" \\
1391nboctants\_ & {\param \_nboctants} & yes & 8 \\
1392type\_ & {\param \_type} & yes & 1 \\
1393\hline
1394\end{tabular}
1395\end{center}
1396
1397It also offers :
1398
1399\begin{itemize}
1400\item some accessors
1401\item implementation of inherited virtual functions concerning ellipsoids :
1402\begin{lstlisting}
1403virtual void computeMB();
1404virtual std::vector<const Point*> boundNodes() const;
1405virtual std::vector<Point*> nodes();
1406virtual std::vector<const Point*> nodes() const;
1407virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1408virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1409\end{lstlisting}
1410\end{itemize}
1411
1412\subsubsection{The {\class Ball} class}
1413
1414To define a ball, you do the same way as for an ellipsoid.
1415
1416\begin{figure}[H]
1417\begin{center}
1418\inputTikZ{geo-ball}
1419\end{center}
1420\end{figure}
1421
1422It offers constructors, taking from 2 to 9 {\class Parameter}. Authorized keys and what they do are:
1423
1424\begin{center}
1425\begin{tabular}{|c|c|c|c|}
1426\hline
1427attribute & key to use to define it & optional & default value \\
1428\hline
1429\hline
1430p\_ & ({\param \_center}, {\param \_v1}, {\param \_v2}, {\param \_v6}) & no & no \\
1431 & or ({\param \_center}, {\param \_length}) & & \\
1432 & or ({\param \_origin}, {\param \_length}) & & \\
1433n\_ & {\param \_nnodes} & yes & 2 \\
1434h\_ & {\param \_hsteps} & yes & no \\
1435domName\_ & {\param \_domain\_name} & yes & "" \\
1436sideNames\_ & {\param \_side\_names} & yes & "" \\
1437nboctants\_ & {\param \_nboctants} & yes & 8 \\
1438type\_ & {\param \_type} & yes & 1 \\
1439\hline
1440\end{tabular}
1441\end{center}
1442
1443It also offers :
1444
1445\begin{itemize}
1446\item some accessors
1447\item implementation of inherited virtual functions concerning balls :
1448\begin{lstlisting}
1449virtual void computeMB();
1450virtual std::vector<const Point*> boundNodes() const;
1451virtual std::vector<Point*> nodes();
1452virtual std::vector<const Point*> nodes() const;
1453virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1454virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1455\end{lstlisting}
1456\end{itemize}
1457
1458\subsubsection{The {\class Trunk} class}
1459
1460A trunk is a generalized truncated cone. To define a trunk, you need to give a surface, namely a polygonal surface or a elliptical surface. To define the other surface, you just need to give a point of this surface ($origin$), and the scale factor according to the first surface.
1461
1462For a trunk with polygonal basis, $origin$ is the equivalent of the first vertex of the surface you give, as you can see on the following figure of a trunk with triangular basis. The triangle being defined by its vertices $p_1$, $p_2$ and $p_3$, $origin$ is the equivalent of $p_1$:
1463
1464\begin{figure}[H]
1465\begin{center}
1466\inputTikZ{geo-trunk-triangle}
1467\end{center}
1468\end{figure}
1469
1470For a trunk with elliptical basis, $origin$ is the center of the second basis, as you can see on the following figure of a trunk with elliptical basis.
1471
1472\begin{figure}[H]
1473\begin{center}
1474\inputTikZ{geo-trunk-ellipse}
1475\end{center}
1476\end{figure}
1477
1478\begin{lstlisting}
1479class Trunk : public Volume
1480{
1481  protected:
1482    Surface * basis_;
1483    real_t scale_;
1484    std::vector<Point> p_;
1485    std::vector<number_t> n_;
1486    std::vector<real_t> h_;
1487    Point origin_, center1_, p1_, p2_;
1488    bool isElliptical_, isN_;
1489\end{lstlisting}
1490
1491It offers constructors, taking from 3 to 8 {\class Parameter}. Authorized keys and what they do are:
1492
1493\begin{center}
1494\begin{tabular}{|c|c|c|c|}
1495\hline
1496attribute & key to use to define it & optional & default value \\
1497\hline
1498\hline
1499p\_ & ({\param \_center1}, {\param \_v1}, {\param \_v2}, {\param \_center2}, {\param \_scale}) & no & no \\
1500 & or ({\param \_basis}, {\param \_origin}, {\param \_scale}) & & \\
1501n\_ & {\param \_nnodes} & yes & 2 \\
1502h\_ & {\param \_hsteps} & yes & no \\
1503domName\_ & {\param \_domain\_name} & yes & "" \\
1504sideNames\_ & {\param \_side\_names} & yes & "" \\
1505\hline
1506\end{tabular}
1507\end{center}
1508
1509It also offers :
1510
1511\begin{itemize}
1512\item some accessors
1513\item implementation of inherited virtual functions concerning trunks :
1514\begin{lstlisting}
1515virtual void computeMB();
1516virtual std::vector<const Point*> boundNodes() const;
1517virtual std::vector<Point*> nodes();
1518virtual std::vector<const Point*> nodes() const;
1519virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1520virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1521\end{lstlisting}
1522\end{itemize}
1523
1524\subsubsection{The {\class Cylinder} class}
1525
1526A cylinder is a truncated cone whose apex is at infinite distance. So it is the geometry defined by the extrusion of a surface by translation.
1527
1528\begin{lstlisting}
1529class Cylinder : public Trunk
1530{
1531  protected:
1532    //! direction vector
1533    Point dir_;
1534\end{lstlisting}
1535
1536\begin{figure}[H]
1537\begin{center}
1538\inputTikZ{geo-cylinder}
1539\end{center}
1540\end{figure}
1541
1542It offers constructors, taking from 2 to 7 {\class Parameter}. Authorized keys and what they do are:
1543
1544\begin{center}
1545\begin{tabular}{|c|c|c|c|}
1546\hline
1547attribute & key to use to define it & optional & default value \\
1548\hline
1549\hline
1550p\_ & ({\param \_center1}, {\param \_v1}, {\param \_v2}, {\param \_center2}) & no & no \\
1551 & or ({\param \_basis}, {\param \_direction}) & & \\
1552n\_ & {\param \_nnodes} & yes & 2 \\
1553h\_ & {\param \_hsteps} & yes & no \\
1554domName\_ & {\param \_domain\_name} & yes & "" \\
1555sideNames\_ & {\param \_side\_names} & yes & "" \\
1556\hline
1557\end{tabular}
1558\end{center}
1559
1560It also offers :
1561
1562\begin{itemize}
1563\item some accessors
1564\item implementation of inherited virtual functions concerning cylinders :
1565\begin{lstlisting}
1566virtual void computeMB();
1567virtual std::vector<const Point*> boundNodes() const;
1568virtual std::vector<Point*> nodes();
1569virtual std::vector<const Point*> nodes() const;
1570virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1571virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1572\end{lstlisting}
1573\end{itemize}
1574
1575\subsubsection{The {\class Prism} class}
1576
1577A prism is by definition a cylinder whose basis is a polygonal surface. Often a prism refers to a cylinder with triangular basis (as the finite element cell).
1578
1579\begin{lstlisting}
1580class Prism : public Cylinder
1581{
1582  private:
1583    bool isTriangular_;
1584    Point p3_;
1585\end{lstlisting}
1586
1587\begin{figure}[H]
1588\begin{center}
1589\inputTikZ{geo-prism}
1590\end{center}
1591\end{figure}
1592
1593It offers constructors, taking from 2 to 7 {\class Parameter}. Authorized keys and what they do are:
1594
1595\begin{center}
1596\begin{tabular}{|c|c|c|c|}
1597\hline
1598attribute & key to use to define it & optional & default value \\
1599\hline
1600\hline
1601p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v3}, {\param \_direction}) & no & no \\
1602 & or ({\param \_basis}, {\param \_direction}) & & \\
1603n\_ & {\param \_nnodes} & yes & 2 \\
1604h\_ & {\param \_hsteps} & yes & no \\
1605domName\_ & {\param \_domain\_name} & yes & "" \\
1606sideNames\_ & {\param \_side\_names} & yes & "" \\
1607\hline
1608\end{tabular}
1609\end{center}
1610
1611It also offers :
1612
1613\begin{itemize}
1614\item some accessors
1615\item implementation of inherited virtual functions concerning prisms :
1616\begin{lstlisting}
1617virtual void computeMB();
1618virtual std::vector<const Point*> boundNodes() const;
1619virtual std::vector<Point*> nodes();
1620virtual std::vector<const Point*> nodes() const;
1621virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1622virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1623\end{lstlisting}
1624\end{itemize}
1625
1626\subsubsection{The {\class Cone} class}
1627
1628A cone is defined by a surface and an apex. As for trunks and cylinders, you can also define directly a cone with elliptical basis.
1629
1630\begin{figure}[H]
1631\begin{center}
1632\inputTikZ{geo-cone}
1633\end{center}
1634\end{figure}
1635
1636It offers constructors, taking from 3 to 7 {\class Parameter}. Authorized keys and what they do are:
1637
1638\begin{center}
1639\begin{tabular}{|c|c|c|c|}
1640\hline
1641attribute & key to use to define it & optional & default value \\
1642\hline
1643\hline
1644p\_ & ({\param \_center1}, {\param \_v1}, {\param \_v2}, {\param \_apex}) & no & no \\
1645 & or ({\param \_basis}, {\param \_apex}) & & \\
1646n\_ & {\param \_nnodes} & yes & 2 \\
1647h\_ & {\param \_hsteps} & yes & no \\
1648domName\_ & {\param \_domain\_name} & yes & "" \\
1649sideNames\_ & {\param \_side\_names} & yes & "" \\
1650\hline
1651\end{tabular}
1652\end{center}
1653
1654It also offers :
1655
1656\begin{itemize}
1657\item some accessors
1658\item implementation of inherited virtual functions concerning cones :
1659\begin{lstlisting}
1660virtual void computeMB();
1661virtual std::vector<const Point*> boundNodes() const;
1662virtual std::vector<Point*> nodes();
1663virtual std::vector<const Point*> nodes() const;
1664virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1665virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1666\end{lstlisting}
1667\end{itemize}
1668
1669\subsubsection{The {\class Pyramid} class}
1670
1671A pyramid is a cone with a polygonal basis. Often a pyramid refers to a cone with quadrangular basis (as the finite element cell).
1672
1673\begin{lstlisting}
1674class Pyramid : public Cone
1675{
1676  private:
1677    bool isQuadrangular_;
1678    Point p3_, p4_;
1679\end{lstlisting}
1680
1681\begin{figure}[H]
1682\begin{center}
1683\inputTikZ{geo-pyramid}
1684\end{center}
1685\end{figure}
1686
1687It offers constructors, taking from 2 to 8 {\class Parameter}. Authorized keys and what they do are:
1688
1689\begin{center}
1690\begin{tabular}{|c|c|c|c|}
1691\hline
1692attribute & key to use to define it & optional & default value \\
1693\hline
1694\hline
1695p\_ & ({\param \_v1}, {\param \_v2}, {\param \_v3}, {\param \_v4}, {\param \_apex}) & no & no \\
1696 & or ({\param \_basis}, {\param \_apex}) & & \\
1697n\_ & {\param \_nnodes} & yes & 2 \\
1698h\_ & {\param \_hsteps} & yes & no \\
1699domName\_ & {\param \_domain\_name} & yes & "" \\
1700sideNames\_ & {\param \_side\_names} & yes & "" \\
1701\hline
1702\end{tabular}
1703\end{center}
1704
1705It also offers :
1706
1707\begin{itemize}
1708\item some accessors
1709\item implementation of inherited virtual functions concerning pyramids :
1710\begin{lstlisting}
1711virtual void computeMB();
1712virtual std::vector<const Point*> boundNodes() const;
1713virtual std::vector<Point*> nodes();
1714virtual std::vector<const Point*> nodes() const;
1715virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1716virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1717\end{lstlisting}
1718\end{itemize}
1719
1720\subsubsection{The {\class RevTrunk} class}
1721
1722A revolution trunk is a right trunk with circular basis. {\class RevTrunk} offers you more geometry abilities. Indeed, you can decide to add extensions at ends of  the revolution trunk. Extensions can be : none, flat, ellipsoid, or cone.
1723
1724
1725\begin{lstlisting}
1726class RevTrunk : public Trunk
1727{
1728  protected:
1729    real_t radius1_, radius2_;
1730    number_t nbSubDomains_;
1731    GeometricEndShape endShape1_;
1732    real_t distance1_;
1733    GeometricEndShape endShape2_;
1734    real_t distance2_;
1735    dimen_t type_;
1736\end{lstlisting}
1737
1738\begin{figure}[H]
1739\begin{center}
1740\inputTikZ{geo-revtrunk}
1741\end{center}
1742\end{figure}
1743
1744It offers constructors, taking from 4 to 13 {\class Parameter}. Authorized keys and what they do are:
1745
1746\begin{center}
1747\begin{tabular}{|c|c|c|c|}
1748\hline
1749attribute & key to use to define it & optional & default value \\
1750\hline
1751\hline
1752p\_ & ({\param \_center1}, {\param \_center2}, {\param \_radius1}, {\param \_radius2}) & no & no \\
1753endShape1\_ & {\param \_end1\_shape} &  yes & \_gesflat\\
1754endShape2\_ & {\param \_end2\_shape} & yes & \_gesFlat \\
1755distance1\_ & {\param \_end1\_distance} &  yes & 0. \\
1756distance2\_ & {\param \_end2\_distance} & yes & 0. \\
1757n\_ & {\param \_nnodes} & yes & 2 \\
1758h\_ & {\param \_hsteps} & yes & no \\
1759domName\_ & {\param \_domain\_name} & yes & "" \\
1760sideNames\_ & {\param \_side\_names} & yes & "" \\
1761nbSubDomains\_ & {\param \_nbsubdomains} & yes & 1 \\
1762type\_ & {\param \_type} & yes & 1 \\
1763\hline
1764\end{tabular}
1765\end{center}
1766
1767It also offers :
1768
1769\begin{itemize}
1770\item some accessors
1771\item implementation of inherited virtual functions concerning revolution trunks :
1772\begin{lstlisting}
1773virtual void computeMB();
1774virtual std::vector<const Point*> boundNodes() const;
1775virtual std::vector<Point*> nodes();
1776virtual std::vector<const Point*> nodes() const;
1777virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1778virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1779\end{lstlisting}
1780\end{itemize}
1781
1782\subsubsection{The {\class RevCylinder} class}
1783
1784A revolution cylinder is a revolution trunk where both radiuses are equal. So, we need centers of both bases, and the radius. As {\class RevTrunk}, {\class RevCylinder} offers you the ability to add extensions at ends of  the revolution cylinder.
1785
1786\begin{figure}[H]
1787\begin{center}
1788\inputTikZ{geo-revcylinder}
1789\end{center}
1790\end{figure}
1791
1792It offers constructors, taking from 3 to 12 {\class Parameter}. Authorized keys and what they do are:
1793
1794\begin{center}
1795\begin{tabular}{|c|c|c|c|}
1796\hline
1797attribute & key to use to define it & optional & default value \\
1798\hline
1799\hline
1800p\_ & ({\param \_center1}, {\param \_center2}, {\param \_radius}) & no & no \\
1801endShape1\_ & {\param \_end1\_shape} &  yes & \_gesflat\\
1802endShape2\_ & {\param \_end2\_shape} & yes & \_gesFlat \\
1803distance1\_ & {\param \_end1\_distance} &  yes & 0. \\
1804distance2\_ & {\param \_end2\_distance} & yes & 0. \\
1805n\_ & {\param \_nnodes} & yes & 2 \\
1806h\_ & {\param \_hsteps} & yes & no \\
1807domName\_ & {\param \_domain\_name} & yes & "" \\
1808sideNames\_ & {\param \_side\_names} & yes & "" \\
1809nbSubDomains\_ & {\param \_nbsubdomains} & yes & 1 \\
1810type\_ & {\param \_type} & yes & 1 \\
1811\hline
1812\end{tabular}
1813\end{center}
1814
1815It also offers :
1816
1817\begin{itemize}
1818\item some accessors
1819\item implementation of inherited virtual functions concerning revolution cylinders :
1820\begin{lstlisting}
1821virtual void computeMB();
1822virtual std::vector<const Point*> boundNodes() const;
1823virtual std::vector<Point*> nodes();
1824virtual std::vector<const Point*> nodes() const;
1825virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1826virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1827\end{lstlisting}
1828\end{itemize}
1829
1830\subsubsection{The {\class RevCone} class}
1831
1832A revolution cone is a revolution trunk where second radius is equal to 0. As {\class RevTrunk},  {\class RevCone} offers you more the ability to add an extension to the basis of a revolution cone.
1833
1834\begin{figure}[H]
1835\begin{center}
1836\inputTikZ{geo-revcone}
1837\end{center}
1838\end{figure}
1839
1840It offers constructors, taking from 3 to 10 {\class Parameter}. Authorized keys and what they do are:
1841
1842\begin{center}
1843\begin{tabular}{|c|c|c|c|}
1844\hline
1845attribute & key to use to define it & optional & default value \\
1846\hline
1847\hline
1848p\_ & ({\param \_center}, {\param \_radius}, {\param \_apex}) & no & no \\
1849endShape1\_ & {\param \_end\_shape} &  yes & \_gesflat\\
1850distance1\_ & {\param \_end\_distance} &  yes & 0. \\
1851n\_ & {\param \_nnodes} & yes & 2 \\
1852h\_ & {\param \_hsteps} & yes & no \\
1853domName\_ & {\param \_domain\_name} & yes & "" \\
1854sideNames\_ & {\param \_side\_names} & yes & "" \\
1855nbSubDomains\_ & {\param \_nbsubdomains} & yes & 1 \\
1856type\_ & {\param \_type} & yes & 1 \\
1857\hline
1858\end{tabular}
1859\end{center}
1860
1861It also offers :
1862
1863\begin{itemize}
1864\item some accessors
1865\item implementation of inherited virtual functions concerning revolution cones :
1866\begin{lstlisting}
1867virtual void computeMB();
1868virtual std::vector<const Point*> boundNodes() const;
1869virtual std::vector<Point*> nodes();
1870virtual std::vector<const Point*> nodes() const;
1871virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > curves() const;
1872virtual std::vector<std::pair<ShapeType,std::vector<const Point*> > > surfs() const;
1873\end{lstlisting}
1874\end{itemize}
1875
1876\subsubsection{The {\class SetOfElems} class}
1877
1878