PolyFEM
Loading...
Searching...
No Matches
polyfem::mesh::Remesher Class Referenceabstract

#include <Remesher.hpp>

Inheritance diagram for polyfem::mesh::Remesher:
[legend]
Collaboration diagram for polyfem::mesh::Remesher:
[legend]

Classes

struct  GlobalProjectionCache
 

Public Types

template<typename T >
using EdgeMap = std::unordered_map< std::array< size_t, 2 >, T, polyfem::utils::HashUnorderedArray< size_t, 2 >, polyfem::utils::EqualUnorderedArray< size_t, 2 > >
 Map from a (sorted) edge to an integer (ID)
 
template<typename T >
using FaceMap = std::unordered_map< std::array< size_t, 3 >, T, polyfem::utils::HashUnorderedArray< size_t, 3 >, polyfem::utils::EqualUnorderedArray< size_t, 3 > >
 Map from a (sorted) edge to an integer (ID)
 
template<typename T >
using TetMap = std::unordered_map< std::array< size_t, 4 >, T, polyfem::utils::HashUnorderedArray< size_t, 4 >, polyfem::utils::EqualUnorderedArray< size_t, 4 > >
 
template<typename T >
using BoundaryMap = std::variant< EdgeMap< T >, FaceMap< T > >
 

Public Member Functions

 Remesher (const State &state, const Eigen::MatrixXd &obstacle_displacements, const Eigen::MatrixXd &obstacle_quantities, const double current_time, const double starting_energy)
 Construct a new Remesher object.
 
virtual ~Remesher ()=default
 
virtual void init (const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXd &positions, const Eigen::MatrixXi &elements, const Eigen::MatrixXd &projection_quantities, const BoundaryMap< int > &boundary_to_id, const std::vector< int > &body_ids, const EdgeMap< double > &elastic_energy, const EdgeMap< double > &contact_energy)
 Initialize the mesh.
 
virtual bool execute ()=0
 Execute the remeshing.
 
virtual int dim () const =0
 Dimension of the mesh.
 
virtual bool is_volume () const
 Is the mesh a volumetric mesh.
 
virtual Eigen::MatrixXd rest_positions () const =0
 Exports rest positions of the stored mesh.
 
virtual Eigen::MatrixXd displacements () const =0
 Exports positions of the stored mesh.
 
virtual Eigen::MatrixXd positions () const =0
 Exports displacements of the stored mesh.
 
virtual Eigen::MatrixXi edges () const =0
 Exports edges of the stored mesh.
 
virtual Eigen::MatrixXi elements () const =0
 Exports elements of the stored mesh.
 
virtual Eigen::MatrixXi boundary_edges () const =0
 Exports boundary edges of the stored mesh.
 
virtual Eigen::MatrixXi boundary_faces () const =0
 Exports boundary faces of the stored mesh.
 
virtual Eigen::MatrixXd projection_quantities () const =0
 Exports projected quantities of the stored mesh.
 
virtual BoundaryMap< int > boundary_ids () const =0
 Exports boundary ids of the stored mesh.
 
virtual std::vector< int > body_ids () const =0
 Exports body ids of the stored mesh.
 
virtual std::vector< int > boundary_nodes (const Eigen::VectorXi &vertex_to_basis) const =0
 Get the boundary nodes of the stored mesh.
 
virtual int n_quantities () const =0
 Number of projection quantities (not including the position)
 
const Obstacleobstacle () const
 Get a reference to the collision obstacles.
 
const Eigen::MatrixXd & obstacle_displacements () const
 Get a reference to the collision obstacles' displacements.
 
const Eigen::MatrixXd & obstacle_quantities () const
 Get a reference to the collision obstacles' extra quantities.
 
virtual void set_rest_positions (const Eigen::MatrixXd &positions)=0
 Set rest positions of the stored mesh.
 
virtual void set_positions (const Eigen::MatrixXd &positions)=0
 Set deformed positions of the stored mesh.
 
virtual void set_projection_quantities (const Eigen::MatrixXd &projection_quantities)=0
 Set projected quantities of the stored mesh.
 
virtual void set_fixed (const std::vector< bool > &fixed)=0
 Set if a vertex is fixed.
 
