PolyFEM
Loading...
Searching...
No Matches
StateSolveNonlinear.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
5
14
24
25#include <ipc/ipc.hpp>
26
27namespace polyfem
28{
29 using namespace mesh;
30 using namespace solver;
31 using namespace time_integrator;
32 using namespace io;
33 using namespace utils;
34
35 std::shared_ptr<polysolve::nonlinear::Solver> State::make_nl_solver(bool for_al) const
36 {
37 return polysolve::nonlinear::Solver::create(for_al ? args["solver"]["augmented_lagrangian"]["nonlinear"] : args["solver"]["nonlinear"], args["solver"]["linear"], units.characteristic_length(), logger());
38 }
39
40 void State::solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
41 {
42 init_nonlinear_tensor_solve(sol, t0 + dt);
43
44 // Write the total energy to a CSV file
45 int save_i = 0;
46 EnergyCSVWriter energy_csv(resolve_output_path("energy.csv"), solve_data);
47 RuntimeStatsCSVWriter stats_csv(resolve_output_path("stats.csv"), *this, t0, dt);
48 const bool remesh_enabled = args["space"]["remesh"]["enabled"];
49#ifndef POLYFEM_WITH_REMESHING
50 if (remesh_enabled)
51 log_and_throw_error("Remeshing is not enabled in this build! Set POLYFEM_WITH_REMESHING=ON in CMake to enable it.");
52#endif
53 // const double save_dt = remesh_enabled ? (dt / 3) : dt;
54
55 // Save the initial solution
56 energy_csv.write(save_i, sol);
57 save_timestep(t0, save_i, t0, dt, sol, Eigen::MatrixXd()); // no pressure
58 save_i++;
59
61 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
62
63 for (int t = 1; t <= time_steps; ++t)
64 {
65 double forward_solve_time = 0, remeshing_time = 0, global_relaxation_time = 0;
66
67 {
68 POLYFEM_SCOPED_TIMER(forward_solve_time);
70 }
71
72#ifdef POLYFEM_WITH_REMESHING
73 if (remesh_enabled)
74 {
75 energy_csv.write(save_i, sol);
76 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
77 save_i++;
78
79 bool remesh_success;
80 {
81 POLYFEM_SCOPED_TIMER(remeshing_time);
82 remesh_success = this->remesh(t0 + dt * t, dt, sol);
83 }
84
85 // Save the solution after remeshing
86 energy_csv.write(save_i, sol);
87 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
88 save_i++;
89
90 // Only do global relaxation if remeshing was successful
91 if (remesh_success)
92 {
93 POLYFEM_SCOPED_TIMER(global_relaxation_time);
94 solve_tensor_nonlinear(sol, t, false); // solve the scene again after remeshing
95 }
96 }
97#endif
98 // Always save the solution for consistency
99 energy_csv.write(save_i, sol);
100 save_timestep(t0 + dt * t, t, t0, dt, sol, Eigen::MatrixXd()); // no pressure
101 save_i++;
102
104 {
105 cache_transient_adjoint_quantities(t, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
106 }
107
108 {
109 POLYFEM_SCOPED_TIMER("Update quantities");
110
111 solve_data.time_integrator->update_quantities(sol);
112
113 solve_data.nl_problem->update_quantities(t0 + (t + 1) * dt, sol);
114
117 }
118
119 logger().info("{}/{} t={}", t, time_steps, t0 + dt * t);
120
121 const std::string rest_mesh_path = args["output"]["data"]["rest_mesh"].get<std::string>();
122 if (!rest_mesh_path.empty())
123 {
124 Eigen::MatrixXd V;
125 Eigen::MatrixXi F;
128 resolve_output_path(fmt::format(args["output"]["data"]["rest_mesh"], t)),
129 V, F, mesh->get_body_ids(), mesh->is_volume(), /*binary=*/true);
130 }
131
132 const std::string &state_path = resolve_output_path(fmt::format(args["output"]["data"]["state"], t));
133 if (!state_path.empty())
134 solve_data.time_integrator->save_state(state_path);
135
136 // save restart file
137 save_restart_json(t0, dt, t);
138 stats_csv.write(t, forward_solve_time, remeshing_time, global_relaxation_time, sol);
139 }
140 }
141
142 void State::init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t, const bool init_time_integrator)
143 {
144 assert(sol.cols() == 1);
145 assert(!assembler->is_linear() || is_contact_enabled()); // non-linear
146 assert(!problem->is_scalar()); // tensor
147 assert(mixed_assembler == nullptr);
148
150 {
151 if (initial_sol_update.size() == ndof())
152 sol = initial_sol_update;
153 else
154 initial_sol_update = sol;
155 }
156
157 // --------------------------------------------------------------------
158 // Check for initial intersections
159 if (is_contact_enabled())
160 {
161 POLYFEM_SCOPED_TIMER("Check for initial intersections");
162
163 const Eigen::MatrixXd displaced = collision_mesh.displace_vertices(
164 utils::unflatten(sol, mesh->dimension()));
165
166 if (ipc::has_intersections(collision_mesh, displaced, args["solver"]["contact"]["CCD"]["broad_phase"]))
167 {
169 resolve_output_path("intersection.obj"), displaced,
170 collision_mesh.edges(), collision_mesh.faces());
171 log_and_throw_error("Unable to solve, initial solution has intersections!");
172 }
173 }
174
175 // --------------------------------------------------------------------
176 // Initialize time integrator
177 if (problem->is_time_dependent())
178 {
179 if (init_time_integrator)
180 {
181 POLYFEM_SCOPED_TIMER("Initialize time integrator");
183
184 Eigen::MatrixXd solution, velocity, acceleration;
185 initial_solution(solution); // Reload this because we need all previous solutions
186 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
187 assert(solution.rows() == sol.size());
188 initial_velocity(velocity);
189 assert(velocity.rows() == sol.size());
190 initial_acceleration(acceleration);
191 assert(acceleration.rows() == sol.size());
192
194 {
195 if (initial_vel_update.size() == ndof())
196 velocity = initial_vel_update;
197 else
198 initial_vel_update = velocity;
199 }
200
201 const double dt = args["time"]["dt"];
202 solve_data.time_integrator->init(solution, velocity, acceleration, dt);
203 }
204 assert(solve_data.time_integrator != nullptr);
205 }
206 else
207 {
208 solve_data.time_integrator = nullptr;
209 }
210
211 // --------------------------------------------------------------------
212 // Initialize forms
213
214 damping_assembler = std::make_shared<assembler::ViscousDamping>();
216
218
219 // for backward solve
220 damping_prev_assembler = std::make_shared<assembler::ViscousDampingPrev>();
222
223 const std::vector<std::shared_ptr<Form>> forms = solve_data.init_forms(
224 // General
225 units,
226 mesh->dimension(), t,
227 // Elastic form
229 // Body form
232 n_boundary_samples(), rhs, sol, mass_matrix_assembler->density(),
233 // Pressure form
235 // Inertia form
236 args.value("/time/quasistatic"_json_pointer, true), mass,
237 damping_assembler->is_valid() ? damping_assembler : nullptr,
238 // Lagged regularization form
239 args["solver"]["advanced"]["lagged_regularization_weight"],
240 args["solver"]["advanced"]["lagged_regularization_iterations"],
241 // Augmented lagrangian form
242 obstacle.ndof(),
243 // Contact form
244 args["contact"]["enabled"], args["contact"]["periodic"].get<bool>() ? periodic_collision_mesh : collision_mesh, args["contact"]["dhat"],
245 avg_mass, args["contact"]["use_convergent_formulation"],
246 args["solver"]["contact"]["barrier_stiffness"],
247 args["solver"]["contact"]["CCD"]["broad_phase"],
248 args["solver"]["contact"]["CCD"]["tolerance"],
249 args["solver"]["contact"]["CCD"]["max_iterations"],
251 // Homogenization
253 // Periodic contact
254 args["contact"]["periodic"], periodic_collision_mesh_to_basis,
255 // Friction form
256 args["contact"]["friction_coefficient"],
257 args["contact"]["epsv"],
258 args["solver"]["contact"]["friction_iterations"],
259 // Rayleigh damping form
260 args["solver"]["rayleigh_damping"]);
261
262 for (const auto &form : forms)
263 form->set_output_dir(output_dir);
264
265 if (solve_data.contact_form != nullptr)
266 solve_data.contact_form->save_ccd_debug_meshes = args["output"]["advanced"]["save_ccd_debug_meshes"];
267
268 // --------------------------------------------------------------------
269 // Initialize nonlinear problems
270
271 const int ndof = n_bases * mesh->dimension();
272 solve_data.nl_problem = std::make_shared<NLProblem>(
275 solve_data.nl_problem->init(sol);
276 solve_data.nl_problem->update_quantities(t, sol);
277 // --------------------------------------------------------------------
278
279 stats.solver_info = json::array();
280 }
281
282 void State::solve_tensor_nonlinear(Eigen::MatrixXd &sol, const int t, const bool init_lagging)
283 {
284 assert(solve_data.nl_problem != nullptr);
285 NLProblem &nl_problem = *(solve_data.nl_problem);
286
287 assert(sol.size() == rhs.size());
288
289 if (nl_problem.uses_lagging())
290 {
291 if (init_lagging)
292 {
293 POLYFEM_SCOPED_TIMER("Initializing lagging");
294 nl_problem.init_lagging(sol); // TODO: this should be u_prev projected
295 }
296 logger().info("Lagging iteration 1:");
297 }
298
299 // ---------------------------------------------------------------------
300
301 // Save the subsolve sequence for debugging
302 int subsolve_count = 0;
303 save_subsolve(subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
304
305 // ---------------------------------------------------------------------
306
307 std::shared_ptr<polysolve::nonlinear::Solver> nl_solver = make_nl_solver(true);
308
309 ALSolver al_solver(
311 args["solver"]["augmented_lagrangian"]["initial_weight"],
312 args["solver"]["augmented_lagrangian"]["scaling"],
313 args["solver"]["augmented_lagrangian"]["max_weight"],
314 args["solver"]["augmented_lagrangian"]["eta"],
315 [&](const Eigen::VectorXd &x) {
316 this->solve_data.update_barrier_stiffness(sol);
317 });
318
319 al_solver.post_subsolve = [&](const double al_weight) {
320 stats.solver_info.push_back(
321 {{"type", al_weight > 0 ? "al" : "rc"},
322 {"t", t}, // TODO: null if static?
323 {"info", nl_solver->info()}});
324 if (al_weight > 0)
325 stats.solver_info.back()["weight"] = al_weight;
326 save_subsolve(++subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
327 };
328
329 Eigen::MatrixXd prev_sol = sol;
330 al_solver.solve_al(nl_solver, nl_problem, sol);
331
332 nl_solver = make_nl_solver(false);
333 al_solver.solve_reduced(nl_solver, nl_problem, sol);
334
335 // ---------------------------------------------------------------------
336
337 // TODO: Make this more general
338 const double lagging_tol = args["solver"]["contact"].value("friction_convergence_tol", 1e-2) * units.characteristic_length();
339
341 {
342 // Lagging loop (start at 1 because we already did an iteration above)
343 bool lagging_converged = !nl_problem.uses_lagging();
344 for (int lag_i = 1; !lagging_converged; lag_i++)
345 {
346 Eigen::VectorXd tmp_sol = nl_problem.full_to_reduced(sol);
347
348 // Update the lagging before checking for convergence
349 nl_problem.update_lagging(tmp_sol, lag_i);
350
351 // Check if lagging converged
352 Eigen::VectorXd grad;
353 nl_problem.gradient(tmp_sol, grad);
354 const double delta_x_norm = (prev_sol - sol).lpNorm<Eigen::Infinity>();
355 logger().debug("Lagging convergence grad_norm={:g} tol={:g} (||Δx||={:g})", grad.norm(), lagging_tol, delta_x_norm);
356 if (grad.norm() <= lagging_tol)
357 {
358 logger().info(
359 "Lagging converged in {:d} iteration(s) (grad_norm={:g} tol={:g})",
360 lag_i, grad.norm(), lagging_tol);
361 lagging_converged = true;
362 break;
363 }
364
365 if (delta_x_norm <= 1e-12)
366 {
367 logger().warn(
368 "Lagging produced tiny update between iterations {:d} and {:d} (grad_norm={:g} grad_tol={:g} ||Δx||={:g} Δx_tol={:g}); stopping early",
369 lag_i - 1, lag_i, grad.norm(), lagging_tol, delta_x_norm, 1e-6);
370 lagging_converged = false;
371 break;
372 }
373
374 // Check for convergence first before checking if we can continue
375 if (lag_i >= nl_problem.max_lagging_iterations())
376 {
377 logger().warn(
378 "Lagging failed to converge with {:d} iteration(s) (grad_norm={:g} tol={:g})",
379 lag_i, grad.norm(), lagging_tol);
380 lagging_converged = false;
381 break;
382 }
383
384 // Solve the problem with the updated lagging
385 logger().info("Lagging iteration {:d}:", lag_i + 1);
386 nl_problem.init(sol);
388 nl_solver->minimize(nl_problem, tmp_sol);
389 prev_sol = sol;
390 sol = nl_problem.reduced_to_full(tmp_sol);
391
392 // Save the subsolve sequence for debugging and info
393 stats.solver_info.push_back(
394 {{"type", "rc"},
395 {"t", t}, // TODO: null if static?
396 {"lag_i", lag_i},
397 {"info", nl_solver->info()}});
398 save_subsolve(++subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
399 }
400 }
401 }
402} // namespace polyfem
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
Eigen::MatrixXd initial_vel_update
Definition State.hpp:717
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:725
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:389
int n_bases
number of bases
Definition State.hpp:178
void cache_transient_adjoint_quantities(const int current_step, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad)
int n_pressure_bases
number of pressure bases
Definition State.hpp:180
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
int n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:263
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:437
std::shared_ptr< assembler::ViscousDamping > damping_assembler
Definition State.hpp:164
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:468
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
ipc::CollisionMesh periodic_collision_mesh
IPC collision mesh under periodic BC.
Definition State.hpp:541
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:538
std::string output_dir
Directory for output files.
Definition State.hpp:593
solver::CacheLevel optimization_enabled
Definition State.hpp:667
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
void solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
solves transient tensor nonlinear problem
void initial_solution(Eigen::MatrixXd &solution) const
Load or compute the initial solution.
std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(bool for_al) const
factory to create the nl solver depending on input
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:202
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
json args
main input arguments containing all defaults
Definition State.hpp:101
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_velocity(Eigen::MatrixXd &velocity) const
Load or compute the initial velocity.
io::OutStatsData stats
Other statistics.
Definition State.hpp:605
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:543
void initial_acceleration(Eigen::MatrixXd &acceleration) const
Load or compute the initial acceleration.
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
std::vector< mesh::LocalBoundary > local_boundary
mapping from elements to nodes for dirichlet boundary conditions
Definition State.hpp:433
double avg_mass
average system mass, used for contact with IPC
Definition State.hpp:204
int ndof() const
Definition State.hpp:673
assembler::AssemblyValsCache mass_ass_vals_cache
Definition State.hpp:197
void solve_tensor_nonlinear(Eigen::MatrixXd &sol, const int t=0, const bool init_lagging=true)
solves nonlinear problems
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:427
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t=1.0, const bool init_time_integrator=true)
initialize the nonlinear solver
std::shared_ptr< assembler::ViscousDampingPrev > damping_prev_assembler
Definition State.hpp:165
bool is_contact_enabled() const
does the simulation has contact
Definition State.hpp:574
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:435
std::shared_ptr< assembler::PressureAssembler > build_pressure_assembler() const
Definition State.hpp:256
std::shared_ptr< assembler::PressureAssembler > elasticity_pressure_assembler
Definition State.hpp:162
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.
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:159
Eigen::MatrixXd initial_sol_update
Definition State.hpp:717
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
std::unordered_map< int, std::vector< mesh::LocalBoundary > > local_pressure_cavity
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:439
double characteristic_length() const
Definition Units.hpp:22
void write(const int i, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2927
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:398
void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2953
std::function< void(const double)> post_subsolve
Definition ALSolver.hpp:33
void solve_reduced(std::shared_ptr< NLSolver > nl_solver, NLProblem &nl_problem, Eigen::MatrixXd &sol)
Definition ALSolver.cpp:93
void solve_al(std::shared_ptr< NLSolver > nl_solver, NLProblem &nl_problem, Eigen::MatrixXd &sol)
Definition ALSolver.cpp:25
virtual void init(const TVector &x0) override
virtual TVector full_to_reduced(const TVector &full) const
virtual void gradient(const TVector &x, TVector &gradv) override
void init_lagging(const TVector &x) override
Definition NLProblem.cpp:62
virtual TVector reduced_to_full(const TVector &reduced) const
void update_lagging(const TVector &x, const int iter_num) override
Definition NLProblem.cpp:67
void update_dt()
updates the dt inside the different forms
std::shared_ptr< solver::BCPenaltyForm > al_pen_form
std::shared_ptr< solver::NLProblem > nl_problem
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 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::RhsAssembler > rhs_assembler
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.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71