PolyFEM
|
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc More...
#include <Assembler.hpp>
Public Member Functions | |
LinearAssembler () | |
virtual | ~LinearAssembler ()=default |
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 override |
assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by the overloaded assemble (see below) function that the subclass (eg Laplacian) defines sets stiffness and modifies cache if it has not already been computed | |
virtual bool | is_linear () const override |
virtual Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > | assemble (const LinearAssemblerData &data) const =0 |
local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value | |
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 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 |
virtual 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 |
virtual 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 |
virtual 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 |
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, ParamFunc > | parameters () 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 ¶ms, 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 |
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 |
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc
Definition at line 208 of file Assembler.hpp.
polyfem::assembler::LinearAssembler::LinearAssembler | ( | ) |
Definition at line 153 of file Assembler.cpp.
|
virtualdefault |
|
overridevirtual |
assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by the overloaded assemble (see below) function that the subclass (eg Laplacian) defines sets stiffness and modifies cache if it has not already been computed
Reimplemented from polyfem::assembler::Assembler.
Reimplemented in polyfem::assembler::LinearElasticity, polyfem::assembler::Mass, polyfem::assembler::StokesVelocity, and polyfem::assembler::StokesPressure.
Referenced by polyfem::basis::RBFWithQuadratic::compute_constraints_matrix_2d(), polyfem::basis::RBFWithQuadraticLagrange::compute_constraints_matrix_2d(), and polyfem::basis::PolygonalBasis2d::compute_integral_constraints().
|
pure virtual |
local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value
Implemented in polyfem::assembler::BilaplacianMain, polyfem::assembler::BilaplacianAux, polyfem::assembler::Helmholtz, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::IncompressibleLinearElasticityDispacement, polyfem::assembler::IncompressibleLinearElasticityPressure, polyfem::assembler::Laplacian, polyfem::assembler::LinearElasticity, polyfem::assembler::Mass, polyfem::assembler::StokesVelocity, polyfem::assembler::StokesPressure, polyfem::assembler::BilaplacianMain, polyfem::assembler::BilaplacianAux, polyfem::assembler::Helmholtz, polyfem::assembler::HookeLinearElasticity, polyfem::assembler::IncompressibleLinearElasticityDispacement, polyfem::assembler::IncompressibleLinearElasticityPressure, polyfem::assembler::Laplacian, polyfem::assembler::LinearElasticity, polyfem::assembler::Mass, polyfem::assembler::StokesVelocity, and polyfem::assembler::StokesPressure.
|
inlineoverridevirtual |
Implements polyfem::assembler::Assembler.
Reimplemented in polyfem::assembler::HookeLinearElasticity, and polyfem::assembler::LinearElasticity.
Definition at line 230 of file Assembler.hpp.