virtual void set_boundary_ids (const BoundaryMap< int > &boundary_to_id)=0
 Set the boundary IDs of all edges.
 
virtual void set_body_ids (const std::vector< int > &body_ids)=0
 Set the body IDs of all elements.
 
void write_mesh (const std::string &path) const
 Writes a VTU mesh file.
 
void init_assembler (const std::vector< int > &body_ids) const
 Create an assembler object.
 

Static Public Member Functions

static Eigen::MatrixXd combine_time_integrator_quantities (const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > &time_integrator)
 Combine the quantities of a time integrator into a single matrix (one column per quantity)
 
static void split_time_integrator_quantities (const Eigen::MatrixXd &quantities, const int dim, Eigen::MatrixXd &x_prevs, Eigen::MatrixXd &v_prevs, Eigen::MatrixXd &a_prevs)
 Split the quantities of a time integrator into separate vectors.
 
static int build_bases (const Mesh &mesh, const std::string &assembler_formulation, std::vector< polyfem::basis::ElementBases > &bases, std::vector< LocalBoundary > &local_boundary, Eigen::VectorXi &vertex_to_basis)
 Build bases for a given mesh (V, F)
 
static void log_timings ()
 

Public Attributes

const Statestate
 Reference to the simulation state.
 
int max_op_attempts = 1
 

Static Public Attributes

static std::unordered_map< std::string, utils::Timingtimings
 Timings for the remeshing operations.
 
static double total_time = 0
 
static size_t num_solves = 0
 
static size_t total_ndofs = 0
 

Protected Member Functions

virtual void init_attributes_and_connectivity (const size_t num_vertices, const Eigen::MatrixXi &elements)=0
 Create an internal mesh representation and associate attributes.
 
void project_quantities ()
 Update the mesh positions and other projection quantities.
 
void cache_before ()
 Cache quantities before applying an operation.
 

Protected Attributes

GlobalProjectionCache global_projection_cache
 
const json args
 Copy of remesh args.
 
const Eigen::MatrixXd m_obstacle_displacements
 Collision obstacles' displacements.
 
Eigen::MatrixXd m_obstacle_quantities
 Collision obstacles' extra quantities.
 
const double current_time
 Current time.
 
const double starting_energy
 Starting energy.
 

Detailed Description

Definition at line 19 of file Remesher.hpp.

Member Typedef Documentation

◆ BoundaryMap

template<typename T >
using polyfem::mesh::Remesher::BoundaryMap = std::variant<EdgeMap<T>, FaceMap<T> >

Definition at line 42 of file Remesher.hpp.

◆ EdgeMap

template<typename T >
using polyfem::mesh::Remesher::EdgeMap = std::unordered_map< std::array<size_t, 2>, T, polyfem::utils::HashUnorderedArray<size_t, 2>, polyfem::utils::EqualUnorderedArray<size_t, 2> >

Map from a (sorted) edge to an integer (ID)

Definition at line 26 of file Remesher.hpp.

◆ FaceMap

template<typename T >
using polyfem::mesh::Remesher::FaceMap = std::unordered_map< std::array<size_t, 3>, T, polyfem::utils::HashUnorderedArray<size_t, 3>, polyfem::utils::EqualUnorderedArray<size_t, 3> >

Map from a (sorted) edge to an integer (ID)

Definition at line 32 of file Remesher.hpp.

◆ TetMap

template<typename T >
using polyfem::mesh::Remesher::TetMap = std::unordered_map< std::array<size_t, 4>, T, polyfem::utils::HashUnorderedArray<size_t, 4>, polyfem::utils::EqualUnorderedArray<size_t, 4> >

Definition at line 37 of file Remesher.hpp.

Constructor & Destructor Documentation

◆ Remesher()

polyfem::mesh::Remesher::Remesher ( const State state,
const Eigen::MatrixXd &  obstacle_displacements,
const Eigen::MatrixXd &  obstacle_quantities,
const double  current_time,
const double  starting_energy 
)

Construct a new Remesher object.

Parameters
stateSimulation current state

Definition at line 21 of file Remesher.cpp.

◆ ~Remesher()

virtual polyfem::mesh::Remesher::~Remesher ( )
virtualdefault

Member Function Documentation

◆ body_ids()

virtual std::vector< int > polyfem::mesh::Remesher::body_ids ( ) const
pure virtual

