PolyFEM
Loading...
Searching...
No Matches
SolveData.cpp
Go to the documentation of this file.
1#include "SolveData.hpp"
2
22
23namespace polyfem::solver
24{
25 using namespace polyfem::time_integrator;
26
27 std::vector<std::shared_ptr<Form>> SolveData::init_forms(
28 // General
29 const Units &units,
30 const int dim,
31 const double t,
32
33 // Elastic form
34 const int n_bases,
35 const std::vector<basis::ElementBases> &bases,
36 const std::vector<basis::ElementBases> &geom_bases,
37 const assembler::Assembler &assembler,
38 const assembler::AssemblyValsCache &ass_vals_cache,
39 const assembler::AssemblyValsCache &mass_ass_vals_cache,
40
41 // Body form
42 const int n_pressure_bases,
43 const std::vector<int> &boundary_nodes,
44 const std::vector<mesh::LocalBoundary> &local_boundary,
45 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
46 const int n_boundary_samples,
47 const Eigen::MatrixXd &rhs,
48 const Eigen::MatrixXd &sol,
49 const assembler::Density &density,
50
51 // Pressure form
52 const std::vector<mesh::LocalBoundary> &local_pressure_boundary,
53 const std::unordered_map<int, std::vector<mesh::LocalBoundary>> &local_pressure_cavity,
54 const std::shared_ptr<assembler::PressureAssembler> pressure_assembler,
55
56 // Inertia form
57 const bool ignore_inertia,
58 const StiffnessMatrix &mass,
59 const std::shared_ptr<assembler::ViscousDamping> damping_assembler,
60
61 // Lagged regularization form
62 const double lagged_regularization_weight,
63 const int lagged_regularization_iterations,
64
65 // Augemented lagrangian form
66 // const std::vector<int> &boundary_nodes,
67 // const std::vector<mesh::LocalBoundary> &local_boundary,
68 // const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
69 // const int n_boundary_samples,
70 // const StiffnessMatrix &mass,
71 const size_t obstacle_ndof,
72
73 // Contact form
74 const bool contact_enabled,
75 const ipc::CollisionMesh &collision_mesh,
76 const double dhat,
77 const double avg_mass,
78 const bool use_convergent_contact_formulation,
79 const json &barrier_stiffness,
80 const ipc::BroadPhaseMethod broad_phase,
81 const double ccd_tolerance,
82 const long ccd_max_iterations,
83 const bool enable_shape_derivatives,
84
85 // Homogenization
86 const assembler::MacroStrainValue &macro_strain_constraint,
87
88 // Periodic contact
89 const bool periodic_contact,
90 const Eigen::VectorXi &tiled_to_single,
91 const std::shared_ptr<utils::PeriodicBoundary> &periodic_bc,
92
93 // Friction form
94 const double friction_coefficient,
95 const double epsv,
96 const int friction_iterations,
97
98 // Rayleigh damping form
99 const json &rayleigh_damping)
100 {
101 const bool is_time_dependent = time_integrator != nullptr;
102 assert(!is_time_dependent || time_integrator != nullptr);
103 const double dt = is_time_dependent ? time_integrator->dt() : 0.0;
104 const int ndof = n_bases * dim;
105 // if (is_formulation_mixed) // mixed not supported
106 // ndof_ += n_pressure_bases; // pressure is a scalar
107 const bool is_volume = dim == 3;
108
109 std::vector<std::shared_ptr<Form>> forms;
110 al_form.clear();
111
112 elastic_form = std::make_shared<ElasticForm>(
113 n_bases, bases, geom_bases, assembler, ass_vals_cache,
114 t, dt, is_volume);
115 forms.push_back(elastic_form);
116
117 if (rhs_assembler != nullptr)
118 {
119 body_form = std::make_shared<BodyForm>(
120 ndof, n_pressure_bases, boundary_nodes, local_boundary,
121 local_neumann_boundary, n_boundary_samples, rhs, *rhs_assembler,
122 density, /*is_formulation_mixed=*/false,
123 is_time_dependent);
124 body_form->update_quantities(t, sol);
125 forms.push_back(body_form);
126 }
127
128 if (pressure_assembler != nullptr)
129 {
130 pressure_form = std::make_shared<PressureForm>(
131 ndof,
132 local_pressure_boundary,
133 local_pressure_cavity,
134 boundary_nodes,
135 n_boundary_samples, *pressure_assembler,
136 is_time_dependent);
137 pressure_form->update_quantities(t, sol);
138 forms.push_back(pressure_form);
139 }
140
141 inertia_form = nullptr;
142 damping_form = nullptr;
143 if (is_time_dependent)
144 {
145 if (!ignore_inertia)
146 {
147 assert(time_integrator != nullptr);
148 inertia_form = std::make_shared<InertiaForm>(mass, *time_integrator);
149 forms.push_back(inertia_form);
150 }
151
152 if (damping_assembler != nullptr)
153 {
154 damping_form = std::make_shared<ElasticForm>(
155 n_bases, bases, geom_bases, *damping_assembler, ass_vals_cache, t, dt, is_volume);
156 forms.push_back(damping_form);
157 }
158 }
159 else
160 {
161 if (lagged_regularization_weight > 0)
162 {
163 forms.push_back(std::make_shared<LaggedRegForm>(lagged_regularization_iterations));
164 forms.back()->set_weight(lagged_regularization_weight);
165 }
166 }
167
168 if (rhs_assembler != nullptr)
169 {
170 // assembler::Mass mass_mat_assembler;
171 // mass_mat_assembler.set_size(dim);
172 StiffnessMatrix mass_tmp = mass;
173 // mass_mat_assembler.assemble(dim == 3, n_bases, bases, geom_bases, mass_ass_vals_cache, mass_tmp, true);
174 // assert(mass_tmp.rows() == mass.rows() && mass_tmp.cols() == mass.cols());
175
176 al_form.push_back(std::make_shared<BCLagrangianForm>(
177 ndof, boundary_nodes, local_boundary, local_neumann_boundary,
178 n_boundary_samples, mass_tmp, *rhs_assembler, obstacle_ndof, is_time_dependent, t, periodic_bc));
179 forms.push_back(al_form.back());
180 }
181
182 if (macro_strain_constraint.is_active())
183 {
184 // don't push these two into forms because they take a different input x
185 strain_al_lagr_form = std::make_shared<MacroStrainLagrangianForm>(macro_strain_constraint);
186 }
187
188 contact_form = nullptr;
189 periodic_contact_form = nullptr;
190 friction_form = nullptr;
191 if (contact_enabled)
192 {
193 const bool use_adaptive_barrier_stiffness = !barrier_stiffness.is_number();
194
195 if (periodic_contact)
196 {
197 periodic_contact_form = std::make_shared<PeriodicContactForm>(
198 collision_mesh, tiled_to_single, dhat, avg_mass, use_convergent_contact_formulation,
199 use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase, ccd_tolerance,
200 ccd_max_iterations);
201
202 if (use_adaptive_barrier_stiffness)
203 {
204 periodic_contact_form->set_barrier_stiffness(1);
205 // logger().debug("Using adaptive barrier stiffness");
206 }
207 else
208 {
209 assert(barrier_stiffness.is_number());
210 assert(barrier_stiffness.get<double>() > 0);
211 periodic_contact_form->set_barrier_stiffness(barrier_stiffness);
212 // logger().debug("Using fixed barrier stiffness of {}", contact_form->barrier_stiffness());
213 }
214 }
215 else
216 {
217 contact_form = std::make_shared<ContactForm>(
218 collision_mesh, dhat, avg_mass, use_convergent_contact_formulation,
219 use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase, ccd_tolerance * units.characteristic_length(),
220 ccd_max_iterations);
221
222 if (use_adaptive_barrier_stiffness)
223 {
224 contact_form->set_barrier_stiffness(1);
225 // logger().debug("Using adaptive barrier stiffness");
226 }
227 else
228 {
229 assert(barrier_stiffness.is_number());
230 assert(barrier_stiffness.get<double>() > 0);
231 contact_form->set_barrier_stiffness(barrier_stiffness);
232 // logger().debug("Using fixed barrier stiffness of {}", contact_form->barrier_stiffness());
233 }
234
235 if (contact_form)
236 forms.push_back(contact_form);
237
238 // ----------------------------------------------------------------
239 }
240
241 if (friction_coefficient != 0)
242 {
243 friction_form = std::make_shared<FrictionForm>(
244 collision_mesh, time_integrator, epsv, friction_coefficient,
245 broad_phase, *contact_form, friction_iterations);
246 friction_form->init_lagging(sol);
247 forms.push_back(friction_form);
248 }
249 }
250
251 const std::vector<json> rayleigh_damping_jsons = utils::json_as_array(rayleigh_damping);
252 if (is_time_dependent)
253 {
254 // Map from form name to form so RayleighDampingForm::create can get the correct form to damp
255 const std::unordered_map<std::string, std::shared_ptr<Form>> possible_forms_to_damp = {
256 {"elasticity", elastic_form},
257 {"contact", contact_form},
258 };
259
260 for (const json &params : rayleigh_damping_jsons)
261 {
262 forms.push_back(RayleighDampingForm::create(
263 params, possible_forms_to_damp,
265 }
266 }
267 else if (rayleigh_damping_jsons.size() > 0)
268 {
269 log_and_throw_adjoint_error("Rayleigh damping is only supported for time-dependent problems");
270 }
271
272 update_dt();
273
274 return forms;
275 }
276
277 void SolveData::update_barrier_stiffness(const Eigen::VectorXd &x)
278 {
279 if (contact_form == nullptr || !contact_form->use_adaptive_barrier_stiffness())
280 return;
281
282 Eigen::VectorXd grad_energy = Eigen::VectorXd::Zero(x.size());
283 const std::array<std::shared_ptr<Form>, 4> energy_forms{
285 for (const std::shared_ptr<Form> &form : energy_forms)
286 {
287 if (form == nullptr || !form->enabled())
288 continue;
289
290 Eigen::VectorXd grad_form;
291 form->first_derivative(x, grad_form);
292 grad_energy += grad_form;
293 }
294
295 contact_form->update_barrier_stiffness(x, grad_energy);
296 }
297
299 {
300 if (time_integrator == nullptr) // if is not time dependent
301 return;
302
303 const std::array<std::shared_ptr<Form>, 6> energy_forms{
305 for (const std::shared_ptr<Form> &form : energy_forms)
306 {
307 if (form == nullptr)
308 continue;
309 form->set_weight(time_integrator->acceleration_scaling());
310 }
311 }
312
313 std::vector<std::pair<std::string, std::shared_ptr<solver::Form>>> SolveData::named_forms() const
314 {
315 std::vector<std::pair<std::string, std::shared_ptr<solver::Form>>> res{
316 {"elastic", elastic_form},
317 {"inertia", inertia_form},
318 {"body", body_form},
319 {"contact", contact_form},
320 {"friction", friction_form},
321 {"damping", damping_form},
322 {"pressure", pressure_form},
323 {"strain_augmented_lagrangian_lagr", strain_al_lagr_form},
324 {"periodic_contact", periodic_contact_form},
325 };
326
327 for (const auto &form : al_form)
328 res.push_back({"augmented_lagrangian", form});
329
330 return res;
331 }
332} // namespace polyfem::solver
int x
double characteristic_length() const
Definition Units.hpp:22
Caches basis evaluation and geometric mapping at every element.
static std::shared_ptr< RayleighDampingForm > create(const json &params, const std::unordered_map< std::string, std::shared_ptr< Form > > &forms, const time_integrator::ImplicitTimeIntegrator &time_integrator)
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::MacroStrainLagrangianForm > strain_al_lagr_form
std::shared_ptr< solver::ContactForm > contact_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< solver::AugmentedLagrangianForm > > al_form
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< assembler::PressureAssembler > pressure_assembler
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 std::shared_ptr< utils::PeriodicBoundary > &periodic_bc, 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:27
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
std::vector< T > json_as_array(const json &j)
Return the value of a json object as an array.
Definition JSONUtils.hpp:38
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22