PolyFEM
Loading...
Searching...
No Matches
SolveData.hpp
Go to the documentation of this file.
1#pragma once
2
10
11#include <ipc/collision_mesh.hpp>
12#include <ipc/broad_phase/broad_phase.hpp>
13
14#include <Eigen/Core>
15
16#include <memory>
17#include <string>
18#include <unordered_map>
19
21{
22 class ImplicitTimeIntegrator;
23} // namespace polyfem::time_integrator
24
25namespace polyfem::assembler
26{
27 class ViscousDamping;
28 class MacroStrainValue;
29 class PressureAssembler;
30} // namespace polyfem::assembler
31
32namespace polyfem::solver
33{
34 class NLProblem;
35 class Form;
36 class ContactForm;
37 class PeriodicContactForm;
38 class MacroStrainALForm;
39 class FrictionForm;
40 class BodyForm;
41 class BCLagrangianForm;
42 class BCPenaltyForm;
43 class MacroStrainLagrangianForm;
44 class MacroStrainALForm;
45 class InertiaForm;
46 class ElasticForm;
47 class PressureForm;
48
51 {
52 public:
55 std::vector<std::shared_ptr<Form>> init_forms(
56 // General
57 const Units &units,
58 const int dim,
59 const double t,
60
61 // Elastic form
62 const int n_bases,
63 const std::vector<basis::ElementBases> &bases,
64 const std::vector<basis::ElementBases> &geom_bases,
65 const assembler::Assembler &assembler,
66 const assembler::AssemblyValsCache &ass_vals_cache,
67 const assembler::AssemblyValsCache &mass_ass_vals_cache,
68
69 // Body form
70 const int n_pressure_bases,
71 const std::vector<int> &boundary_nodes,
72 const std::vector<mesh::LocalBoundary> &local_boundary,
73 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
74 const int n_boundary_samples,
75 const Eigen::MatrixXd &rhs,
76 const Eigen::MatrixXd &sol,
77 const assembler::Density &density,
78
79 // Pressure form
80 const std::vector<mesh::LocalBoundary> &local_pressure_boundary,
81 const std::unordered_map<int, std::vector<mesh::LocalBoundary>> &local_pressure_cavity,
82 const std::shared_ptr<assembler::PressureAssembler> pressure_assembler,
83
84 // Inertia form
85 const bool ignore_inertia,
86 const StiffnessMatrix &mass,
87 const std::shared_ptr<assembler::ViscousDamping> damping_assembler,
88
89 // Lagged regularization form
90 const double lagged_regularization_weight,
91 const int lagged_regularization_iterations,
92
93 // Augemented lagrangian form
94 // const std::vector<int> &boundary_nodes,
95 // const std::vector<mesh::LocalBoundary> &local_boundary,
96 // const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
97 // const int n_boundary_samples,
98 // const StiffnessMatrix &mass,
99 const size_t obstacle_ndof,
100
101 // Contact form
102 const bool contact_enabled,
103 const ipc::CollisionMesh &collision_mesh,
104 const double dhat,
105 const double avg_mass,
106 const bool use_convergent_contact_formulation,
107 const json &barrier_stiffness,
108 const ipc::BroadPhaseMethod broad_phase,
109 const double ccd_tolerance,
110 const long ccd_max_iterations,
111 const bool enable_shape_derivatives,
112
113 // Homogenization
114 const assembler::MacroStrainValue &macro_strain_constraint,
115
116 // Periodic contact
117 const bool periodic_contact,
118 const Eigen::VectorXi &tiled_to_single,
119
120 // Friction form
121 const double friction_coefficient,
122 const double epsv,
123 const int friction_iterations,
124
125 // Rayleigh damping form
126 const json &rayleigh_damping);
127
130 void update_barrier_stiffness(const Eigen::VectorXd &x);
131
133 void update_dt();
134
135 std::vector<std::pair<std::string, std::shared_ptr<solver::Form>>> named_forms() const;
136
137 public:
138 std::shared_ptr<assembler::RhsAssembler> rhs_assembler;
139 std::shared_ptr<assembler::PressureAssembler> pressure_assembler;
140 std::shared_ptr<solver::NLProblem> nl_problem;
141
142 std::shared_ptr<solver::BCLagrangianForm> al_lagr_form;
143 std::shared_ptr<solver::BCPenaltyForm> al_pen_form;
144 std::shared_ptr<solver::MacroStrainLagrangianForm> strain_al_lagr_form;
145 std::shared_ptr<solver::MacroStrainALForm> strain_al_pen_form;
146 std::shared_ptr<solver::BodyForm> body_form;
147 std::shared_ptr<solver::ContactForm> contact_form;
148 std::shared_ptr<solver::ElasticForm> damping_form;
149 std::shared_ptr<solver::ElasticForm> elastic_form;
150 std::shared_ptr<solver::FrictionForm> friction_form;
151 std::shared_ptr<solver::InertiaForm> inertia_form;
152 std::shared_ptr<solver::PressureForm> pressure_form;
153
154 std::shared_ptr<solver::PeriodicContactForm> periodic_contact_form;
155
156 std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator;
157 };
158} // namespace polyfem::solver
int x
Caches basis evaluation and geometric mapping at every element.
class to store time stepping data
Definition SolveData.hpp:51
std::shared_ptr< solver::FrictionForm > friction_form
std::shared_ptr< solver::InertiaForm > inertia_form
std::shared_ptr< solver::PeriodicContactForm > periodic_contact_form
void update_dt()
updates the dt inside the different forms
std::shared_ptr< solver::PressureForm > pressure_form
std::shared_ptr< solver::BodyForm > body_form
std::shared_ptr< solver::BCPenaltyForm > al_pen_form
std::shared_ptr< solver::NLProblem > nl_problem
std::shared_ptr< solver::MacroStrainLagrangianForm > strain_al_lagr_form
std::shared_ptr< solver::ContactForm > contact_form
std::shared_ptr< solver::MacroStrainALForm > strain_al_pen_form
std::vector< std::pair< std::string, std::shared_ptr< solver::Form > > > named_forms() const
std::shared_ptr< solver::ElasticForm > damping_form
std::shared_ptr< solver::ElasticForm > elastic_form
std::vector< std::shared_ptr< Form > > init_forms(const Units &units, const int dim, const double t, 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 assembler::AssemblyValsCache &mass_ass_vals_cache, const int n_pressure_bases, const std::vector< int > &boundary_nodes, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const int n_boundary_samples, const Eigen::MatrixXd &rhs, const Eigen::MatrixXd &sol, const assembler::Density &density, const std::vector< mesh::LocalBoundary > &local_pressure_boundary, const std::unordered_map< int, std::vector< mesh::LocalBoundary > > &local_pressure_cavity, const std::shared_ptr< assembler::PressureAssembler > pressure_assembler, const bool ignore_inertia, const StiffnessMatrix &mass, const std::shared_ptr< assembler::ViscousDamping > damping_assembler, const double lagged_regularization_weight, const int lagged_regularization_iterations, const size_t obstacle_ndof, const bool contact_enabled, const ipc::CollisionMesh &collision_mesh, const double dhat, const double avg_mass, const bool use_convergent_contact_formulation, const json &barrier_stiffness, const ipc::BroadPhaseMethod broad_phase, const double ccd_tolerance, const long ccd_max_iterations, const bool enable_shape_derivatives, const assembler::MacroStrainValue &macro_strain_constraint, const bool periodic_contact, const Eigen::VectorXi &tiled_to_single, const double friction_coefficient, const double epsv, const int friction_iterations, const json &rayleigh_damping)
Initialize the forms and return a vector of pointers to them.
Definition SolveData.cpp:29
std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator
void update_barrier_stiffness(const Eigen::VectorXd &x)
update the barrier stiffness for the forms
std::shared_ptr< solver::BCLagrangianForm > al_lagr_form
std::shared_ptr< assembler::PressureAssembler > pressure_assembler
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
Used for test only.
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20