PolyFEM
|
#include <PressureAssembler.hpp>
Public Member Functions | |
PressureAssembler (const Assembler &assembler, const mesh::Mesh &mesh, const mesh::Obstacle &obstacle, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const std::unordered_map< int, std::vector< mesh::LocalBoundary > > &local_pressure_cavity, const std::vector< int > &dirichlet_nodes, const std::vector< int > &primitive_to_nodes, const std::vector< int > &node_to_primitives, const int n_basis, const int size, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Problem &problem) | |
double | compute_energy (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const int resolution, const double t) const |
void | compute_energy_grad (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, const double t, Eigen::VectorXd &grad) const |
void | compute_energy_hess (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, const double t, const bool project_to_psd, StiffnessMatrix &hess) const |
double | compute_cavity_energy (const Eigen::MatrixXd &displacement, const std::unordered_map< int, std::vector< mesh::LocalBoundary > > &local_pressure_cavity, const int resolution, const double t) const |
void | compute_cavity_energy_grad (const Eigen::MatrixXd &displacement, const std::unordered_map< int, std::vector< mesh::LocalBoundary > > &local_pressure_cavity, const std::vector< int > &dirichlet_nodes, const int resolution, const double t, Eigen::VectorXd &grad) const |
void | compute_cavity_energy_hess (const Eigen::MatrixXd &displacement, const std::unordered_map< int, std::vector< mesh::LocalBoundary > > &local_pressure_cavity, const std::vector< int > &dirichlet_nodes, const int resolution, const double t, const bool project_to_psd, StiffnessMatrix &hess) const |
void | compute_force_jacobian (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, const double t, const int n_vertices, StiffnessMatrix &hess) const |
const Problem & | problem () const |
const mesh::Mesh & | mesh () const |
const std::vector< basis::ElementBases > & | bases () const |
const std::vector< basis::ElementBases > & | gbases () const |
const Assembler & | assembler () const |
void | compute_grad_volume_id (const Eigen::MatrixXd &displacement, const int boundary_id, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, Eigen::VectorXd &grad, const double t=0, const bool multiply_pressure=false) const |
Private Member Functions | |
double | compute_volume (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_boundary, const int resolution, const double t=0, const bool multiply_pressure=false) const |
void | compute_grad_volume (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, Eigen::VectorXd &grad, const double t=0, const bool multiply_pressure=false) const |
void | compute_hess_volume_3d (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, StiffnessMatrix &hess, const double t=0, const bool multiply_pressure=false) const |
void | compute_hess_volume_2d (const Eigen::MatrixXd &displacement, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &dirichlet_nodes, const int resolution, StiffnessMatrix &hess, const double t=0, const bool multiply_pressure=false) const |
bool | is_closed_or_boundary_fixed (const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &dirichlet_nodes) const |
Private Attributes | |
const Assembler & | assembler_ |
const mesh::Mesh & | mesh_ |
const mesh::Obstacle & | obstacle_ |
const int | n_basis_ |
const int | size_ |
const std::vector< basis::ElementBases > & | bases_ |
const std::vector< basis::ElementBases > & | gbases_ |
const Problem & | problem_ |
std::unordered_map< int, double > | starting_volumes_ |
std::unique_ptr< ThermodynamicProcess > | cavity_thermodynamics_ |
const std::vector< int > | primitive_to_nodes_ |
const std::vector< int > | node_to_primitives_ |
std::set< int > | relevant_pressure_nodes_ |
Definition at line 74 of file PressureAssembler.hpp.
polyfem::assembler::PressureAssembler::PressureAssembler | ( | const Assembler & | assembler, |
const mesh::Mesh & | mesh, | ||
const mesh::Obstacle & | obstacle, | ||
const std::vector< mesh::LocalBoundary > & | local_pressure_boundary, | ||
const std::unordered_map< int, std::vector< mesh::LocalBoundary > > & | local_pressure_cavity, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const std::vector< int > & | primitive_to_nodes, | ||
const std::vector< int > & | node_to_primitives, | ||
const int | n_basis, | ||
const int | size, | ||
const std::vector< basis::ElementBases > & | bases, | ||
const std::vector< basis::ElementBases > & | gbases, | ||
const Problem & | problem | ||
) |
Definition at line 896 of file PressureAssembler.cpp.
References cavity_thermodynamics_, compute_volume(), is_closed_or_boundary_fixed(), polyfem::logger(), and starting_volumes_.
|
inline |
Definition at line 143 of file PressureAssembler.hpp.
References assembler_.
|
inline |
Definition at line 141 of file PressureAssembler.hpp.
References bases_.
double polyfem::assembler::PressureAssembler::compute_cavity_energy | ( | const Eigen::MatrixXd & | displacement, |
const std::unordered_map< int, std::vector< mesh::LocalBoundary > > & | local_pressure_cavity, | ||
const int | resolution, | ||
const double | t | ||
) | const |
Definition at line 930 of file PressureAssembler.cpp.
References cavity_thermodynamics_, compute_volume(), polyfem::assembler::Problem::pressure_cavity_bc(), problem_, and starting_volumes_.
Referenced by polyfem::solver::PressureForm::value_unweighted().
void polyfem::assembler::PressureAssembler::compute_cavity_energy_grad | ( | const Eigen::MatrixXd & | displacement, |
const std::unordered_map< int, std::vector< mesh::LocalBoundary > > & | local_pressure_cavity, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
const double | t, | ||
Eigen::VectorXd & | grad | ||
) | const |
Definition at line 948 of file PressureAssembler.cpp.
References cavity_thermodynamics_, compute_grad_volume(), compute_volume(), polyfem::assembler::Problem::pressure_cavity_bc(), problem_, and starting_volumes_.
Referenced by polyfem::solver::PressureForm::first_derivative_unweighted().
void polyfem::assembler::PressureAssembler::compute_cavity_energy_hess | ( | const Eigen::MatrixXd & | displacement, |
const std::unordered_map< int, std::vector< mesh::LocalBoundary > > & | local_pressure_cavity, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
const double | t, | ||
const bool | project_to_psd, | ||
StiffnessMatrix & | hess | ||
) | const |
Definition at line 973 of file PressureAssembler.cpp.
References cavity_thermodynamics_, compute_grad_volume(), compute_hess_volume_2d(), compute_hess_volume_3d(), compute_volume(), polyfem::log_and_throw_error(), polyfem::assembler::Problem::pressure_cavity_bc(), problem_, size_, and starting_volumes_.
Referenced by polyfem::solver::PressureForm::second_derivative_unweighted().
double polyfem::assembler::PressureAssembler::compute_energy | ( | const Eigen::MatrixXd & | displacement, |
const std::vector< mesh::LocalBoundary > & | local_pressure_boundary, | ||
const int | resolution, | ||
const double | t | ||
) | const |
Definition at line 1012 of file PressureAssembler.cpp.
References compute_volume().
Referenced by polyfem::solver::PressureForm::value_unweighted().
void polyfem::assembler::PressureAssembler::compute_energy_grad | ( | const Eigen::MatrixXd & | displacement, |
const std::vector< mesh::LocalBoundary > & | local_pressure_boundary, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
const double | t, | ||
Eigen::VectorXd & | grad | ||
) | const |
Definition at line 1021 of file PressureAssembler.cpp.
References compute_grad_volume().
Referenced by polyfem::solver::PressureForm::first_derivative_unweighted().
void polyfem::assembler::PressureAssembler::compute_energy_hess | ( | const Eigen::MatrixXd & | displacement, |
const std::vector< mesh::LocalBoundary > & | local_pressure_boundary, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
const double | t, | ||
const bool | project_to_psd, | ||
StiffnessMatrix & | hess | ||
) | const |
Definition at line 1032 of file PressureAssembler.cpp.
References compute_hess_volume_2d(), compute_hess_volume_3d(), and size_.
Referenced by compute_force_jacobian(), and polyfem::solver::PressureForm::second_derivative_unweighted().
void polyfem::assembler::PressureAssembler::compute_force_jacobian | ( | const Eigen::MatrixXd & | displacement, |
const std::vector< mesh::LocalBoundary > & | local_pressure_boundary, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
const double | t, | ||
const int | n_vertices, | ||
StiffnessMatrix & | hess | ||
) | const |
Definition at line 1047 of file PressureAssembler.cpp.
References compute_energy_hess().
Referenced by polyfem::solver::PressureForm::force_shape_derivative().
|
private |
Definition at line 210 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), polyfem::utils::maybe_parallel_for(), n_basis_, and size_.
Referenced by compute_cavity_energy_grad(), compute_cavity_energy_hess(), and compute_energy_grad().
void polyfem::assembler::PressureAssembler::compute_grad_volume_id | ( | const Eigen::MatrixXd & | displacement, |
const int | boundary_id, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary, | ||
const std::vector< int > & | dirichlet_nodes, | ||
const int | resolution, | ||
Eigen::VectorXd & | grad, | ||
const double | t = 0 , |
||
const bool | multiply_pressure = false |
||
) | const |
Definition at line 335 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), polyfem::utils::maybe_parallel_for(), n_basis_, and size_.
Referenced by polyfem::solver::PressureForm::force_pressure_derivative().
|
private |
Definition at line 652 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), polyfem::utils::maybe_parallel_for(), n_basis_, and size_.
Referenced by compute_cavity_energy_hess(), and compute_energy_hess().
|
private |
Definition at line 465 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), polyfem::utils::maybe_parallel_for(), n_basis_, and size_.
Referenced by compute_cavity_energy_hess(), and compute_energy_hess().
|
private |
Definition at line 90 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), and polyfem::utils::maybe_parallel_for().
Referenced by compute_cavity_energy(), compute_cavity_energy_grad(), compute_cavity_energy_hess(), compute_energy(), and PressureAssembler().
|
inline |
Definition at line 142 of file PressureAssembler.hpp.
References gbases_.
|
private |
Definition at line 821 of file PressureAssembler.cpp.
References polyfem::utils::create_thread_storage(), and polyfem::utils::maybe_parallel_for().
Referenced by PressureAssembler().
|
inline |
Definition at line 140 of file PressureAssembler.hpp.
References mesh_.
|
inline |
Definition at line 139 of file PressureAssembler.hpp.
References problem_.
|
private |
Definition at line 190 of file PressureAssembler.hpp.
Referenced by assembler().
|
private |
Definition at line 195 of file PressureAssembler.hpp.
Referenced by bases().
|
private |
Definition at line 200 of file PressureAssembler.hpp.
Referenced by compute_cavity_energy(), compute_cavity_energy_grad(), compute_cavity_energy_hess(), and PressureAssembler().
|
private |
Definition at line 196 of file PressureAssembler.hpp.
Referenced by gbases().
|
private |
Definition at line 191 of file PressureAssembler.hpp.
Referenced by mesh().
|
private |
Definition at line 193 of file PressureAssembler.hpp.
Referenced by compute_grad_volume(), compute_grad_volume_id(), compute_hess_volume_2d(), and compute_hess_volume_3d().
|
private |
Definition at line 203 of file PressureAssembler.hpp.
|
private |
Definition at line 192 of file PressureAssembler.hpp.
|
private |
Definition at line 202 of file PressureAssembler.hpp.
|
private |
Definition at line 197 of file PressureAssembler.hpp.
Referenced by compute_cavity_energy(), compute_cavity_energy_grad(), compute_cavity_energy_hess(), and problem().
|
private |
Definition at line 204 of file PressureAssembler.hpp.
|
private |
Definition at line 194 of file PressureAssembler.hpp.
Referenced by compute_cavity_energy_hess(), compute_energy_hess(), compute_grad_volume(), compute_grad_volume_id(), compute_hess_volume_2d(), and compute_hess_volume_3d().
|
private |
Definition at line 199 of file PressureAssembler.hpp.
Referenced by compute_cavity_energy(), compute_cavity_energy_grad(), compute_cavity_energy_hess(), and PressureAssembler().