PolyFEM
Loading...
Searching...
No Matches
polyfem::assembler::NLAssembler Class Referenceabstract

#include <Assembler.hpp>

Inheritance diagram for polyfem::assembler::NLAssembler:
[legend]
Collaboration diagram for polyfem::assembler::NLAssembler:
[legend]

Public Member Functions

virtual ~NLAssembler ()=default
 
double assemble_energy (const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
 
Eigen::VectorXd assemble_energy_per_element (const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
 
void assemble_gradient (const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const override
 
void assemble_hessian (const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const override
 
virtual bool is_linear () const override
 
- Public Member Functions inherited from polyfem::assembler::Assembler
virtual ~Assembler ()=default
 
virtual std::string name () const =0
 
int size () const
 
virtual void set_size (const int size)
 
virtual void assemble (const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const
 
virtual void compute_scalar_value (const OutputData &data, std::vector< NamedMatrix > &result) const
 
virtual void compute_tensor_value (const OutputData &data, std::vector< NamedMatrix > &result) const
 
virtual void compute_stiffness_value (const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const
 
virtual void compute_dstress_dmu_dlambda (const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const
 
virtual void compute_stress_grad_multiply_mat (const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad_multiply_stress (const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad_multiply_vect (const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_grad (const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
 
virtual void compute_stress_prev_grad (const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &result) const
 
virtual std::map< std::string, ParamFuncparameters () const =0
 
virtual VectorNd compute_rhs (const AutodiffHessianPt &pt) const
 
virtual Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > kernel (const int dim, const AutodiffGradPt &rvect, const AutodiffScalarGrad &r) const
 
void set_materials (const std::vector< int > &body_ids, const json &body_params, const Units &units)
 
virtual void add_multimaterial (const int index, const json &params, const Units &units)
 
virtual void update_lame_params (const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus)
 
virtual bool is_solution_displacement () const
 
virtual bool is_fluid () const
 
virtual bool is_tensor () const
 

Protected Member Functions

virtual double compute_energy (const NonLinearAssemblerData &data) const =0
 
virtual Eigen::VectorXd assemble_gradient (const NonLinearAssemblerData &data) const =0
 
virtual Eigen::MatrixXd assemble_hessian (const NonLinearAssemblerData &data) const =0
 

Additional Inherited Members

- Public Types inherited from polyfem::assembler::Assembler
typedef std::pair< std::string, Eigen::MatrixXd > NamedMatrix
 
typedef std::function< double(const RowVectorNd &, const RowVectorNd &, double, int)> ParamFunc
 
- Protected Attributes inherited from polyfem::assembler::Assembler
int size_ = -1
 

Detailed Description

Definition at line 238 of file Assembler.hpp.

Constructor & Destructor Documentation

◆ ~NLAssembler()

virtual polyfem::assembler::NLAssembler::~NLAssembler ( )
virtualdefault

Member Function Documentation

◆ assemble_energy()

double polyfem::assembler::NLAssembler::assemble_energy ( const bool  is_volume,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  gbases,
const AssemblyValsCache cache,
const double  t,
const double  dt,
const Eigen::MatrixXd &  displacement,
const Eigen::MatrixXd &  displacement_prev 
) const
overridevirtual

Reimplemented from polyfem::assembler::Assembler.

Reimplemented in polyfem::assembler::SaintVenantElasticity, polyfem::assembler::ViscousDamping, and polyfem::assembler::ViscousDampingPrev.

Definition at line 495 of file Assembler.cpp.

References cache, compute_energy(), polyfem::utils::create_thread_storage(), polyfem::utils::get_local_thread_storage(), polyfem::utils::maybe_parallel_for(), quadrature, val, and vals.

Referenced by assemble_energy_per_element().

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

◆ assemble_energy_per_element()

Eigen::VectorXd polyfem::assembler::NLAssembler::assemble_energy_per_element ( const bool  is_volume,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  gbases,
const AssemblyValsCache cache,
const double  t,
const double  dt,
const Eigen::MatrixXd &  displacement,
const Eigen::MatrixXd &  displacement_prev 
) const
overridevirtual

Reimplemented from polyfem::assembler::Assembler.

Definition at line 533 of file Assembler.cpp.

References assemble_energy(), cache, compute_energy(), polyfem::utils::create_thread_storage(), polyfem::utils::get_local_thread_storage(), polyfem::utils::maybe_parallel_for(), quadrature, val, and vals.

Here is the call graph for this function:

◆ assemble_gradient() [1/2]

void polyfem::assembler::NLAssembler::assemble_gradient ( const bool  is_volume,
const int  n_basis,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  gbases,
const AssemblyValsCache cache,
const double  t,
const double  dt,
const Eigen::MatrixXd &  displacement,
const Eigen::MatrixXd &  displacement_prev,
Eigen::MatrixXd &  rhs 
) const
overridevirtual

◆ assemble_gradient() [2/2]

virtual Eigen::VectorXd polyfem::assembler::NLAssembler::assemble_gradient ( const NonLinearAssemblerData data) const
protectedpure virtual

Implemented in polyfem::assembler::GenericElastic< AMIPSEnergyAutodiff >, polyfem::assembler::GenericElastic< IncompressibleOgdenElasticity >, polyfem::assembler::GenericElastic< MooneyRivlin3ParamElasticity >, polyfem::assembler::GenericElastic< MooneyRivlinElasticity >, polyfem::assembler::GenericElastic< NeoHookeanAutodiff >, polyfem::assembler::GenericElastic< UnconstrainedOgdenElasticity >, polyfem::assembler::AMIPSEnergy, polyfem::assembler::FixedCorotational, polyfem::assembler::GenericElastic< Derived >, polyfem::assembler::GenericElastic< AMIPSEnergyAutodiff >, polyfem::assembler::GenericElastic< IncompressibleOgdenElasticity >, polyfem::assembler::GenericElastic< MooneyRivlin3ParamElasticity >, polyfem::assembler::GenericElastic< MooneyRivlinElasticity >, polyfem::assembler::GenericElastic< NeoHookeanAutodiff >, polyfem::assembler::GenericElastic< UnconstrainedOgdenElasticity >, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::LinearElasticity, polyfem::assembler::MooneyRivlin3ParamSymbolic, polyfem::assembler::MultiModel, polyfem::assembler::NavierStokesVelocity, polyfem::assembler::NeoHookeanElasticity, polyfem::assembler::SaintVenantElasticity, polyfem::assembler::ViscousDamping, polyfem::assembler::ViscousDampingPrev, polyfem::assembler::AMIPSEnergy, polyfem::assembler::FixedCorotational, polyfem::assembler::GenericElastic< Derived >, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::LinearElasticity, polyfem::assembler::MooneyRivlin3ParamSymbolic, polyfem::assembler::MultiModel, polyfem::assembler::NavierStokesVelocity, polyfem::assembler::NeoHookeanElasticity, polyfem::assembler::SaintVenantElasticity, polyfem::assembler::ViscousDamping, and polyfem::assembler::ViscousDampingPrev.

◆ assemble_hessian() [1/2]

void polyfem::assembler::NLAssembler::assemble_hessian ( const bool  is_volume,
const int  n_basis,
const bool  project_to_psd,
const std::vector< basis::ElementBases > &  bases,
const std::vector< basis::ElementBases > &  gbases,
const AssemblyValsCache cache,
const double  t,
const double  dt,
const Eigen::MatrixXd &  displacement,
const Eigen::MatrixXd &  displacement_prev,
utils::MatrixCache mat_cache,
StiffnessMatrix grad 
) const
overridevirtual

◆ assemble_hessian() [2/2]

virtual Eigen::MatrixXd polyfem::assembler::NLAssembler::assemble_hessian ( const NonLinearAssemblerData data) const
protectedpure virtual

Implemented in polyfem::assembler::GenericElastic< AMIPSEnergyAutodiff >, polyfem::assembler::GenericElastic< IncompressibleOgdenElasticity >, polyfem::assembler::GenericElastic< MooneyRivlin3ParamElasticity >, polyfem::assembler::GenericElastic< MooneyRivlinElasticity >, polyfem::assembler::GenericElastic< NeoHookeanAutodiff >, polyfem::assembler::GenericElastic< UnconstrainedOgdenElasticity >, polyfem::assembler::AMIPSEnergy, polyfem::assembler::FixedCorotational, polyfem::assembler::GenericElastic< Derived >, polyfem::assembler::GenericElastic< AMIPSEnergyAutodiff >, polyfem::assembler::GenericElastic< IncompressibleOgdenElasticity >, polyfem::assembler::GenericElastic< MooneyRivlin3ParamElasticity >, polyfem::assembler::GenericElastic< MooneyRivlinElasticity >, polyfem::assembler::GenericElastic< NeoHookeanAutodiff >, polyfem::assembler::GenericElastic< UnconstrainedOgdenElasticity >, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::LinearElasticity, polyfem::assembler::MooneyRivlin3ParamSymbolic, polyfem::assembler::MultiModel, polyfem::assembler::NavierStokesVelocity, polyfem::assembler::NeoHookeanElasticity, polyfem::assembler::SaintVenantElasticity, polyfem::assembler::ViscousDamping, polyfem::assembler::ViscousDampingPrev, polyfem::assembler::AMIPSEnergy, polyfem::assembler::FixedCorotational, polyfem::assembler::GenericElastic< Derived >, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::LinearElasticity, polyfem::assembler::MooneyRivlin3ParamSymbolic, polyfem::assembler::MultiModel, polyfem::assembler::NavierStokesVelocity, polyfem::assembler::NeoHookeanElasticity, polyfem::assembler::SaintVenantElasticity, polyfem::assembler::ViscousDamping, and polyfem::assembler::ViscousDampingPrev.

◆ compute_energy()

◆ is_linear()

virtual bool polyfem::assembler::NLAssembler::is_linear ( ) const
inlineoverridevirtual

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