PolyFEM
|
#include <Electrostatics.hpp>
Public Member Functions | |
Electrostatics () | |
std::string | name () const override |
std::map< std::string, ParamFunc > | parameters () const override |
void | add_multimaterial (const int index, const json ¶ms, const Units &units) override |
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > | assemble (const LinearAssemblerData &data) const override |
computes local stiffness matrix (1x1) for bases i,j where i,j is passed in through data ie integral of grad(phi_i) dot grad(phi_j) | |
void | compute_stress_grad_multiply_mat (const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override |
void | compute_stiffness_value (const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const override |
double | compute_stored_energy (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 Eigen::MatrixXd &solution) |
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 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 | |
![]() | |
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 | ~Assembler ()=default |
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_dstress_dmu_dlambda (const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) 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 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 | 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 |
Private Attributes | |
GenericMatParam | epsilon_ |
Additional Inherited Members | |
![]() | |
typedef std::pair< std::string, Eigen::MatrixXd > | NamedMatrix |
typedef std::function< double(const RowVectorNd &, const RowVectorNd &, double, int)> | ParamFunc |
![]() | |
int | size_ = -1 |
Definition at line 12 of file Electrostatics.hpp.
|
inline |
Definition at line 15 of file Electrostatics.hpp.
|
overridevirtual |
Reimplemented from polyfem::assembler::Assembler.
Definition at line 13 of file Electrostatics.cpp.
References polyfem::assembler::GenericMatParam::add_multimaterial(), epsilon_, polyfem::Units::permittivity(), and polyfem::assembler::Assembler::size().
|
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.
|
overridevirtual |
computes local stiffness matrix (1x1) for bases i,j where i,j is passed in through data ie integral of grad(phi_i) dot grad(phi_j)
Implements polyfem::assembler::LinearAssembler.
Definition at line 20 of file Electrostatics.cpp.
References polyfem::assembler::ElementAssemblyValues::basis_values, polyfem::assembler::LinearAssemblerData::da, polyfem::assembler::ElementAssemblyValues::element_id, epsilon_, polyfem::assembler::LinearAssemblerData::i, polyfem::assembler::LinearAssemblerData::j, polyfem::assembler::LinearAssemblerData::t, polyfem::assembler::ElementAssemblyValues::val, and polyfem::assembler::LinearAssemblerData::vals.
Referenced by compute_stored_energy().
|
virtual |
local assembly function that defines the bilinear form (LHS) computes and returns a single local stiffness value
Implements polyfem::assembler::LinearAssembler.
|
overridevirtual |
Reimplemented from polyfem::assembler::Assembler.
Definition at line 46 of file Electrostatics.cpp.
double polyfem::assembler::Electrostatics::compute_stored_energy | ( | 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 Eigen::MatrixXd & | solution | ||
) |
Definition at line 78 of file Electrostatics.cpp.
References assemble(), and cache.
|
overridevirtual |
Reimplemented from polyfem::assembler::Assembler.
Definition at line 37 of file Electrostatics.cpp.
References polyfem::assembler::OptAssemblerData::grad_u_i.
|
inlineoverridevirtual |
Implements polyfem::assembler::Assembler.
Definition at line 19 of file Electrostatics.hpp.
|
overridevirtual |
Implements polyfem::assembler::Assembler.
Definition at line 67 of file Electrostatics.cpp.
References epsilon_.
|
private |
Definition at line 51 of file Electrostatics.hpp.
Referenced by add_multimaterial(), assemble(), and parameters().