Function

CGAL::lloyd_optimize_mesh_3

Definition

The function lloyd_optimize_mesh_3 is a mesh optimization process based on the minimization of a global energy function.

In lloyd_optimize_mesh_3, the minimized global energy may be interpreted as the L1-norm of the error achieved when the function x2 is interpolated on the mesh domain using a piecewise linear function which is linear in each cell of the Voronoi diagram of the mesh vertices.

The optimizer lloyd_optimize_mesh_3 works in iterative steps. At each iteration, mesh vertices are moved into positions that bring to zero the energy gradient and the Delaunay triangulation is updated. Vertices on the mesh boundaries are handled in a special way so as to preserve an accurate representation of the domain boundaries.

#include <CGAL/lloyd_optimize_mesh_3.h>

template<typename C3T3, typename MeshDomain_3>
Mesh_optimization_return_code
lloyd_optimize_mesh_3 ( C3T3& c3t3,
MeshDomain_3 domain,
double parameters::time_limit=0,
std::size_t parameters::max_iteration_number=0,
double parameters::convergence=0.02,
double parameters::freeze_bound = 0.01)

Precondition

time_limit 0 and 0 convergence 1 and 0 freeze_bound 1

Parameters

Parameter C3T3 is required to be a model of the concept MeshComplex_3InTriangulation_3. The argument c3t3, passed by reference, provides the initial mesh and is modified by the algorithm to represent the final optimized mesh.

Parameter MeshDomain_3 is required to be a model of the concept MeshDomain_3. The argument domain must be the MeshDomain_3 object used to create the c3t3 parameter.

The function has four optional parameters which are named parameters (we use the Boost.Parameter library). Therefore, when calling the function, the parameters can be provided in any order provided that the names of the parameters are used (see example at the bottom of this page).

Return Values

The function lloyd_optimize_mesh_3 returns a value of type Mesh_optimization_return_code which is:

Example

// Lloyd-smoothing until convergence reaches 0.01, freezing vertices which
// move less than 0.001*shortest_incident_edge_length
lloyd_optimize_mesh_3(c3t3, domain, parameters::convergence=0.01,
                      parameters::freeze_bound=0.001);

See Also

CGAL::Mesh_optimization_return_code
CGAL::make_mesh_3
CGAL::refine_mesh_3
CGAL::exude_mesh_3
CGAL::perturb_mesh_3
CGAL::odt_optimize_mesh_3