clean lloyd_optimize_mesh_3

This commit is contained in:
Sébastien Loriot 2022-09-20 14:27:30 +02:00
parent 928ef0b447
commit e23d77bbf7
1 changed files with 113 additions and 117 deletions

View File

@ -32,123 +32,119 @@
namespace CGAL {
/*!
\ingroup PkgMesh3Functions
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 \f$ L^1\f$-norm of the error achieved
when the function \f$ x^2\f$ 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.
\pre `time_limit` \f$ \geq\f$ 0 and 0 \f$ \leq\f$ `convergence` \f$ \leq\f$ 1 and 0 \f$ \leq\f$ `freeze_bound` \f$ \leq\f$ 1
\tparam 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.
\tparam MD is required to be a model of the concept
`MeshDomain_3`. The argument `domain` must be the `MD`
object used to create the `c3t3` parameter.
\tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
@param c3t3 the initial mesh that will be modified by the algorithm to represent the final optimized mesh.
@param domain the domain to be discretized
@param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below:
\cgalNamedParamsBegin
\cgalParamNBegin{time_limit}
\cgalParamDescription{to set up, in seconds, a CPU time limit after which the optimization process is stopped.
This time is measured using `CGAL::Real_timer`. 0 means that there is no time limit.}
\cgalParamType{`double`}
\cgalParamExtra{\pre `time_limit` \f$ \geq\f$ 0}
\cgalParamDefault{0}
\cgalParamNBegin{max_iteration_number}
\cgalParamDescription{limit on the number of performed iterations. 0 means that there is
no limit on the number of performed iterations.}
\cgalParamExtra{\pre `max_iteration_number >=0`}
\cgalParamType{`int`}
\cgalParamDefault{0}
\cgalParamNBegin{freeze_bound}
\cgalParamDescription{designed to reduce running time of each optimization iteration.
Any vertex that has a displacement less than a given fraction of the length
of its shortest incident edge, is frozen (i.e.\ is not relocated).
The parameter `freeze_bound` gives the threshold ratio.
If it is set to 0, freezing of vertices is disabled.}
\cgalParamExtra{\pre `0<= freeze_bound <=1}
\cgalParamType{`double`}
\cgalParamDefault{0.001}
\cgalParamNBegin{convergence}
\cgalParamDescription{threshold ratio of stopping criterion based on convergence: the optimization process is stopped
when at the last iteration the displacement of any vertex is less than
a given fraction of the length of the shortest edge incident to that vertex.}
\cgalParamExtra{\pre `0 <=convergence <= 1`}
\cgalParamType{`double`}
\cgalParamDefault{0.001}
\cgalParamNBegin{do_freeze}
\cgalParamDescription{completes the `freeze_bound` parameter. If it is set to `true` (default value),
frozen vertices will not move anymore in next iterations. Otherwise, at each iteration, any vertex that
moves, unfreezes all its incident vertices.}
\cgalParamType{`bool`}
\cgalParamDefault{true}
\cgalNamedParamsEnd
\return
The function `lloyd_optimize_mesh_3()` returns a value of type `CGAL::Mesh_optimization_return_code`
which is:
<UL>
<LI>`CGAL::TIME_LIMIT_REACHED` when the time limit is reached.
<LI>`CGAL::MAX_ITERATION_NUMBER_REACHED` when `lloyd_optimize_mesh_3()` stops because it has performed `max_iteration_number` iterations.
<LI>`CGAL::CONVERGENCE_REACHED` when `lloyd_optimize_mesh_3()` stops because the convergence criterion
is achieved.
<LI>`CGAL::ALL_VERTICES_FROZEN` when all vertices have been frozen, when the
`do_freeze` parameter is set to true.
<LI>`CGAL::CANT_IMPROVE_ANYMORE` when `lloyd_optimize_mesh_3()` stops because
most vertices have been frozen, and no better convergence can be reached.
</UL>
\cgalHeading{Example}
\code{.cpp}
// 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,
parameters::do_freeze=true);
\endcode
\sa `CGAL::Mesh_optimization_return_code`
\sa `CGAL::make_mesh_3()`
\sa `CGAL::refine_mesh_3()`
\sa `CGAL::exude_mesh_3()`
\sa `CGAL::perturb_mesh_3()`
\sa `CGAL::odt_optimize_mesh_3()`
\note This function requires the \ref thirdpartyEigen library.
*/
* \ingroup PkgMesh3Functions
*
* 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 \f$ L^1\f$-norm of the error achieved
* when the function \f$ x^2\f$ 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.
*
* \tparam 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.
*
* \tparam MD is required to be a model of the concept
* `MeshDomain_3`. The argument `domain` must be the `MD`
* object used to create the `c3t3` parameter.
*
* \tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
*
* @param c3t3 the initial mesh that will be modified by the algorithm to represent the final optimized mesh.
* @param domain the domain to be discretized
* @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below:
*
* \cgalNamedParamsBegin
* \cgalParamNBegin{time_limit}
* \cgalParamDescription{to set up, in seconds, a CPU time limit after which the optimization process is stopped.
* This time is measured using `CGAL::Real_timer`. 0 means that there is no time limit.}
* \cgalParamType{`double`}
* \cgalParamExtra{\pre `time_limit >= 0`}
* \cgalParamDefault{0}
* \cgalParamNEnd
* \cgalParamNBegin{max_iteration_number}
* \cgalParamDescription{limit on the number of performed iterations. 0 means that there is
* no limit on the number of performed iterations.}
* \cgalParamExtra{\pre `max_iteration_number >=0`}
* \cgalParamType{`int`}
* \cgalParamDefault{0}
* \cgalParamNEnd
* \cgalParamNBegin{freeze_bound}
* \cgalParamDescription{designed to reduce running time of each optimization iteration.
* Any vertex that has a displacement less than a given fraction of the length
* of its shortest incident edge, is frozen (i.e.\ is not relocated).
* The parameter `freeze_bound` gives the threshold ratio.
* If it is set to 0, freezing of vertices is disabled.}
* \cgalParamExtra{\pre `0<= freeze_bound <=1`}
* \cgalParamType{`double`}
* \cgalParamDefault{0.001}
* \cgalParamNEnd
* \cgalParamNBegin{convergence}
* \cgalParamDescription{threshold ratio of stopping criterion based on convergence: the optimization process is stopped
* when at the last iteration the displacement of any vertex is less than
* a given fraction of the length of the shortest edge incident to that vertex.}
* \cgalParamExtra{\pre `0 <=convergence <= 1`}
* \cgalParamType{`double`}
* \cgalParamDefault{0.02}
* \cgalParamNEnd
* \cgalParamNBegin{do_freeze}
* \cgalParamDescription{completes the `freeze_bound` parameter. If it is set to `true` (default value),
* frozen vertices will not move anymore in next iterations. Otherwise, at each iteration, any vertex that
* moves, unfreezes all its incident vertices.}
* \cgalParamType{`bool`}
* \cgalParamDefault{true}
* \cgalParamNEnd
* \cgalNamedParamsEnd
*
* \return a value of type `CGAL::Mesh_optimization_return_code` which is:
* <UL>
* <LI>`CGAL::TIME_LIMIT_REACHED` when the time limit is reached.
* <LI>`CGAL::MAX_ITERATION_NUMBER_REACHED` when `lloyd_optimize_mesh_3()` stops because it has performed `max_iteration_number` iterations.
* <LI>`CGAL::CONVERGENCE_REACHED` when `lloyd_optimize_mesh_3()` stops because the convergence criterion
* is achieved.
* <LI>`CGAL::ALL_VERTICES_FROZEN` when all vertices have been frozen, when the
* `do_freeze` parameter is set to true.
* <LI>`CGAL::CANT_IMPROVE_ANYMORE` when `lloyd_optimize_mesh_3()` stops because
* most vertices have been frozen, and no better convergence can be reached.
* </UL>
*
* \cgalHeading{Example}
*
*
* \code{.cpp}
* // 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,
* parameters::do_freeze=true);
*
* \endcode
*
* \sa `CGAL::Mesh_optimization_return_code`
* \sa `CGAL::make_mesh_3()`
* \sa `CGAL::refine_mesh_3()`
* \sa `CGAL::exude_mesh_3()`
* \sa `CGAL::perturb_mesh_3()`
* \sa `CGAL::odt_optimize_mesh_3()`
*
* \note This function requires the \ref thirdpartyEigen library.
*/
template<typename C3T3, typename MeshDomain, typename CGAL_NP_TEMPLATE_PARAMETERS>
Mesh_optimization_return_code lloyd_optimize_mesh_3(C3T3& c3t3, MeshDomain& domain,const CGAL_NP_CLASS& np = parameters::default_values())
{