New API support and doc for odt_optimize_mesh_3.h

This commit is contained in:
GYuvanShankar 2022-06-15 13:31:45 +05:30
parent 11341ed33d
commit 89c07e7718
9 changed files with 259 additions and 49 deletions

View File

@ -7,7 +7,8 @@ INPUT += \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/Polyhedral_complex_mesh_domain_3.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/Mesh_domain_with_polyline_features_3.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/Mesh_3/generate_label_weights.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/exude_mesh_3.h
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/exude_mesh_3.h \
${CGAL_PACKAGE_INCLUDE_DIR}/CGAL/odt_optimize_mesh_3.h
PROJECT_NAME = "CGAL ${CGAL_DOC_VERSION} - 3D Mesh Generation"
HTML_EXTRA_FILES = ${CGAL_PACKAGE_DOC_DIR}/fig/implicit_domain_3.jpg \
${CGAL_PACKAGE_DOC_DIR}/fig/implicit_domain_4.jpg \

View File

@ -59,7 +59,7 @@ to represent the final optimized mesh.
\cgalParamDescription{is used to set up, in seconds, a CPU time limit after which the optimization process
is stopped. This time is measured using the Real_timer class. The default value is
0 and means that there is no time limit.}
\cgalParamType{'double'}
\cgalParamType{`double`}
\cgalParamDefault{0}
\cgalParamNBegin{sliver_bound_new}
\cgalParamDescription{is designed to give, in degrees, a targeted lower bound on dihedral angles of mesh cells.
@ -68,7 +68,7 @@ to represent the final optimized mesh.
The optimization process stops when every cell in the mesh achieves this quality. The
default value is 0 and means that there is no targeted bound: the exuder then runs as long
as it can improve the smallest dihedral angles of the set of cells incident to some vertices.}
\cgalParamType{'double'}
\cgalParamType{`double`}
\cgalParamDefault{0}
\cgalNamedParamsEnd
\return

View File

