Geometric element#
A mesh consists mainly of a list of nodes and a list of geometric elements. XLiFE++ makes a distinction between two types of geometric element:
the complete element with all geometric information (
MeshElement
class)the side element that refers to a side of an other geometric element (
GeomElement
class)
Actually, complete elements are GeomElement
objects that handle a MeshElement
allocated pointer whereas side elements are
GeomElement
objects that handle a list of parent elements (a pair of GeomElement
* and side number). By default, MeshElement
pointer is not allocated for side element,
but it can be allocated using the GeomElement
::buildSideMeshElement
function.

In advanced use, users may have to manipulate GeomElement
/MeshElement
objects.
|
unique id number of geometric element |
|
|
|
list of parents ( |
|
material id ( |
|
domain id ( |
|
useful to mark an element ( |
|
local curvatures angles (x-y plane, x-z plane); ( |
|
twin |
|
a free pointer (void*) allowing users to attach data of any type to the Element |
|
a free void pointers vector ( |
Note
userData_p
and userData_ps
are free pointers dedicated to users to do specific calculations that are not managed by XLiFE++;
the casting and memory management associated with these pointers are the responsibility of users.
|
a subset of |
|
a subset of mesh node numbers (direct access to node numbers, |
|
a subset of mesh node numbers ( |
|
length or area or volume of element and of its sides ( |
|
characteristic size of element ( |
|
centroid of the element ( |
|
element orientation (sign of jacobian determinant, 0 means unset orientation) |
|
pointer to the reference element ( |
|
pointer to a geometric calculation structure ( |
|
true if the geometric map \(\varphi_\ell\) from reference element is linear (for instance first order simplicial element) |
The GeomElement
class collects all the member functions related to geometric element being either a mesh element or a side element.
returns element number as |
|
return the dimension of the element and the nodes dimansion as |
|
|
return the shape type of element (s=0) and its sides s>0 (_segment, _triangle, …) |
return |
|
create a |
|
return number of parents if a side element |
|
return return i-th parent and its side (no check) or access to |
|
|
returns number of nodes, of vertices, of sides, of sides of side of element |
|
return the global number of vertex/node i (i=1,…) |
|
return the global numbers of vertices/nodes if s=0, of side s else (as a |
|
return as |
return the adjacent element by side s or vertex i as a |
|
return the center of element (as a |
|
|
return outward normal/tangent vector on side s>0, if s=0 normal/tangent vector to element |
|
return true if mesh element contains a point pt given as a |
|
projection of a point onto mesh element (return a |
|
display information related to current mesh element |
Most of member functions of GeomElement
call member functions of MeshElement
.
But there are specific functions attached to MeshElement
:
|
return reference element ( as a |
|
build |
|
clear |
compute measure of element and its side |
|
compute orientation of element: sign(det(jacobian)) |
MeshElement
manages a pointer to a GeomMapData
object that handles data and tools attached
to the map \(\varphi_\ell\) from reference element to geometric element:
|
pointer to |
|
point in reference element used for calculation |
|
jacobian matrix of mapping from reference element to current one, inverse jacobian matrix (when needed) |
|
jacobian determinant and differential element for elementary integrals |
|
unit outward normal/tangent/bi-tangent vector at a given boundary point |
|
symetric metric tensor (t_i|t_j) and metric_tensor determinant |
|
mapping from ref. elt onto geom. elt (return a |
|
inverse mapping from elt onto ref. elt. |
|
Piola mapping from ref. elt onto geom. elt. Element (shape functions known) |
|
return the covariant/contraariant Piola map matrix at current point ( |
|
compute jacobian matrix at a point (return a |
|
compute jacobian determinant (return a |
|
compute inverse of jacobian matrix (jacobian already computed) |
|
compute differential element assuming jacobian matrix is update |
|
computes (oriented) normal vector at a point of element, assuming jacobian matrix is update |
|
compute tangent (and bitangent) vector at currernr point |
|
compute metric tensor for surface differential geometry |
|
compute surface gradient gradS from 2D reference gradient (grad0, grad1) |
|
display information |