PolyFEM
Loading...
Searching...
No Matches
StateSolveNonlinear.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
5#include <polyfem/Common.hpp>
6
15
27
28#include <Eigen/Core>
29
30#include <ipc/ipc.hpp>
31
32#include <spdlog/fmt/fmt.h>
33
34#include <cassert>
35#include <memory>
36#include <string>
37#include <vector>
38
39namespace polyfem
40{
41 using namespace mesh;
42 using namespace solver;
43 using namespace time_integrator;
44 using namespace io;
45 using namespace utils;
46
47 std::shared_ptr<polysolve::nonlinear::Solver> State::make_nl_solver(bool for_al) const
48 {
49 return polysolve::nonlinear::Solver::create(for_al ? args["solver"]["augmented_lagrangian"]["nonlinear"] : args["solver"]["nonlinear"], args["solver"]["linear"], units.characteristic_length(), logger());
50 }
51
52 void State::solve_transient_tensor_nonlinear(const int time_steps,
53 const double t0,
54 const double dt,
55 Eigen::MatrixXd &sol,
56 UserPostStepCallback user_post_step,
57 const InitialConditionOverride *ic_override)
58 {
59 init_nonlinear_tensor_solve(sol, t0 + dt, true, ic_override);
60
61 // File index offset for continuing file numbering from a previous run.
62 // Set via restart JSON; defaults to 0 for fresh simulations.
63 const int t_offset = args["output"]["data"]["file_index_offset"].get<int>();
64
65 // Write the total energy to a CSV file
66 int save_i = 0;
67
68 std::unique_ptr<EnergyCSVWriter> energy_csv = nullptr;
69 std::unique_ptr<RuntimeStatsCSVWriter> stats_csv = nullptr;
70
71 if (args["output"]["stats"])
72 {
73 logger().debug("Saving nl stats to {} and {}", resolve_output_path("energy.csv"), resolve_output_path("stats.csv"));
74 energy_csv = std::make_unique<EnergyCSVWriter>(resolve_output_path("energy.csv"), solve_data);
75 stats_csv = std::make_unique<RuntimeStatsCSVWriter>(resolve_output_path("stats.csv"), *this, t0, dt);
76 }
77 const bool remesh_enabled = args["space"]["remesh"]["enabled"];
78 // const double save_dt = remesh_enabled ? (dt / 3) : dt;
79
80 // Save the initial solution
81 if (energy_csv)
82 energy_csv->write(save_i, sol);
83 save_timestep(t0, t_offset, t0, dt, sol, Eigen::MatrixXd()); // no pressure
84 save_i++;
85
86 // Step 0.
87 if (user_post_step)
88 {
89 user_post_step(0, *this, sol, nullptr, nullptr);
90 }
91
92 for (int t = 1; t <= time_steps; ++t)
93 {
94 double forward_solve_time = 0, remeshing_time = 0, global_relaxation_time = 0;
95
96 {
97 POLYFEM_SCOPED_TIMER(forward_solve_time);
98 solve_tensor_nonlinear(t, sol, true);
99 }
100
101 if (remesh_enabled)
102 {
103 if (energy_csv)
104 energy_csv->write(save_i, sol);
105 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
106 save_i++;
107
108 bool remesh_success;
109 {
110 POLYFEM_SCOPED_TIMER(remeshing_time);
111 remesh_success = this->remesh(t0 + dt * t, dt, sol);
112 }
113
114 // Save the solution after remeshing
115 if (energy_csv)
116 energy_csv->write(save_i, sol);
117 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
118 save_i++;
119
120 // Only do global relaxation if remeshing was successful
121 if (remesh_success)
122 {
123 POLYFEM_SCOPED_TIMER(global_relaxation_time);
124 solve_tensor_nonlinear(t, sol, false); // solve the scene again after remeshing
125 }
126 }
127
128 // Always save the solution for consistency
129 if (energy_csv)
130 energy_csv->write(save_i, sol);
131 save_timestep(t0 + dt * t, t + t_offset, t0, dt, sol, Eigen::MatrixXd()); // no pressure
132 save_i++;
133
134 if (user_post_step)
135 {
136 user_post_step(t, *this, sol, nullptr, nullptr);
137 }
138
139 {
140 POLYFEM_SCOPED_TIMER("Update quantities");
141
142 solve_data.time_integrator->update_quantities(sol);
143
144 solve_data.nl_problem->update_quantities(t0 + (t + 1) * dt, sol);
145
148 }
149
150 logger().info("{}/{} t={}", t, time_steps, t0 + dt * t);
151 if (time_callback)
152 time_callback(t, time_steps, t0 + dt * t, t0 + dt * time_steps);
153
154 const std::string rest_mesh_path = args["output"]["data"]["rest_mesh"].get<std::string>();
155 if (!rest_mesh_path.empty())
156 {
157 Eigen::MatrixXd V;
158 Eigen::MatrixXi F;
161 resolve_output_path(fmt::format(args["output"]["data"]["rest_mesh"], t + t_offset)),
162 V, F, mesh->get_body_ids(), mesh->is_volume(), /*binary=*/true);
163 }
164
165 const std::string &state_path = resolve_output_path(fmt::format(args["output"]["data"]["state"], t + t_offset));
166 if (!state_path.empty())
167 solve_data.time_integrator->save_state(state_path);
168
169 // save restart file
170 save_restart_json(t0, dt, t);
171 if (remesh_enabled && stats_csv)
172 stats_csv->write(t, forward_solve_time, remeshing_time, global_relaxation_time, sol);
173 }
174 }
175
176 void State::init_nonlinear_tensor_solve(Eigen::MatrixXd &sol,
177 const double t,
178 const bool init_time_integrator,
179 const InitialConditionOverride *ic_override)
180 {
181 assert(sol.cols() == 1);
182 assert(!is_problem_linear()); // non-linear
183 assert(!problem->is_scalar()); // tensor
184 assert(mixed_assembler == nullptr);
185
186 // --------------------------------------------------------------------
187 // Check for initial intersections
188 if (is_contact_enabled())
189 {
190 POLYFEM_SCOPED_TIMER("Check for initial intersections");
191
192 const Eigen::MatrixXd displaced = collision_mesh.displace_vertices(
193 utils::unflatten(sol, mesh->dimension()));
194
195 if (ipc::has_intersections(collision_mesh, displaced, ipc::create_broad_phase(args["solver"]["contact"]["CCD"]["broad_phase"]).get()))
196 {
198 resolve_output_path("intersection.obj"), displaced,
199 collision_mesh.edges(), collision_mesh.faces());
200 log_and_throw_error("Unable to solve, initial solution has intersections!");
201 }
202 }
203
204 // --------------------------------------------------------------------
205 // Initialize time integrator
206 if (problem->is_time_dependent())
207 {
208 if (init_time_integrator)
209 {
210 POLYFEM_SCOPED_TIMER("Initialize time integrator");
212
213 Eigen::MatrixXd solution, velocity, acceleration;
214 initial_solution(solution, ic_override); // Reload this because we need all previous solutions
215 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
216 assert(solution.rows() == sol.size());
217 initial_velocity(velocity, ic_override);
218 assert(velocity.rows() == sol.size());
219 initial_acceleration(acceleration, ic_override);
220 assert(acceleration.rows() == sol.size());
221
222 if (solution.cols() != velocity.cols() || solution.cols() != acceleration.cols())
223 {
225 "Incompatible initial-condition history for transient solve: "
226 "solution has {} columns, velocity has {}, acceleration has {}.",
227 solution.cols(), velocity.cols(), acceleration.cols());
228 }
229
230 const double dt = args["time"]["dt"];
231 solve_data.time_integrator->init(solution, velocity, acceleration, dt);
232 }
233 assert(solve_data.time_integrator != nullptr);
234 }
235 else
236 {
237 solve_data.time_integrator = nullptr;
238 }
239
240 // --------------------------------------------------------------------
241 // Initialize forms
242
243 damping_assembler = std::make_shared<assembler::ViscousDamping>();
245
247
248 // for backward solve
249 damping_prev_assembler = std::make_shared<assembler::ViscousDampingPrev>();
251
252 const ElementInversionCheck check_inversion = args["solver"]["advanced"]["check_inversion"];
253 const std::vector<std::shared_ptr<Form>> forms = solve_data.init_forms(
254 // General
255 units,
256 mesh->dimension(), t, in_node_to_node,
257 // Elastic form
258 n_bases, bases, geom_bases(), *assembler, ass_vals_cache, mass_ass_vals_cache, args["solver"]["advanced"]["jacobian_threshold"], check_inversion,
259 // Body form
262 n_boundary_samples(), rhs, sol, mass_matrix_assembler->density(),
263 // Pressure form
265 // Inertia form
266 args.value("/time/quasistatic"_json_pointer, true), mass,
267 damping_assembler->is_valid() ? damping_assembler : nullptr,
268 // Lagged regularization form
269 args["solver"]["advanced"]["lagged_regularization_weight"],
270 args["solver"]["advanced"]["lagged_regularization_iterations"],
271 // Augmented lagrangian form
272 obstacle.ndof(), args["constraints"]["hard"], args["constraints"]["soft"],
273 // Contact form
274 args["contact"]["enabled"], args["contact"]["periodic"].get<bool>() ? periodic_collision_mesh : collision_mesh, args["contact"]["dhat"],
275 avg_mass, args["contact"]["use_convergent_formulation"] ? bool(args["contact"]["use_area_weighting"]) : false,
276 args["contact"]["use_convergent_formulation"] ? bool(args["contact"]["use_improved_max_operator"]) : false,
277 args["contact"]["use_convergent_formulation"] ? bool(args["contact"]["use_physical_barrier"]) : false,
278 args["solver"]["contact"]["barrier_stiffness"],
279 args["solver"]["contact"]["initial_barrier_stiffness"],
280 args["solver"]["contact"]["CCD"]["broad_phase"],
281 args["solver"]["contact"]["CCD"]["tolerance"],
282 args["solver"]["contact"]["CCD"]["max_iterations"],
284 // Smooth Contact Form
285 args["contact"]["use_gcp_formulation"],
286 args["contact"]["alpha_t"],
287 args["contact"]["alpha_n"],
288 args["contact"]["use_adaptive_dhat"],
289 args["contact"]["min_distance_ratio"],
290 // Normal Adhesion Form
291 args["contact"]["adhesion"]["adhesion_enabled"],
292 args["contact"]["adhesion"]["dhat_p"],
293 args["contact"]["adhesion"]["dhat_a"],
294 args["contact"]["adhesion"]["adhesion_strength"],
295 // Tangential Adhesion Form
296 args["contact"]["adhesion"]["tangential_adhesion_coefficient"],
297 args["contact"]["adhesion"]["epsa"],
298 args["solver"]["contact"]["tangential_adhesion_iterations"],
299 // Homogenization
301 // Periodic contact
302 args["contact"]["periodic"], periodic_collision_mesh_to_basis, periodic_bc,
303 // Friction form
304 args["contact"]["friction_coefficient"],
305 args["contact"]["epsv"],
306 args["solver"]["contact"]["friction_iterations"],
307 // Rayleigh damping form
308 args["solver"]["rayleigh_damping"]);
309
310 for (const auto &form : forms)
311 form->set_output_dir(output_dir);
312
313 if (solve_data.contact_form != nullptr)
314 solve_data.contact_form->save_ccd_debug_meshes = args["output"]["advanced"]["save_ccd_debug_meshes"];
315
316 // --------------------------------------------------------------------
317 // Initialize nonlinear problems
318
319 const int ndof = n_bases * mesh->dimension();
320 solve_data.nl_problem = std::make_shared<NLProblem>(
322 polysolve::linear::Solver::create(args["solver"]["linear"], logger()));
323 solve_data.nl_problem->init(sol);
324 solve_data.nl_problem->update_quantities(t, sol);
325 // --------------------------------------------------------------------
326
327 stats.solver_info = json::array();
328 }
329
330 void State::solve_tensor_nonlinear(int step, Eigen::MatrixXd &sol, const bool init_lagging, UserPostStepCallback user_post_step)
331 {
332 assert(solve_data.nl_problem != nullptr);
333 NLProblem &nl_problem = *(solve_data.nl_problem);
334
335 assert(sol.size() == rhs.size());
336
337 if (nl_problem.uses_lagging())
338 {
339 if (init_lagging)
340 {
341 POLYFEM_SCOPED_TIMER("Initializing lagging");
342 nl_problem.init_lagging(sol); // TODO: this should be u_prev projected
343 }
344 logger().info("Lagging iteration 1:");
345 }
346
347 // ---------------------------------------------------------------------
348
349 // Save the subsolve sequence for debugging
350 int subsolve_count = 0;
351 save_subsolve(subsolve_count, step, sol, Eigen::MatrixXd()); // no pressure
352
353 // ---------------------------------------------------------------------
354
355 std::shared_ptr<polysolve::nonlinear::Solver> nl_solver = make_nl_solver(true);
356
357 ALSolver al_solver(
359 args["solver"]["augmented_lagrangian"]["initial_weight"],
360 args["solver"]["augmented_lagrangian"]["scaling"],
361 args["solver"]["augmented_lagrangian"]["max_weight"],
362 args["solver"]["augmented_lagrangian"]["eta"],
363 [&](const Eigen::VectorXd &x) {
364 this->solve_data.update_barrier_stiffness(sol);
365 });
366
367 al_solver.post_subsolve = [&](const double al_weight) {
368 stats.solver_info.push_back(
369 {{"type", al_weight > 0 ? "al" : "rc"},
370 {"t", step}, // TODO: null if static?
371 {"info", nl_solver->info()}});
372 if (al_weight > 0)
373 stats.solver_info.back()["weight"] = al_weight;
374 save_subsolve(++subsolve_count, step, sol, Eigen::MatrixXd()); // no pressure
375 };
376
377 Eigen::MatrixXd prev_sol = sol;
378 al_solver.solve_al(nl_problem, sol,
379 args["solver"]["augmented_lagrangian"]["nonlinear"], args["solver"]["linear"], units.characteristic_length());
380
381 al_solver.solve_reduced(nl_problem, sol,
382 args["solver"]["nonlinear"], args["solver"]["linear"], units.characteristic_length());
383
384 if (args["space"]["advanced"]["count_flipped_els_continuous"])
385 {
386 const auto invalidList = count_invalid(mesh->dimension(), bases, geom_bases(), sol);
387 logger().debug("Flipped elements (cnt {}) : {}", invalidList.size(), invalidList);
388 }
389
390 // ---------------------------------------------------------------------
391
392 // TODO: Make this more general
393 const double lagging_tol = args["solver"]["contact"].value("friction_convergence_tol", 1e-2) * units.characteristic_length();
394
396 {
397 // Lagging loop (start at 1 because we already did an iteration above)
398 bool lagging_converged = !nl_problem.uses_lagging();
399 for (int lag_i = 1; !lagging_converged; lag_i++)
400 {
401 Eigen::VectorXd tmp_sol = nl_problem.full_to_reduced(sol);
402
403 // Update the lagging before checking for convergence
404 nl_problem.update_lagging(tmp_sol, lag_i);
405
406 // Check if lagging converged
407 Eigen::VectorXd grad;
408 nl_problem.gradient(tmp_sol, grad);
409 const double delta_x_norm = (prev_sol - sol).lpNorm<Eigen::Infinity>();
410 logger().debug("Lagging convergence grad_norm={:g} tol={:g} (||Δx||={:g})", grad.norm(), lagging_tol, delta_x_norm);
411 if (grad.norm() <= lagging_tol)
412 {
413 logger().info(
414 "Lagging converged in {:d} iteration(s) (grad_norm={:g} tol={:g})",
415 lag_i, grad.norm(), lagging_tol);
416 lagging_converged = true;
417 break;
418 }
419
420 if (delta_x_norm <= 1e-12)
421 {
422 logger().warn(
423 "Lagging produced tiny update between iterations {:d} and {:d} (grad_norm={:g} grad_tol={:g} ||Δx||={:g} Δx_tol={:g}); stopping early",
424 lag_i - 1, lag_i, grad.norm(), lagging_tol, delta_x_norm, 1e-6);
425 lagging_converged = false;
426 break;
427 }
428
429 // Check for convergence first before checking if we can continue
430 if (lag_i >= nl_problem.max_lagging_iterations())
431 {
432 logger().warn(
433 "Lagging failed to converge with {:d} iteration(s) (grad_norm={:g} tol={:g})",
434 lag_i, grad.norm(), lagging_tol);
435 lagging_converged = false;
436 break;
437 }
438
439 // Solve the problem with the updated lagging
440 logger().info("Lagging iteration {:d}:", lag_i + 1);
441 nl_problem.init(sol);
443 nl_problem.normalize_forms();
444 nl_solver->minimize(nl_problem, tmp_sol);
445 nl_problem.finish();
446 prev_sol = sol;
447 sol = nl_problem.reduced_to_full(tmp_sol);
448
449 // Save the subsolve sequence for debugging and info
450 stats.solver_info.push_back(
451 {{"type", "rc"},
452 {"t", step}, // TODO: null if static?
453 {"lag_i", lag_i},
454 {"info", nl_solver->info()}});
455 save_subsolve(++subsolve_count, step, sol, Eigen::MatrixXd()); // no pressure
456 }
457 }
458
459 if (user_post_step)
460 {
461 user_post_step(step, *this, sol, nullptr, nullptr);
462 }
463 }
464} // namespace polyfem
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
Runtime override for initial-condition histories.
Definition State.hpp:90
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:805
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:499
int n_bases
number of bases
Definition State.hpp:212
int n_pressure_bases
number of pressure bases
Definition State.hpp:214
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:230
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:257
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:558
std::shared_ptr< assembler::ViscousDamping > damping_assembler
Definition State.hpp:198
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
Definition State.hpp:571
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:589
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:189
void initial_velocity(Eigen::MatrixXd &velocity, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial velocity.
ipc::CollisionMesh periodic_collision_mesh
IPC collision mesh under periodic BC.
Definition State.hpp:666
void save_subsolve(const int i, const int t, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a subsolve when save_solve_sequence_debug is true
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:663
std::string output_dir
Directory for output files.
Definition State.hpp:729
void set_materials(std::vector< std::shared_ptr< assembler::Assembler > > &assemblers) const
set the material and the problem dimension
void save_timestep(const double time, const int t, const double t0, const double dt, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a timestep
std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(bool for_al) const
factory to create the nl solver depending on input
void solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves transient tensor nonlinear problem
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:191
bool remesh(const double time, const double dt, Eigen::MatrixXd &sol)
Remesh the FE space and update solution(s).
void solve_tensor_nonlinear(int step, Eigen::MatrixXd &sol, const bool init_lagging=true, UserPostStepCallback user_post_step={})
solves nonlinear problems
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:587
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:236
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:202
json args
main input arguments containing all defaults
Definition State.hpp:135
void save_restart_json(const double t0, const double dt, const int t) const
Save a JSON sim file for restarting the simulation at time t.
void initial_acceleration(Eigen::MatrixXd &acceleration, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial acceleration.
io::OutStatsData stats
Other statistics.
Definition State.hpp:735
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:668
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:205
void initial_solution(Eigen::MatrixXd &solution, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial solution.
std::vector< mesh::LocalBoundary > local_boundary
mapping from elements to nodes for dirichlet boundary conditions
Definition State.hpp:554
double avg_mass
average system mass, used for contact with IPC
Definition State.hpp:238
bool optimization_enabled
Definition State.hpp:799
std::function< void(int, int, double, double)> time_callback
Definition State.hpp:740
QuadratureOrders n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:297
int ndof() const
Definition State.hpp:309
assembler::AssemblyValsCache mass_ass_vals_cache
Definition State.hpp:231
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:548
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:375
void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t=1.0, const bool init_time_integrator=true, const InitialConditionOverride *ic_override=nullptr)
initialize the nonlinear solver
std::shared_ptr< assembler::ViscousDampingPrev > damping_prev_assembler
Definition State.hpp:199
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:699
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:556
std::shared_ptr< assembler::PressureAssembler > build_pressure_assembler() const
Definition State.hpp:290
std::shared_ptr< assembler::PressureAssembler > elasticity_pressure_assembler
Definition State.hpp:196
void build_mesh_matrices(Eigen::MatrixXd &V, Eigen::MatrixXi &F)
Build the mesh matrices (vertices and elements) from the mesh using the bases node ordering.
bool is_problem_linear() const
Returns whether the system is linear. Collisions and pressure add nonlinearity to the problem.
Definition State.hpp:523
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:193
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:241
std::unordered_map< int, std::vector< mesh::LocalBoundary > > local_pressure_cavity
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:560
double characteristic_length() const
Definition Units.hpp:22
static void write(const std::string &path, const mesh::Mesh &mesh, const bool binary)
saves the mesh
Definition MshWriter.cpp:7
static bool write(const std::string &path, const Eigen::MatrixXd &v, const Eigen::MatrixXi &e, const Eigen::MatrixXi &f)
Definition OBJWriter.cpp:18
json solver_info
information of the solver, eg num iteration, time, errors, etc the informations varies depending on t...
Definition OutData.hpp:382
void solve_reduced(NLProblem &nl_problem, Eigen::MatrixXd &sol, std::shared_ptr< polysolve::nonlinear::Solver > nl_solver)
Definition ALSolver.hpp:41
std::function< void(const double)> post_subsolve
Definition ALSolver.hpp:53
void solve_al(NLProblem &nl_problem, Eigen::MatrixXd &sol, std::shared_ptr< polysolve::nonlinear::Solver > nl_solver)
Definition ALSolver.hpp:29
virtual void init(const TVector &x0) override
double normalize_forms() override
TVector full_to_reduced(const TVector &full) const
virtual void gradient(const TVector &x, TVector &gradv) override
void init_lagging(const TVector &x) override
TVector reduced_to_full(const TVector &reduced) const
void update_lagging(const TVector &x, const int iter_num) override
void update_dt()
updates the dt inside the different forms
std::shared_ptr< solver::NLProblem > nl_problem
std::shared_ptr< solver::ContactForm > contact_form
std::vector< std::shared_ptr< solver::AugmentedLagrangianForm > > al_form
std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator
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 QuadratureOrders &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
void update_barrier_stiffness(const Eigen::VectorXd &x)
update the barrier stiffness for the forms
static std::shared_ptr< ImplicitTimeIntegrator > construct_time_integrator(const json &params)
Factory method for constructing implicit time integrators from the name of the integrator.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
std::vector< int > count_invalid(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u)
Definition Jacobian.cpp:212
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
std::function< void(int step, State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd *disp_grad, const Eigen::MatrixXd *pressure)> UserPostStepCallback
User callback at the end of every solver step.
Definition State.hpp:86
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73