@ -21,7 +21,7 @@
#include <CGAL/disable_warnings.h>
#include <CGAL/boost/parameter.h>
#include <CGAL/Named_function_parameters.h>
#include <CGAL/Mesh_3/Mesh_global_optimizer.h>
#include <CGAL/Mesh_3/Odt_move.h>
#include <CGAL/Mesh_3/Mesh_sizing_field.h>
@ -29,43 +29,136 @@
#include <CGAL/Mesh_3/parameters_defaults.h>
#include <CGAL/Mesh_3/internal/check_weights.h>
#include <boost/parameter/preprocessor.hpp>
namespace CGAL {
/*!
@ingroup PkgMesh3Functions
#if defined(BOOST_MSVC)
# pragma warning(push)
# pragma warning(disable:4003) // not enough actual parameters for macro
#endif
The function `odt_optimize_mesh_3()` is a mesh optimization process
based on the minimization of a global energy function.
// see <CGAL/config.h>
CGAL_PRAGMA_DIAG_PUSH
// see <CGAL/boost/parameter.h>
CGAL_IGNORE_BOOST_PARAMETER_NAME_WARNINGS
In `odt_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 mesh cell.
BOOST_PARAMETER_FUNCTION(
(Mesh_optimization_return_code),
odt_optimize_mesh_3,
parameters::tag,
(required (in_out(c3t3),*) (domain,*) )
(optional
(time_limit_, *, 0 )
(max_iteration_number_, *, 0 )
(convergence_, *, parameters::default_values_for_mesh_3::odt_convergence_ratio )
(freeze_bound_, *, parameters::default_values_for_mesh_3::odt_freeze_ratio )
(do_freeze_, *, parameters::default_values_for_mesh_3::do_freeze ))
)
The optimizer `odt_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 iis 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 MeshDomain 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 ...
@param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below:
\cgalNamedParamsBegin
\cgalParamNBegin{time_limit_new}
\cgalParamDescription{is used to set up, in seconds,
a CPU time limit after which the optimization process is stopped. This time is
measured using `Real_timer`.
The default value is 0 and means that there is no time limit.}
\cgalParamType{`double`}
\cgalParamDefault{0}
\cgalParamNBegin{max_iteration_number_new}
\cgalParamDescription{sets a limit on the number of performed iterations.
The default value of 0 means that there is
no limit on the number of performed iterations.}
\cgalParamType{`std::size_t`}
\cgalParamDefault{0}
\cgalParamNBegin{convergence_new}
\cgalParamDescription{is a 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 percentage of the length
the shortest edge incident to that vertex.
The parameter `convergence` gives the threshold ratio.}
\cgalParamType{`double`}
\cgalParamDefault{0.02}
\cgalParamNBegin{freeze_bound_new}
\cgalParamDescription{is designed to reduce running time of each optimization iteration. Any vertex
that has a displacement less than a given percentage of the length of its shortest incident edge, is frozen (i.e.\ is
not relocated). The parameter `freeze_bound` gives the threshold ratio.}
\cgalParamType{`double`}
\cgalParamDefault{0.01}
\cgalParamNBegin{do_freeze_new}
\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 `odt_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 `odt_optimize_mesh_3()` stops because it has performed `max_iteration_number` iterations.
<LI>`CGAL::CONVERGENCE_REACHED` when `odt_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 `odt_optimize_mesh_3()` stops because
most vertices have been frozen, and no better convergence can be reached.
</UL>
\cgalHeading{Example}
\code{.cpp}
// 100 iterations of ODT-smoothing
odt_optimize_mesh_3(c3t3,
domain,
parameters::max_iteration_number = 100,
parameters::convergence = 0);
\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::lloyd_optimize_mesh_3()`
*/
template<typename C3T3,typename MeshDomain,typename CGAL_NP_TEMPLATE_PARAMETERS>
Mesh_optimization_return_code odt_optimize_mesh_3(C3T3& c3t3, MeshDomain& domain, const CGAL_NP_CLASS& np = parameters::default_values())
{
return odt_optimize_mesh_3_impl(c3t3, domain,
time_limit_, max_iteration_number_,
convergence_, freeze_bound_
, do_freeze_ );
using parameters::choose_parameter;
using parameters::get_parameter;
double time_limit=choose_parameter(get_parameter(np,internal_np::maximum_running_time),0);
std::size_t max_iteration_number=choose_parameter(get_parameter(np,internal_np::number_of_iterations),0);
double convergence=choose_parameter(get_parameter(np,internal_np::convergence_ratio),0.02);
double freeze_bound=choose_parameter(get_parameter(np,internal_np::vertex_freeze_bound),0.01);
bool do_freeze=choose_parameter(get_parameter(np,internal_np::freeze),true);
return odt_optimize_mesh_3_impl(c3t3, domain, time_limit, max_iteration_number, convergence, freeze_bound, do_freeze);
}
CGAL_PRAGMA_DIAG_POP
#if defined(BOOST_MSVC)
# pragma warning(pop)
#endif
#ifndef DOXYGEN_RUNNING
#ifndef CGAL_NO_DEPRECATED_CODE
template<typename C3T3, typename MeshDomain, typename... NP_Pack>
Mesh_optimization_return_code odt_optimize_mesh_3(C3T3& c3t3, MeshDomain& domain, const NP_Pack& ...nps)
{
return odt_optimize_mesh_3(c3t3, domain, internal_np::combine_named_parameters(nps...));
}
#endif //CGAL_NO_DEPRECATED_CODE
template <typename C3T3, typename MeshDomain>
Mesh_optimization_return_code
@ -106,8 +199,122 @@ odt_optimize_mesh_3_impl(C3T3& c3t3,
return opt(static_cast<int>(max_iteration_number));
}
#else
namespace CGAL {
/*!
\ingroup PkgMesh3Functions
\deprecated This function is deprecated since \cgal 5.5, the overload using `NamedParameters` must be used instead.
The function `odt_optimize_mesh_3()` is a mesh optimization process
based on the minimization of a global energy function.
In `odt_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 mesh cell.
The optimizer `odt_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.
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).
\cgalHeading{Named Parameters}
- <b>`parameters::time_limit`</b>
is used to set up, in seconds,
a CPU time limit after which the optimization process is stopped. This time is
measured using `Real_timer`.
The default value is 0 and means that there is no time limit.
- <b>`parameters::max_iteration_number`</b> sets a limit on the
number of performed iterations. The default value of 0 means that there is
no limit on the number of performed iterations.
- <b>`parameters::convergence`</b> is a 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 percentage of the length
the shortest edge incident to that vertex.
The parameter `convergence` gives the threshold ratio.
- <b>`parameters::freeze_bound`</b> is designed to reduce running time of each optimization iteration. Any vertex
that has a displacement less than a given percentage of the length (the of its shortest incident edge, is frozen (i.e.\ is
not relocated). The parameter `freeze_bound` gives the threshold ratio.
- <b>`parameters::do_freeze`</b> 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.
\return
The function `odt_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 `odt_optimize_mesh_3()` stops because it has performed `max_iteration_number` iterations.
<LI>`CGAL::CONVERGENCE_REACHED` when `odt_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 `odt_optimize_mesh_3()` stops because
most vertices have been frozen, and no better convergence can be reached.
</UL>
\cgalHeading{Example}
\code{.cpp}
// 100 iterations of ODT-smoothing
odt_optimize_mesh_3(c3t3,
domain,
parameters::max_iteration_number = 100,
parameters::convergence = 0);
\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::lloyd_optimize_mesh_3()`
*/
template<typename C3T3, typename MD>
Mesh_optimization_return_code
odt_optimize_mesh_3(C3T3& c3t3,
const MD& 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,
bool parameters::do_freeze=true);
} /* namespace CGAL */
#endif //DOXYGEN_RUNNING
} // end namespace CGAL
#include <CGAL/enable_warnings.h>
#endif // CGAL_ODT_OPTIMIZE_MESH_3_H

