PolyFEM
Loading...
Searching...
No Matches
polyfem::solver::ElasticForm Class Reference

Form of the elasticity potential and forces. More...

#include <ElasticForm.hpp>

Inheritance diagram for polyfem::solver::ElasticForm:
[legend]
Collaboration diagram for polyfem::solver::ElasticForm:
[legend]

Public Member Functions

 ElasticForm (const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &geom_bases, const assembler::Assembler &assembler, const assembler::AssemblyValsCache &ass_vals_cache, const double t, const double dt, const bool is_volume)
 Construct a new Elastic Form object.
 
std::string name () const override
 
bool is_step_valid (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
 Determine if a step from solution x0 to solution x1 is allowed.
 
void update_quantities (const double t, const Eigen::VectorXd &x) override
 Update time-dependent fields.
 
void force_material_derivative (const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
 Compute the derivative of the force wrt lame/damping parameters, then multiply the resulting matrix with adjoint_sol.
 
void force_shape_derivative (const double t, const int n_verts, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
 Compute the derivative of the force wrt vertex positions, then multiply the resulting matrix with adjoint_sol.
 
- Public Member Functions inherited from polyfem::solver::Form
virtual ~Form ()
 
virtual void init (const Eigen::VectorXd &x)
 Initialize the form.
 
virtual double value (const Eigen::VectorXd &x) const
 Compute the value of the form multiplied with the weigth.
 
Eigen::VectorXd value_per_element (const Eigen::VectorXd &x) const
 Compute the value of the form multiplied with the weigth.
 
virtual void first_derivative (const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
 Compute the first derivative of the value wrt x multiplied with the weigth.
 
void second_derivative (const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
 Compute the second derivative of the value wrt x multiplied with the weigth.
 
virtual double max_step_size (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
 Determine the maximum step size allowable between the current and next solution.
 
virtual void line_search_begin (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
 Initialize variables used during the line search.
 
virtual void line_search_end ()
 Clear variables used during the line search.
 
virtual void post_step (const polysolve::nonlinear::PostStepData &data)
 Update fields after a step in the optimization.
 
virtual void solution_changed (const Eigen::VectorXd &new_x)
 Update cached fields upon a change in the solution.
 
virtual void init_lagging (const Eigen::VectorXd &x)
 Initialize lagged fields TODO: more than one step.
 
virtual void update_lagging (const Eigen::VectorXd &x, const int iter_num)
 Update lagged fields.
 
virtual int max_lagging_iterations () const
 Get the maximum number of lagging iteration allowable.
 
virtual bool uses_lagging () const
 Does this form require lagging?
 
void set_project_to_psd (bool val)
 Set project to psd.
 
bool is_project_to_psd () const
 Get if the form's second derivative is projected to psd.
 
void enable ()
 Enable the form.
 
void disable ()
 Disable the form.
 
void set_enabled (const bool enabled)
 Set if the form is enabled.
 
bool enabled () const
 Determine if the form is enabled.
 
virtual double weight () const
 Get the form's multiplicative constant weight.
 
void set_weight (const double weight)
 Set the form's multiplicative constant weight.
 
virtual bool is_step_collision_free (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
 Checks if the step is collision free.
 
void set_output_dir (const std::string &output_dir)
 

Protected Member Functions

virtual double value_unweighted (const Eigen::VectorXd &x) const override
 Compute the elastic potential value.
 
Eigen::VectorXd value_per_element_unweighted (const Eigen::VectorXd &x) const override
 Compute the value of the form multiplied with the weigth.
 
virtual void first_derivative_unweighted (const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
 Compute the first derivative of the value wrt x.
 
void second_derivative_unweighted (const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
 Compute the second derivative of the value wrt x.
 
- Protected Member Functions inherited from polyfem::solver::Form
std::string resolve_output_path (const std::string &path) const
 

Private Member Functions

void compute_cached_stiffness ()
 Compute the stiffness matrix (cached)
 

Private Attributes

const int n_bases_
 
const std::vector< basis::ElementBases > & bases_
 
const std::vector< basis::ElementBases > & geom_bases_
 
const assembler::Assemblerassembler_
 Reference to the assembler.
 
const assembler::AssemblyValsCacheass_vals_cache_
 
double t_
 
const double dt_
 
const bool is_volume_
 
StiffnessMatrix cached_stiffness_
 Cached stiffness matrix for linear elasticity.
 
std::unique_ptr< utils::MatrixCachemat_cache_
 Matrix cache (mutable because it is modified in second_derivative_unweighted)
 
Eigen::VectorXd x_prev_
 

Additional Inherited Members

- Protected Attributes inherited from polyfem::solver::Form
bool project_to_psd_ = false
 If true, the form's second derivative is projected to be positive semidefinite.
 
double weight_ = 1
 weight of the form (e.g., AL penalty weight or Δt²)
 
bool enabled_ = true
 If true, the form is enabled.
 
std::string output_dir_
 

Detailed Description

Form of the elasticity potential and forces.

Definition at line 14 of file ElasticForm.hpp.

Constructor & Destructor Documentation

◆ ElasticForm()

polyfem::solver::ElasticForm::ElasticForm ( const int  n_bases,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  geom_bases,
const assembler::Assembler assembler,
const assembler::AssemblyValsCache ass_vals_cache,
const double  t,
const double  dt,
const bool  is_volume 
)

Construct a new Elastic Form object.

Parameters
stateReference to the simulation state

Definition at line 34 of file ElasticForm.cpp.

References assembler_, compute_cached_stiffness(), polyfem::assembler::Assembler::is_linear(), and mat_cache_.

Here is the call graph for this function:

Member Function Documentation

◆ compute_cached_stiffness()

void polyfem::solver::ElasticForm::compute_cached_stiffness ( )
private

Compute the stiffness matrix (cached)

Definition at line 116 of file ElasticForm.cpp.

References ass_vals_cache_, polyfem::assembler::Assembler::assemble(), assembler_, bases_, cached_stiffness_, geom_bases_, polyfem::assembler::Assembler::is_linear(), is_volume_, n_bases_, and t_.

Referenced by ElasticForm().

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

◆ first_derivative_unweighted()

void polyfem::solver::ElasticForm::first_derivative_unweighted ( const Eigen::VectorXd &  x,
Eigen::VectorXd &  gradv 
) const
overrideprotectedvirtual

Compute the first derivative of the value wrt x.

Parameters
[in]xCurrent solution
[out]gradvOutput gradient of the value wrt x

Implements polyfem::solver::Form.

Definition at line 71 of file ElasticForm.cpp.

References ass_vals_cache_, polyfem::assembler::Assembler::assemble_gradient(), assembler_, bases_, dt_, geom_bases_, is_volume_, n_bases_, t_, x, and x_prev_.

Here is the call graph for this function:

◆ force_material_derivative()

void polyfem::solver::ElasticForm::force_material_derivative ( const double  t,
const Eigen::MatrixXd &  x,
const Eigen::MatrixXd &  x_prev,
const Eigen::MatrixXd &  adjoint,
Eigen::VectorXd &  term 
)

Compute the derivative of the force wrt lame/damping parameters, then multiply the resulting matrix with adjoint_sol.

Parameters
tCurrent time
[in]xCurrent solution
[in]adjointCurrent adjoint solution
[out]termDerivative of force multiplied by the adjoint

Definition at line 125 of file ElasticForm.cpp.

References ass_vals_cache_, assembler_, bases_, polyfem::assembler::AssemblyValsCache::compute(), polyfem::assembler::Assembler::compute_dstress_dmu_dlambda(), polyfem::assembler::ViscousDamping::compute_dstress_dpsi_dphi(), polyfem::utils::create_thread_storage(), polyfem::assembler::ElementAssemblyValues::det, dt_, geom_bases_, polyfem::utils::get_local_thread_storage(), polyfem::io::Evaluator::interpolate_at_local_vals(), is_volume_, polyfem::utils::maybe_parallel_for(), polyfem::assembler::Assembler::name(), polyfem::quadrature::Quadrature::points, polyfem::assembler::ElementAssemblyValues::quadrature, quadrature, polyfem::assembler::ElementAssemblyValues::val, vals, polyfem::utils::vector2matrix(), polyfem::quadrature::Quadrature::weights, and x.

Here is the call graph for this function:

◆ force_shape_derivative()

void polyfem::solver::ElasticForm::force_shape_derivative ( const double  t,
const int  n_verts,
const Eigen::MatrixXd &  x,
const Eigen::MatrixXd &  x_prev,
const Eigen::MatrixXd &  adjoint,
Eigen::VectorXd &  term 
)

◆ is_step_valid()

bool polyfem::solver::ElasticForm::is_step_valid ( const Eigen::VectorXd &  x0,
const Eigen::VectorXd &  x1 
) const
overridevirtual

Determine if a step from solution x0 to solution x1 is allowed.

Parameters
x0Current solution
x1Proposed next solution
Returns
True if the step is allowed

Reimplemented from polyfem::solver::Form.

Definition at line 99 of file ElasticForm.cpp.

References polyfem::solver::Form::first_derivative().

Here is the call graph for this function:

◆ name()

std::string polyfem::solver::ElasticForm::name ( ) const
inlineoverridevirtual

Implements polyfem::solver::Form.

Definition at line 27 of file ElasticForm.hpp.

◆ second_derivative_unweighted()

void polyfem::solver::ElasticForm::second_derivative_unweighted ( const Eigen::VectorXd &  x,
StiffnessMatrix hessian 
) const
overrideprotectedvirtual

Compute the second derivative of the value wrt x.

Parameters
[in]xCurrent solution
[out]hessianOutput Hessian of the value wrt x

Implements polyfem::solver::Form.

Definition at line 79 of file ElasticForm.cpp.

References ass_vals_cache_, polyfem::assembler::Assembler::assemble_hessian(), assembler_, bases_, cached_stiffness_, dt_, geom_bases_, polyfem::assembler::Assembler::is_linear(), is_volume_, mat_cache_, n_bases_, POLYFEM_SCOPED_TIMER, polyfem::solver::Form::project_to_psd_, t_, x, and x_prev_.

Here is the call graph for this function:

◆ update_quantities()

void polyfem::solver::ElasticForm::update_quantities ( const double  t,
const Eigen::VectorXd &  x 
)
inlineoverridevirtual

Update time-dependent fields.

Parameters
tCurrent time
xCurrent solution at time t

Reimplemented from polyfem::solver::Form.

Definition at line 60 of file ElasticForm.hpp.

References t_, x, and x_prev_.

◆ value_per_element_unweighted()

Eigen::VectorXd polyfem::solver::ElasticForm::value_per_element_unweighted ( const Eigen::VectorXd &  x) const
overrideprotectedvirtual

Compute the value of the form multiplied with the weigth.

Parameters
xCurrent solution
Returns
Computed value

Reimplemented from polyfem::solver::Form.

Definition at line 63 of file ElasticForm.cpp.

References ass_vals_cache_, polyfem::assembler::Assembler::assemble_energy_per_element(), assembler_, bases_, dt_, geom_bases_, is_volume_, t_, value_unweighted(), x, and x_prev_.

Here is the call graph for this function:

◆ value_unweighted()

double polyfem::solver::ElasticForm::value_unweighted ( const Eigen::VectorXd &  x) const
overrideprotectedvirtual

Compute the elastic potential value.

Parameters
xCurrent solution
Returns
Value of the elastic potential

Implements polyfem::solver::Form.

Definition at line 56 of file ElasticForm.cpp.

References ass_vals_cache_, polyfem::assembler::Assembler::assemble_energy(), assembler_, bases_, dt_, geom_bases_, is_volume_, t_, x, and x_prev_.

Referenced by value_per_element_unweighted().

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

Member Data Documentation

◆ ass_vals_cache_

◆ assembler_

◆ bases_

◆ cached_stiffness_

StiffnessMatrix polyfem::solver::ElasticForm::cached_stiffness_
private

Cached stiffness matrix for linear elasticity.

Definition at line 92 of file ElasticForm.hpp.

Referenced by compute_cached_stiffness(), and second_derivative_unweighted().

◆ dt_

◆ geom_bases_

◆ is_volume_

◆ mat_cache_

std::unique_ptr<utils::MatrixCache> polyfem::solver::ElasticForm::mat_cache_
mutableprivate

Matrix cache (mutable because it is modified in second_derivative_unweighted)

Definition at line 93 of file ElasticForm.hpp.

Referenced by ElasticForm(), and second_derivative_unweighted().

◆ n_bases_

const int polyfem::solver::ElasticForm::n_bases_
private

◆ t_

◆ x_prev_

Eigen::VectorXd polyfem::solver::ElasticForm::x_prev_
private

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