Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SolveData.cpp
Go to the documentation of this file.
1#include "SolveData.hpp"
2
29
30#include <h5pp/h5pp.h>
31
32namespace polyfem::solver
33{
34 using namespace polyfem::time_integrator;
35
36 std::vector<std::shared_ptr<Form>> SolveData::init_forms(
37 // General
38 const Units &units,
39 const int dim,
40 const double t,
41 const Eigen::VectorXi &in_node_to_node,
42
43 // Elastic form
44 const int n_bases,
45 std::vector<basis::ElementBases> &bases,
46 const std::vector<basis::ElementBases> &geom_bases,
47 const assembler::Assembler &assembler,
48 assembler::AssemblyValsCache &ass_vals_cache,
49 const assembler::AssemblyValsCache &mass_ass_vals_cache,
50 const double jacobian_threshold,
51 const ElementInversionCheck check_inversion,
52
53 // Body form
54 const int n_pressure_bases,
55 const std::vector<int> &boundary_nodes,
56 const std::vector<mesh::LocalBoundary> &local_boundary,
57 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
58 const int n_boundary_samples,
59 const Eigen::MatrixXd &rhs,
60 const Eigen::MatrixXd &sol,
61 const assembler::Density &density,
62
63 // Pressure form
64 const std::vector<mesh::LocalBoundary> &local_pressure_boundary,
65 const std::unordered_map<int, std::vector<mesh::LocalBoundary>> &local_pressure_cavity,
66 const std::shared_ptr<assembler::PressureAssembler> pressure_assembler,
67
68 // Inertia form
69 const bool ignore_inertia,
70 const StiffnessMatrix &mass,
71 const std::shared_ptr<assembler::ViscousDamping> damping_assembler,
72
73 // Lagged regularization form
74 const double lagged_regularization_weight,
75 const int lagged_regularization_iterations,
76
77 // Augemented lagrangian form
78 const size_t obstacle_ndof,
79 const std::vector<std::string> &hard_constraint_files,
80 const std::vector<json> &soft_constraint_files,
81
82 // Contact form
83 const bool contact_enabled,
84 const ipc::CollisionMesh &collision_mesh,
85 const double dhat,
86 const double avg_mass,
87 const bool use_area_weighting,
88 const bool use_improved_max_operator,
89 const bool use_physical_barrier,
90 const json &barrier_stiffness,
91 const double initial_barrier_stiffness,
92 const ipc::BroadPhaseMethod broad_phase,
93 const double ccd_tolerance,
94 const long ccd_max_iterations,
95 const bool enable_shape_derivatives,
96
97 // Smooth Contact Form
98 const bool use_gcp_formulation,
99 const double alpha_t,
100 const double alpha_n,
101 const bool use_adaptive_dhat,
102 const double min_distance_ratio,
103
104 // Normal Adhesion Form
105 const bool adhesion_enabled,
106 const double dhat_p,
107 const double dhat_a,
108 const double Y,
109
110 // Tangential Adhesion Form
111 const double tangential_adhesion_coefficient,
112 const double epsa,
113 const int tangential_adhesion_iterations,
114
115 // Homogenization
116 const assembler::MacroStrainValue &macro_strain_constraint,
117
118 // Periodic contact
119 const bool periodic_contact,
120 const Eigen::VectorXi &tiled_to_single,
121 const std::shared_ptr<utils::PeriodicBoundary> &periodic_bc,
122
123 // Friction form
124 const double friction_coefficient,
125 const double epsv,
126 const int friction_iterations,
127
128 // Rayleigh damping form
129 const json &rayleigh_damping)
130 {
131 const bool is_time_dependent = time_integrator != nullptr;
132 assert(!is_time_dependent || time_integrator != nullptr);
133 const double dt = is_time_dependent ? time_integrator->dt() : 0.0;
134 const int ndof = n_bases * dim;
135 // if (is_formulation_mixed) // mixed not supported
136 // ndof_ += n_pressure_bases; // pressure is a scalar
137 const bool is_volume = dim == 3;
138
139 std::vector<std::shared_ptr<Form>> forms;
140 al_form.clear();
141
142 elastic_form = std::make_shared<ElasticForm>(
143 n_bases, bases, geom_bases, assembler, ass_vals_cache,
144 t, dt, is_volume, jacobian_threshold, check_inversion);
145 forms.push_back(elastic_form);
146
147 if (rhs_assembler != nullptr)
148 {
149 body_form = std::make_shared<BodyForm>(
150 ndof, n_pressure_bases, boundary_nodes, local_boundary,
151 local_neumann_boundary, n_boundary_samples, rhs, *rhs_assembler,
152 density, /*is_formulation_mixed=*/false,
153 is_time_dependent);
154 body_form->update_quantities(t, sol);
155 forms.push_back(body_form);
156 }
157
158 if (pressure_assembler != nullptr)
159 {
160 pressure_form = std::make_shared<PressureForm>(
161 ndof,
162 local_pressure_boundary,
163 local_pressure_cavity,
164 boundary_nodes,
165 n_boundary_samples, *pressure_assembler,
166 is_time_dependent);
167 pressure_form->update_quantities(t, sol);
168 forms.push_back(pressure_form);
169 }
170
171 inertia_form = nullptr;
172 damping_form = nullptr;
173 if (is_time_dependent)
174 {
175 if (!ignore_inertia)
176 {
177 assert(time_integrator != nullptr);
178 inertia_form = std::make_shared<InertiaForm>(mass, *time_integrator);
179 forms.push_back(inertia_form);
180 }
181
182 if (damping_assembler != nullptr)
183 {
184 damping_form = std::make_shared<ElasticForm>(
185 n_bases, bases, geom_bases, *damping_assembler, ass_vals_cache, t, dt, is_volume);
186 forms.push_back(damping_form);
187 }
188 }
189 else
190 {
191 if (lagged_regularization_weight > 0)
192 {
193 forms.push_back(std::make_shared<LaggedRegForm>(lagged_regularization_iterations));
194 forms.back()->set_weight(lagged_regularization_weight);
195 }
196 }
197
198 if (rhs_assembler != nullptr)
199 {
200 // assembler::Mass mass_mat_assembler;
201 // mass_mat_assembler.set_size(dim);
202 StiffnessMatrix mass_tmp = mass;
203 // mass_mat_assembler.assemble(dim == 3, n_bases, bases, geom_bases, mass_ass_vals_cache, mass_tmp, true);
204 // assert(mass_tmp.rows() == mass.rows() && mass_tmp.cols() == mass.cols());
205
206 if (!boundary_nodes.empty())
207 al_form.push_back(std::make_shared<BCLagrangianForm>(
208 ndof, boundary_nodes, local_boundary, local_neumann_boundary,
209 n_boundary_samples, mass_tmp, *rhs_assembler, obstacle_ndof, is_time_dependent, t));
210 // forms.push_back(al_form.back());
211 }
212
213 if (periodic_bc != nullptr)
214 {
215 al_form.push_back(std::make_shared<PeriodicLagrangianForm>(ndof, periodic_bc));
216 }
217
218 for (const auto &path : hard_constraint_files)
219 {
220 logger().debug("Setting up hard constraints for {}", path);
221 h5pp::File file(path, h5pp::FileAccess::READONLY);
222 std::vector<int> local2global;
223 if (!file.findDatasets("local2global").empty())
224 local2global = file.readDataset<std::vector<int>>("local2global");
225
226 if (local2global.empty())
227 {
228 local2global.resize(in_node_to_node.size());
229
230 for (int i = 0; i < local2global.size(); ++i)
231 local2global[i] = in_node_to_node[i];
232 }
233 else
234 {
235 for (auto &v : local2global)
236 v = in_node_to_node[v];
237 }
238
239 Eigen::MatrixXd bin = file.readDataset<Eigen::MatrixXd>("b");
240
241 StiffnessMatrix A, A_proj;
242 Eigen::MatrixXd b, b_proj;
243
244 if (!file.findDatasets("A").empty())
245 {
246 Eigen::MatrixXd Ain = file.readDataset<Eigen::MatrixXd>("A");
247 utils::scatter_matrix(ndof, dim, Ain, bin, local2global, A, b);
248
249 if (!file.findDatasets("A_proj").empty())
250 {
251 Eigen::MatrixXd A_proj_in = file.readDataset<Eigen::MatrixXd>("A_proj");
252 if (file.findDatasets("b_proj").empty())
253 log_and_throw_error("Missing b_proj in hard constraint file");
254
255 Eigen::MatrixXd b_proj_in = file.readDataset<Eigen::MatrixXd>("b_proj");
256 utils::scatter_matrix_col(ndof, dim, A_proj_in, b_proj_in, local2global, A_proj, b_proj);
257 }
258 }
259 else
260 {
261 std::vector<double> values = file.readDataset<std::vector<double>>("A_triplets/values");
262 std::vector<int> rows = file.readDataset<std::vector<int>>("A_triplets/rows");
263 std::vector<int> cols = file.readDataset<std::vector<int>>("A_triplets/cols");
264 std::vector<long> shape = file.readDataset<std::vector<long>>("A_triplets/shape");
265 utils::scatter_matrix(ndof, dim, shape, rows, cols, values, bin, local2global, A, b);
266
267 if (!file.findGroups("A_proj_triplets").empty())
268 {
269 if (file.findDatasets("b_proj").empty())
270 log_and_throw_error("Missing b_proj in hard constraint file");
271 if (file.findDatasets("rows", "/A_proj_triplets").empty())
272 log_and_throw_error("Missing A_proj_triplets/rows in hard constraint file");
273 if (file.findDatasets("cols", "/A_proj_triplets").empty())
274 log_and_throw_error("Missing A_proj_triplets/cols in hard constraint file");
275 if (file.findDatasets("values", "/A_proj_triplets").empty())
276 log_and_throw_error("Missing A_proj_triplets/values in hard constraint file");
277
278 std::vector<double> values_proj = file.readDataset<std::vector<double>>("A_proj_triplets/values");
279 std::vector<int> rows_proj = file.readDataset<std::vector<int>>("A_proj_triplets/rows");
280 std::vector<int> cols_proj = file.readDataset<std::vector<int>>("A_proj_triplets/cols");
281 Eigen::MatrixXd b_projin = file.readDataset<Eigen::MatrixXd>("b_proj");
282 std::vector<long> shape_proj = file.readDataset<std::vector<long>>("A_proj_triplets/shape");
283
284 utils::scatter_matrix_col(ndof, dim, shape_proj, rows_proj, cols_proj, values_proj, b_projin, local2global, A_proj, b_proj);
285 }
286 }
287
288 al_form.push_back(std::make_shared<MatrixLagrangianForm>(A, b, A_proj, b_proj));
289 // forms.push_back(al_form.back());
290 }
291
292 for (const auto &j : soft_constraint_files)
293 {
294 const std::string &path = j["data"];
295 double weight = j["weight"];
296
297 logger().debug("Setting up soft constraints for {}", path);
298 h5pp::File file(path, h5pp::FileAccess::READONLY);
299 std::vector<int> local2global;
300 if (!file.findDatasets("local2global").empty())
301 local2global = file.readDataset<std::vector<int>>("local2global");
302
303 if (local2global.empty())
304 {
305 local2global.resize(in_node_to_node.size());
306
307 for (int i = 0; i < local2global.size(); ++i)
308 local2global[i] = in_node_to_node[i];
309 }
310 else
311 {
312 for (auto &v : local2global)
313 v = in_node_to_node[v];
314 }
315
316 Eigen::MatrixXd bin = file.readDataset<Eigen::MatrixXd>("b");
317
319 Eigen::MatrixXd b;
320
321 if (!file.findDatasets("A").empty())
322 {
323 Eigen::MatrixXd Ain = file.readDataset<Eigen::MatrixXd>("A");
324 utils::scatter_matrix(ndof, dim, Ain, bin, local2global, A, b);
325 }
326 else
327 {
328 std::vector<double> values = file.readDataset<std::vector<double>>("A_triplets/values");
329 std::vector<int> rows = file.readDataset<std::vector<int>>("A_triplets/rows");
330 std::vector<int> cols = file.readDataset<std::vector<int>>("A_triplets/cols");
331 std::vector<long> shape = file.readDataset<std::vector<long>>("A_triplets/shape");
332
333 utils::scatter_matrix(ndof, dim, shape, rows, cols, values, bin, local2global, A, b);
334 }
335
336 forms.push_back(std::make_shared<QuadraticPenaltyForm>(A, b, weight));
337 }
338
339 if (macro_strain_constraint.is_active())
340 {
341 // don't push these two into forms because they take a different input x
342 strain_al_lagr_form = std::make_shared<MacroStrainLagrangianForm>(macro_strain_constraint);
343 }
344
345 contact_form = nullptr;
346 periodic_contact_form = nullptr;
347 friction_form = nullptr;
348 if (contact_enabled)
349 {
350 const bool use_adaptive_barrier_stiffness = !barrier_stiffness.is_number();
351
352 if (periodic_contact)
353 {
354 periodic_contact_form = std::make_shared<PeriodicContactForm>(
355 collision_mesh, tiled_to_single, dhat, avg_mass, use_area_weighting, use_improved_max_operator, use_physical_barrier,
356 use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase, ccd_tolerance,
357 ccd_max_iterations);
358
359 if (use_adaptive_barrier_stiffness)
360 {
361 periodic_contact_form->set_barrier_stiffness(1);
362 // logger().debug("Using adaptive barrier stiffness");
363 }
364 else
365 {
366 assert(barrier_stiffness.is_number());
367 assert(barrier_stiffness.get<double>() > 0);
368 periodic_contact_form->set_barrier_stiffness(barrier_stiffness);
369 // logger().debug("Using fixed barrier stiffness of {}", contact_form->barrier_stiffness());
370 }
371
372 // periodic_contact_form is not pushed into forms since it takes different input vectors.
373 }
374 else
375 {
376 if (use_gcp_formulation)
377 {
378 contact_form = std::make_shared<SmoothContactForm>(
379 collision_mesh, dhat, avg_mass, alpha_t, alpha_n, use_adaptive_dhat, min_distance_ratio,
380 use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase,
381 ccd_tolerance * units.characteristic_length(), ccd_max_iterations);
382 }
383 else
384 {
385 contact_form = std::make_shared<BarrierContactForm>(
386 collision_mesh, dhat, avg_mass, use_area_weighting, use_improved_max_operator, use_physical_barrier,
387 use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase, ccd_tolerance * units.characteristic_length(),
388 ccd_max_iterations);
389 }
390
391 if (use_adaptive_barrier_stiffness)
392 {
393 contact_form->set_barrier_stiffness(initial_barrier_stiffness);
394 // logger().debug("Using adaptive barrier stiffness");
395 }
396 else
397 {
398 assert(barrier_stiffness.is_number());
399 assert(barrier_stiffness.get<double>() > 0);
400 contact_form->set_barrier_stiffness(barrier_stiffness);
401 // logger().debug("Using fixed barrier stiffness of {}", contact_form->barrier_stiffness());
402 }
403
404 forms.push_back(contact_form);
405 }
406
407 if (friction_coefficient != 0)
408 {
409 friction_form = std::make_shared<FrictionForm>(
410 collision_mesh, time_integrator, epsv, friction_coefficient,
411 broad_phase, *contact_form, friction_iterations);
412 friction_form->init_lagging(sol);
413 forms.push_back(friction_form);
414 }
415
416 if (adhesion_enabled)
417 {
418 normal_adhesion_form = std::make_shared<NormalAdhesionForm>(
419 collision_mesh, dhat_p, dhat_a, Y, is_time_dependent, enable_shape_derivatives,
420 broad_phase, ccd_tolerance * units.characteristic_length(), ccd_max_iterations
421 );
422 forms.push_back(normal_adhesion_form);
423
424 if (tangential_adhesion_coefficient != 0)
425 {
426 tangential_adhesion_form = std::make_shared<TangentialAdhesionForm>(
427 collision_mesh, time_integrator, epsa, tangential_adhesion_coefficient,
428 broad_phase, *normal_adhesion_form, tangential_adhesion_iterations
429 );
430 forms.push_back(tangential_adhesion_form);
431 }
432 }
433
434
435 }
436
437 const std::vector<json> rayleigh_damping_jsons = utils::json_as_array(rayleigh_damping);
438 if (is_time_dependent)
439 {
440 // Map from form name to form so RayleighDampingForm::create can get the correct form to damp
441 const std::unordered_map<std::string, std::shared_ptr<Form>> possible_forms_to_damp = {
442 {"elasticity", elastic_form},
443 {"contact", contact_form},
444 };
445
446 for (const json &params : rayleigh_damping_jsons)
447 {
448 forms.push_back(RayleighDampingForm::create(
449 params, possible_forms_to_damp,
451 }
452 }
453 else if (rayleigh_damping_jsons.size() > 0)
454 {
455 log_and_throw_adjoint_error("Rayleigh damping is only supported for time-dependent problems");
456 }
457
458 update_dt();
459
460 return forms;
461 }
462
463 void SolveData::update_barrier_stiffness(const Eigen::VectorXd &x)
464 {
465 if (contact_form == nullptr || !contact_form->use_adaptive_barrier_stiffness())
466 return;
467
468 Eigen::VectorXd grad_energy = Eigen::VectorXd::Zero(x.size());
469 const std::array<std::shared_ptr<Form>, 4> energy_forms{
471 for (const std::shared_ptr<Form> &form : energy_forms)
472 {
473 if (form == nullptr || !form->enabled())
474 continue;
475
476 Eigen::VectorXd grad_form;
477 form->first_derivative(x, grad_form);
478 grad_energy += grad_form;
479 }
480
481 contact_form->update_barrier_stiffness(x, grad_energy);
482 }
483
485 {
486 if (time_integrator == nullptr) // if is not time dependent
487 return;
488
489 const std::array<std::shared_ptr<Form>, 6> energy_forms{
491 for (const std::shared_ptr<Form> &form : energy_forms)
492 {
493 if (form == nullptr)
494 continue;
495 form->set_weight(time_integrator->acceleration_scaling());
496 }
497 }
498
499 std::vector<std::pair<std::string, std::shared_ptr<solver::Form>>> SolveData::named_forms() const
500 {
501 std::vector<std::pair<std::string, std::shared_ptr<solver::Form>>> res{
502 {"elastic", elastic_form},
503 {"inertia", inertia_form},
504 {"body", body_form},
505 {"contact", contact_form},
506 {"friction", friction_form},
507 {"damping", damping_form},
508 {"pressure", pressure_form},
509 {"strain_augmented_lagrangian_lagr", strain_al_lagr_form},
510 {"periodic_contact", periodic_contact_form},
511 };
512
513 for (const auto &form : al_form)
514 res.push_back({"augmented_lagrangian", form});
515
516 return res;
517 }
518} // 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::NormalAdhesionForm > normal_adhesion_form
std::shared_ptr< solver::ContactForm > contact_form
std::vector< std::shared_ptr< Form > > init_forms(const Units &units, const int dim, const double t, const Eigen::VectorXi &in_node_to_node, const int n_bases, std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &geom_bases, const assembler::Assembler &assembler, assembler::AssemblyValsCache &ass_vals_cache, const assembler::AssemblyValsCache &mass_ass_vals_cache, const double jacobian_threshold, const solver::ElementInversionCheck check_inversion, 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 std::vector< std::string > &hard_constraint_files, const std::vector< json > &soft_constraint_files, const bool contact_enabled, const ipc::CollisionMesh &collision_mesh, const double dhat, const double avg_mass, const bool use_area_weighting, const bool use_improved_max_operator, const bool use_physical_barrier, const json &barrier_stiffness, const double initial_barrier_stiffness, const ipc::BroadPhaseMethod broad_phase, const double ccd_tolerance, const long ccd_max_iterations, const bool enable_shape_derivatives, const bool use_gcp_formulation, const double alpha_t, const double alpha_n, const bool use_adaptive_dhat, const double min_distance_ratio, const bool adhesion_enabled, const double dhat_p, const double dhat_a, const double Y, const double tangential_adhesion_coefficient, const double epsa, const int tangential_adhesion_iterations, 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:36
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< solver::TangentialAdhesionForm > tangential_adhesion_form
std::shared_ptr< assembler::PressureAssembler > pressure_assembler
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
void scatter_matrix(const int n_dofs, const int dim, const Eigen::MatrixXd &A, const Eigen::MatrixXd &b, const std::vector< int > &local_to_global, StiffnessMatrix &Aout, Eigen::MatrixXd &bout)
std::vector< T > json_as_array(const json &j)
Return the value of a json object as an array.
Definition JSONUtils.hpp:38
void scatter_matrix_col(const int n_dofs, const int dim, const Eigen::MatrixXd &A, const Eigen::MatrixXd &b, const std::vector< int > &local_to_global, StiffnessMatrix &Aout, Eigen::MatrixXd &bout)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:79
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22