View File

@ -567,10 +567,10 @@ void refine_mesh_3_impl(C3T3& c3t3,
{
odt_optimize_mesh_3(c3t3,
domain,
parameters::time_limit = odt.time_limit(),
parameters::max_iteration_number = odt.max_iteration_number(),
parameters::convergence = odt.convergence(),
parameters::freeze_bound = odt.bound());
parameters::time_limit_new = odt.time_limit(),
parameters::max_iteration_number_new = odt.max_iteration_number(),
parameters::convergence_new = odt.convergence(),
parameters::freeze_bound_new = odt.bound());
}
// Lloyd

View File

@ -98,7 +98,7 @@ void test()
oss.clear();
//ODT (2)
CGAL::odt_optimize_mesh_3(c3t3, domain, max_iteration_number = nb_odt);
CGAL::odt_optimize_mesh_3(c3t3, domain, max_iteration_number_new = nb_odt);
c3t3.output_to_medit(oss);
output_c3t3.push_back(oss.str());//[i*5+2]
oss.clear();
@ -120,7 +120,7 @@ void test()
oss.clear();
//EXUDE (4)
CGAL::exude_mesh_3(c3t3, CGAL::parameters::sliver_bound_new=exude_bound);
CGAL::exude_mesh_3(c3t3, sliver_bound_new=exude_bound);
c3t3.output_to_medit(oss);
output_c3t3.push_back(oss.str());//[i*5+4]
oss.clear();

View File

@ -174,8 +174,8 @@ struct Tester
// Vertex number should not change (obvious)
C3t3 odt_c3t3(c3t3);
std::cerr << "Odt...\n";
CGAL::odt_optimize_mesh_3(odt_c3t3, domain, CGAL::parameters::time_limit=5,
CGAL::parameters::convergence=0.001, CGAL::parameters::freeze_bound=0.0005);
CGAL::odt_optimize_mesh_3(odt_c3t3, domain, CGAL::parameters::time_limit_new=5,
CGAL::parameters::convergence_new=0.001, CGAL::parameters::freeze_bound_new=0.0005);
verify_c3t3(odt_c3t3,domain,domain_type,v,v);
verify_c3t3_volume(odt_c3t3, volume*0.95, volume*1.05);
verify_c3t3_hausdorff_distance(odt_c3t3, domain, domain_type, hdist);

View File

@ -254,10 +254,10 @@ void refine_periodic_3_mesh_3_impl(C3T3& c3t3,
if(odt)
{
odt_optimize_mesh_3(c3t3, domain,
parameters::time_limit = odt.time_limit(),
parameters::max_iteration_number = odt.max_iteration_number(),
parameters::convergence = odt.convergence(),
parameters::freeze_bound = odt.bound());
parameters::time_limit_new = odt.time_limit(),
parameters::max_iteration_number_new = odt.max_iteration_number(),
parameters::convergence_new = odt.convergence(),
parameters::freeze_bound_new = odt.bound());
}
// Lloyd

View File

@ -342,6 +342,7 @@ const Boost_parameter_compatibility_wrapper<internal_np::maximum_running_time_t>
const Boost_parameter_compatibility_wrapper<internal_np::i_seed_begin_iterator_t> seeds_begin_new;
const Boost_parameter_compatibility_wrapper<internal_np::i_seed_end_iterator_t> seeds_end_new;
const Boost_parameter_compatibility_wrapper<internal_np::seeds_are_in_domain_t> mark_new;
const Boost_parameter_compatibility_wrapper<internal_np::freeze_t> do_freeze_new;
//Compatibility wrappers for exude_mesh_3.h
const Boost_parameter_compatibility_wrapper<internal_np::lower_sliver_bound_t> sliver_bound_new;
#endif

View File

@ -250,4 +250,5 @@ CGAL_add_named_parameter(i_seed_begin_iterator_t, i_seed_begin_iterator, i_seed_
CGAL_add_named_parameter(i_seed_end_iterator_t, i_seed_end_iterator, i_seed_end_iterator)
//List of named parameters used in exude_mesh_3.h
CGAL_add_named_parameter(lower_sliver_bound_t,lower_sliver_bound,lower_sliver_bound)
CGAL_add_named_parameter(lower_sliver_bound_t,lower_sliver_bound,lower_sliver_bound)
CGAL_add_named_parameter(freeze_t,freeze,freeze)