Exports body ids of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ boundary_edges()

virtual Eigen::MatrixXi polyfem::mesh::Remesher::boundary_edges ( ) const
pure virtual

Exports boundary edges of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by project_quantities().

Here is the caller graph for this function:

◆ boundary_faces()

virtual Eigen::MatrixXi polyfem::mesh::Remesher::boundary_faces ( ) const
pure virtual

Exports boundary faces of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by project_quantities().

Here is the caller graph for this function:

◆ boundary_ids()

virtual BoundaryMap< int > polyfem::mesh::Remesher::boundary_ids ( ) const
pure virtual

◆ boundary_nodes()

virtual std::vector< int > polyfem::mesh::Remesher::boundary_nodes ( const Eigen::VectorXi &  vertex_to_basis) const
pure virtual

Get the boundary nodes of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by project_quantities().

Here is the caller graph for this function:

◆ build_bases()

int polyfem::mesh::Remesher::build_bases ( const Mesh mesh,
const std::string &  assembler_formulation,
std::vector< polyfem::basis::ElementBases > &  bases,
std::vector< LocalBoundary > &  local_boundary,
Eigen::VectorXi &  vertex_to_basis 
)
static

Build bases for a given mesh (V, F)

Parameters
VMatrix of vertex (rest) positions
FMatrix of elements indices
basesOutput element bases
vertex_to_basisMap from vertex to reordered nodes
Returns
Number of bases

Definition at line 254 of file Remesher.cpp.

References polyfem::basis::LagrangeBasis2d::build_bases(), polyfem::basis::LagrangeBasis3d::build_bases(), polyfem::mesh::Mesh::dimension(), polyfem::mesh::Mesh::n_vertices(), and polyfem::mesh::Mesh::point().

Referenced by polyfem::mesh::LocalMesh< M >::build_bases(), polyfem::mesh::LocalRelaxationData< M >::init_bases(), and project_quantities().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cache_before()

void polyfem::mesh::Remesher::cache_before ( )
protected

Cache quantities before applying an operation.

Definition at line 87 of file Remesher.cpp.

References elements(), polyfem::mesh::Remesher::GlobalProjectionCache::elements, global_projection_cache, projection_quantities(), polyfem::mesh::Remesher::GlobalProjectionCache::projection_quantities, rest_positions(), and polyfem::mesh::Remesher::GlobalProjectionCache::rest_positions.

Here is the call graph for this function:

◆ combine_time_integrator_quantities()

Eigen::MatrixXd polyfem::mesh::Remesher::combine_time_integrator_quantities ( const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > &  time_integrator)
static

Combine the quantities of a time integrator into a single matrix (one column per quantity)

Definition at line 375 of file Remesher.cpp.

References projection_quantities(), and x.

Here is the call graph for this function:

◆ dim()

virtual int polyfem::mesh::Remesher::dim ( ) const
pure virtual

Dimension of the mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by is_volume(), project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ displacements()

virtual Eigen::MatrixXd polyfem::mesh::Remesher::displacements ( ) const
pure virtual

Exports positions of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by write_mesh().

Here is the caller graph for this function:

◆ edges()

virtual Eigen::MatrixXi polyfem::mesh::Remesher::edges ( ) const
pure virtual

◆ elements()

virtual Eigen::MatrixXi polyfem::mesh::Remesher::elements ( ) const
pure virtual

Exports elements of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by cache_before(), init(), project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ execute()

virtual bool polyfem::mesh::Remesher::execute ( )
pure virtual

Execute the remeshing.

Returns
True if any operation was performed.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

◆ init()

void polyfem::mesh::Remesher::init ( const Eigen::MatrixXd &  rest_positions,
const Eigen::MatrixXd &  positions,
const Eigen::MatrixXi &  elements,
const Eigen::MatrixXd &  projection_quantities,
const BoundaryMap< int > &  boundary_to_id,
const std::vector< int > &  body_ids,
const EdgeMap< double > &  elastic_energy,
const EdgeMap< double > &  contact_energy 
)
virtual

Initialize the mesh.

