Level-sets, Xfem, fictitious domains
getfem++ offers (since v2.0) a certain number of
functionalities concerning level-sets, support for Xfem and
fictitious domain methods and discontinuous field across a
level-set.
Important : All the tools listed below needs the package
qhull installed on your system.
This package is widely available. It computes convex hull and
delaunay triangulations in arbitrary dimension. Everything here is
considered "work in progress", it is still subject to major changes
if needed.
The program crack.cc on the tests directory of the distribution is a good example of use of these tools.
Representation of level-sets
getfem++ deals with level-set defined by piecewise polynomial function on a mesh. It will be defined as the zero of this function. In the file getfem/getfem_level_set.h a level-set is represented by a function defined on a lagrange fem of a certain degree on a mesh. The constructor to define a new getfem::level_set is the following :
getfem::level_set ls(mesh, degree = 1, with_secondary = false);
where mesh is a valid mesh of type getfem::mesh, degree is the degree of the polynomials (1 is the defaule value), and with_secondary is a boolean whose default value is false. The secondary level-set is used to represent fractures (if p(x) is the primary levelset function and s(x) is the secondary levelset function, the crack is defined by p(x) = 0 and s(x) <= 0 : the role of the secondary is to stop the crack).
Each level-set function is defined by a mesh_fem mf and the dof values over this mesh_fem, in a vector. The object getfem::level_set contains a mesh_fem and the vectors of dof for the corresponding function(s). The method ls.value(0) returns the vector of dof for the primary level-set function, so that these values can be set. The method ls.value(1) returns the dof vector for the secondary level-set function if any. The method ls.get_mesh_fem() returns a reference on the getfem::mesh_fem object.
Mesh cut by level-sets
In order to compute adapted integration methods and finite element methods to represent a field which is discontinuous accross a level-set, a certain number of pre-computations have to be done at the mesh level. The file getfem/getfem_mesh_level_set.h defines the object getfem::mesh_level_set which handles these pre-computations. The constructor of this object is the following:
getfem::mesh_level_set mls(mesh);
where mesh is a valid mesh of
type getfem::mesh. In order to
indicate that the mesh is cut by a level-set, one has to call the
method mls.add_level_set(ls), where
ls is an object of type getfem::level_set. An arbitrary number of
level-sets can be added. To initialize the object or to actualize
it when the value of the level-set function is modified, one has to
call the method mls.adapt().
In particular a subdivision of each element cut by the level-set is made with simplices.
Adapted integration methods
For fields which are discontinuous across a level-set, integration methods have to be adapted. The object getfem::mesh_im_level_set defined in the file getfem/getfem_mesh_im_level_set.h defines a composite integration method for the elements cut by the level-set. The constructor of this object is the following:
getfem::mesh_im_level_set mim(mls, where, regular_im = 0, singular_im = 0);
where mls is an object of type getfem::mesh_level_set, where is an enum for which possible values are
- getfem::mesh_im_level_set::INTEGRATE_INSIDE (integrate over p(x)<0),
- getfem::mesh_im_level_set::INTEGRATE_OUTSIDE (integrate over p(x)>0),
- getfem::mesh_im_level_set::INTEGRATE_ALL,
- getfem::mesh_im_level_set::INTEGRATE_BOUNDARY (integrate over p(x)=0 and s(x) <= 0)
The argument regular_im should be
of type pintegration_method, and
will be the integration method applied on each sub-simplex of the
composite integration for convexes cut by the levelset. The
optional singular_im should be also
of type pintegration_method and is
used for crack singular functions: it is applied to sub-simplices
which share a vertex with the crack tip (the specific integration
method IM_QUASI_POLAR(..) is well
suited for this purpose).
The object getfem::mesh_im_level_set
can be used as a classical getfem::mesh_im object (for instance the
method mim.set_integration_method(...) allows to set
the integration methods for the elements which are not cut by the
level-set).
To initialize the object or to actualize it when the value of
the level-set function is modified, one has to call the method
mim.adapt().
Discontinuous field across some level-sets
The object getfem::mesh_fem_level_set is defined in the file getfem/getfem_mesh_fem_level_set.h. It is derivated from getfem::mesh_fem object and can be used in the same way. It defines a finite element method with discontinuity across the level-sets (it can deal with an arbitrary number of level-sets). The constructor is the following:
getfem::mesh_fem_level_set mfls(mesh, mf);
where mesh is a valid mesh of
type getfem::mesh and mf is the an object of type getfem::mesh_fem which defines the finite
element method used for elements which are not cut by the
level-sets.
To initialize the object or to actualize it when the value of
the level-set function is modified, one has to call the method
mfls.adapt().
To represent discontinuous fields, the finite element method is enriched with discontinuous functions which are the product of a Heaviside function by the base functions of the finite element method represented by mf (see [4] for more details).
Fictitious domain approach with Xfem
An example of a Poisson problem with a Dirichlet condition posed on a boundary independant of the mesh is present on the tests directory of the distribution. See xfem_dirichlet.cc file.
