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