Parameters
rest_positionsRest positions of the mesh (|V| × dim)
positionsCurrent positions of the mesh (|V| × dim)
elementsElements of the mesh (|T| × (dim + 1))
projection_quantitiesQuantities to be projected to the new mesh (dim rows per vertex and 1 column per quantity)
edge_to_boundary_idMap from edge to boundary id (of size |E|)
body_idsBody ids of the mesh (of size |T|)

Reimplemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Definition at line 35 of file Remesher.cpp.

References body_ids(), elements(), init_attributes_and_connectivity(), positions(), projection_quantities(), rest_positions(), set_body_ids(), set_boundary_ids(), set_fixed(), set_positions(), set_projection_quantities(), and set_rest_positions().

Referenced by polyfem::mesh::WildRemesher< WMTKMesh >::init().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_assembler()

void polyfem::mesh::Remesher::init_assembler ( const std::vector< int > &  body_ids) const

Create an assembler object.

Parameters
body_idsOne body ID per element.
Returns
Assembler object

◆ init_attributes_and_connectivity()

virtual void polyfem::mesh::Remesher::init_attributes_and_connectivity ( const size_t  num_vertices,
const Eigen::MatrixXi &  elements 
)
protectedpure virtual

Create an internal mesh representation and associate attributes.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ is_volume()

virtual bool polyfem::mesh::Remesher::is_volume ( ) const
inlinevirtual

Is the mesh a volumetric mesh.

Note
Assumes non-volumetric meshes are 2D

Definition at line 102 of file Remesher.hpp.

References dim().

Referenced by project_quantities(), and polyfem::mesh::PhysicsRemesher< wmtk::TetMesh >::smooth_after().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ log_timings()

void polyfem::mesh::Remesher::log_timings ( )
static

Definition at line 414 of file Remesher.cpp.

References polyfem::logger(), num_solves, timings, total_ndofs, and total_time.

Here is the call graph for this function:

◆ n_quantities()

virtual int polyfem::mesh::Remesher::n_quantities ( ) const
pure virtual

Number of projection quantities (not including the position)

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by project_quantities().

Here is the caller graph for this function:

◆ obstacle()

const Obstacle & polyfem::mesh::Remesher::obstacle ( ) const
inline

Get a reference to the collision obstacles.

Definition at line 131 of file Remesher.hpp.

References polyfem::State::obstacle, and state.

Referenced by project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ obstacle_displacements()

const Eigen::MatrixXd & polyfem::mesh::Remesher::obstacle_displacements ( ) const
inline

Get a reference to the collision obstacles' displacements.

Definition at line 133 of file Remesher.hpp.

References m_obstacle_displacements.

Referenced by write_mesh().

Here is the caller graph for this function:

◆ obstacle_quantities()

const Eigen::MatrixXd & polyfem::mesh::Remesher::obstacle_quantities ( ) const
inline

Get a reference to the collision obstacles' extra quantities.

Definition at line 135 of file Remesher.hpp.

References m_obstacle_quantities.

Referenced by project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ positions()

virtual Eigen::MatrixXd polyfem::mesh::Remesher::positions ( ) const
pure virtual

Exports displacements of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ project_quantities()

◆ projection_quantities()

virtual Eigen::MatrixXd polyfem::mesh::Remesher::projection_quantities ( ) const
pure virtual

Exports projected quantities of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by cache_before(), combine_time_integrator_quantities(), init(), project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ rest_positions()

virtual Eigen::MatrixXd polyfem::mesh::Remesher::rest_positions ( ) const
pure virtual

Exports rest positions of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by cache_before(), init(), project_quantities(), and write_mesh().

Here is the caller graph for this function:

◆ set_body_ids()

virtual void polyfem::mesh::Remesher::set_body_ids ( const std::vector< int > &  body_ids)
pure virtual

Set the body IDs of all elements.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ set_boundary_ids()

virtual void polyfem::mesh::Remesher::set_boundary_ids ( const BoundaryMap< int > &  boundary_to_id)
pure virtual

Set the boundary IDs of all edges.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ set_fixed()

virtual void polyfem::mesh::Remesher::set_fixed ( const std::vector< bool > &  fixed)
pure virtual

Set if a vertex is fixed.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ set_positions()

virtual void polyfem::mesh::Remesher::set_positions ( const Eigen::MatrixXd &  positions)
pure virtual

Set deformed positions of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ set_projection_quantities()

virtual void polyfem::mesh::Remesher::set_projection_quantities ( const Eigen::MatrixXd &  projection_quantities)
pure virtual

Set projected quantities of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init(), and project_quantities().

Here is the caller graph for this function:

◆ set_rest_positions()

virtual void polyfem::mesh::Remesher::set_rest_positions ( const Eigen::MatrixXd &  positions)
pure virtual

Set rest positions of the stored mesh.

Implemented in polyfem::mesh::WildRemesher< WMTKMesh >, polyfem::mesh::WildRemesher< wmtk::TetMesh >, and polyfem::mesh::WildRemesher< wmtk::TriMesh >.

Referenced by init().

Here is the caller graph for this function:

◆ split_time_integrator_quantities()

void polyfem::mesh::Remesher::split_time_integrator_quantities ( const Eigen::MatrixXd &  quantities,
const int  dim,
Eigen::MatrixXd &  x_prevs,
Eigen::MatrixXd &  v_prevs,
Eigen::MatrixXd &  a_prevs 
)
static

Split the quantities of a time integrator into separate vectors.

Definition at line 396 of file Remesher.cpp.

Referenced by polyfem::mesh::LocalRelaxationData< M >::init_solve_data().

Here is the caller graph for this function:

◆ write_mesh()

void polyfem::mesh::Remesher::write_mesh ( const std::string &  path) const

Writes a VTU mesh file.

Parameters
pathOutput path

Definition at line 314 of file Remesher.cpp.

References polyfem::utils::append_rows(), dim(), displacements(), polyfem::mesh::Obstacle::e(), elements(), polyfem::utils::StringUtils::endswith(), polyfem::mesh::Obstacle::f(), polyfem::mesh::Obstacle::n_edges(), obstacle(), obstacle_displacements(), obstacle_quantities(), projection_quantities(), rest_positions(), and polyfem::utils::unflatten().

Here is the call graph for this function:

Member Data Documentation

◆ args

const json polyfem::mesh::Remesher::args
protected

Copy of remesh args.

Definition at line 213 of file Remesher.hpp.

Referenced by polyfem::mesh::PhysicsTriRemesher::swap_edge_before().

◆ current_time

const double polyfem::mesh::Remesher::current_time
protected

Current time.

Definition at line 219 of file Remesher.hpp.

◆ global_projection_cache

GlobalProjectionCache polyfem::mesh::Remesher::global_projection_cache
protected

Definition at line 210 of file Remesher.hpp.

Referenced by cache_before(), and project_quantities().

◆ m_obstacle_displacements

const Eigen::MatrixXd polyfem::mesh::Remesher::m_obstacle_displacements
protected

Collision obstacles' displacements.

Definition at line 215 of file Remesher.hpp.

Referenced by obstacle_displacements().

◆ m_obstacle_quantities

Eigen::MatrixXd polyfem::mesh::Remesher::m_obstacle_quantities
protected

Collision obstacles' extra quantities.

Definition at line 217 of file Remesher.hpp.

Referenced by obstacle_quantities().

◆ max_op_attempts

int polyfem::mesh::Remesher::max_op_attempts = 1

Definition at line 196 of file Remesher.hpp.

◆ num_solves

size_t polyfem::mesh::Remesher::num_solves = 0
static

Definition at line 231 of file Remesher.hpp.

Referenced by log_timings().

◆ starting_energy

const double polyfem::mesh::Remesher::starting_energy
protected

Starting energy.

Definition at line 221 of file Remesher.hpp.

◆ state

const State& polyfem::mesh::Remesher::state

Reference to the simulation state.

Definition at line 191 of file Remesher.hpp.

Referenced by obstacle(), and project_quantities().

◆ timings

decltype(Remesher::timings) polyfem::mesh::Remesher::timings
static

Timings for the remeshing operations.

Definition at line 229 of file Remesher.hpp.

Referenced by log_timings().

◆ total_ndofs

size_t polyfem::mesh::Remesher::total_ndofs = 0
static

Definition at line 232 of file Remesher.hpp.

Referenced by log_timings().

◆ total_time

double polyfem::mesh::Remesher::total_time = 0
static

Definition at line 230 of file Remesher.hpp.

Referenced by log_timings().


The documentation for this class was generated from the following files: