PolyFEM
Loading...
Searching...
No Matches
StateSolveNonlinear.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
5
14
25
26#include <ipc/ipc.hpp>
27
28namespace polyfem
29{
30 using namespace mesh;
31 using namespace solver;
32 using namespace time_integrator;
33 using namespace io;
34 using namespace utils;
35
36 std::shared_ptr<polysolve::nonlinear::Solver> State::make_nl_solver(bool for_al) const
37 {
38 return polysolve::nonlinear::Solver::create(for_al ? args["solver"]["augmented_lagrangian"]["nonlinear"] : args["solver"]["nonlinear"], args["solver"]["linear"], units.characteristic_length(), logger());
39 }
40
41 void State::solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
42 {
43 init_nonlinear_tensor_solve(sol, t0 + dt);
44
45 // Write the total energy to a CSV file
46 int save_i = 0;
47 EnergyCSVWriter energy_csv(resolve_output_path("energy.csv"), solve_data);
48 RuntimeStatsCSVWriter stats_csv(resolve_output_path("stats.csv"), *this, t0, dt);
49 const bool remesh_enabled = args["space"]["remesh"]["enabled"];
50 // const double save_dt = remesh_enabled ? (dt / 3) : dt;
51
52 // Save the initial solution
53 energy_csv.write(save_i, sol);
54 save_timestep(t0, save_i, t0, dt, sol, Eigen::MatrixXd()); // no pressure
55 save_i++;
56
58 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
59
60 for (int t = 1; t <= time_steps; ++t)
61 {
62 double forward_solve_time = 0, remeshing_time = 0, global_relaxation_time = 0;
63
64 {
65 POLYFEM_SCOPED_TIMER(forward_solve_time);
67 }
68
69 if (remesh_enabled)
70 {
71 energy_csv.write(save_i, sol);
72 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
73 save_i++;
74
75 bool remesh_success;
76 {
77 POLYFEM_SCOPED_TIMER(remeshing_time);
78 remesh_success = this->remesh(t0 + dt * t, dt, sol);
79 }
80
81 // Save the solution after remeshing
82 energy_csv.write(save_i, sol);
83 // save_timestep(t0 + dt * t, save_i, t0, save_dt, sol, Eigen::MatrixXd()); // no pressure
84 save_i++;
85
86 // Only do global relaxation if remeshing was successful
87 if (remesh_success)
88 {
89 POLYFEM_SCOPED_TIMER(global_relaxation_time);
90 solve_tensor_nonlinear(sol, t, false); // solve the scene again after remeshing
91 }
92 }
93
94 // Always save the solution for consistency
95 energy_csv.write(save_i, sol);
96 save_timestep(t0 + dt * t, t, t0, dt, sol, Eigen::MatrixXd()); // no pressure
97 save_i++;
98
100 {
101 cache_transient_adjoint_quantities(t, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
102 }
103
104 {
105 POLYFEM_SCOPED_TIMER("Update quantities");
106
107 solve_data.time_integrator->update_quantities(sol);
108
109 solve_data.nl_problem->update_quantities(t0 + (t + 1) * dt, sol);
110
113 }
114
115 logger().info("{}/{} t={}", t, time_steps, t0 + dt * t);
116
117 const std::string rest_mesh_path = args["output"]["data"]["rest_mesh"].get<std::string>();
118 if (!rest_mesh_path.empty())
119 {
120 Eigen::MatrixXd V;
121 Eigen::MatrixXi F;
124 resolve_output_path(fmt::format(args["output"]["data"]["rest_mesh"], t)),
125 V, F, mesh->get_body_ids(), mesh->is_volume(), /*binary=*/true);
126 }
127
128 const std::string &state_path = resolve_output_path(fmt::format(args["output"]["data"]["state"], t));
129 if (!state_path.empty())
130 solve_data.time_integrator->save_state(state_path);
131
132 // save restart file
133 save_restart_json(t0, dt, t);
134 if (remesh_enabled)
135 stats_csv.write(t, forward_solve_time, remeshing_time, global_relaxation_time, sol);
136 }
137 }
138
139 void State::init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t, const bool init_time_integrator)
140 {
141 assert(sol.cols() == 1);
142 assert(!assembler->is_linear() || is_contact_enabled()); // non-linear
143 assert(!problem->is_scalar()); // tensor
144 assert(mixed_assembler == nullptr);
145
147 {
148 if (initial_sol_update.size() == ndof())
149 sol = initial_sol_update;
150 else
151 initial_sol_update = sol;
152 }
153
154 // --------------------------------------------------------------------
155 // Check for initial intersections
156 if (is_contact_enabled())
157 {
158 POLYFEM_SCOPED_TIMER("Check for initial intersections");
159
160 const Eigen::MatrixXd displaced = collision_mesh.displace_vertices(
161 utils::unflatten(sol, mesh->dimension()));
162
163 if (ipc::has_intersections(collision_mesh, displaced, args["solver"]["contact"]["CCD"]["broad_phase"]))
164 {
166 resolve_output_path("intersection.obj"), displaced,
167 collision_mesh.edges(), collision_mesh.faces());
168 log_and_throw_error("Unable to solve, initial solution has intersections!");
169 }
170 }
171
172 // --------------------------------------------------------------------
173 // Initialize time integrator
174 if (problem->is_time_dependent())
175 {
176 if (init_time_integrator)
177 {
178 POLYFEM_SCOPED_TIMER("Initialize time integrator");
180
181 Eigen::MatrixXd solution, velocity, acceleration;
182 initial_solution(solution); // Reload this because we need all previous solutions
183 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
184 assert(solution.rows() == sol.size());
185 initial_velocity(velocity);
186 assert(velocity.rows() == sol.size());
187 initial_acceleration(acceleration);
188 assert(acceleration.rows() == sol.size());
189
191 {
192 if (initial_vel_update.size() == ndof())
193 velocity = initial_vel_update;
194 else
195 initial_vel_update = velocity;
196 }
197
198 const double dt = args["time"]["dt"];
199 solve_data.time_integrator->init(solution, velocity, acceleration, dt);
200 }
201 assert(solve_data.time_integrator != nullptr);
202 }
203 else
204 {
205 solve_data.time_integrator = nullptr;
206 }
207
208 // --------------------------------------------------------------------
209 // Initialize forms
210
211 damping_assembler = std::make_shared<assembler::ViscousDamping>();
213
215
216 // for backward solve
217 damping_prev_assembler = std::make_shared<assembler::ViscousDampingPrev>();
219
220 const std::vector<std::shared_ptr<Form>> forms = solve_data.init_forms(
221 // General
222 units,
223 mesh->dimension(), t,
224 // Elastic form
226 // Body form
229 n_boundary_samples(), rhs, sol, mass_matrix_assembler->density(),
230 // Pressure form
232 // Inertia form
233 args.value("/time/quasistatic"_json_pointer, true), mass,
234 damping_assembler->is_valid() ? damping_assembler : nullptr,
235 // Lagged regularization form
236 args["solver"]["advanced"]["lagged_regularization_weight"],
237 args["solver"]["advanced"]["lagged_regularization_iterations"],
238 // Augmented lagrangian form
239 obstacle.ndof(),
240 // Contact form
241 args["contact"]["enabled"], args["contact"]["periodic"].get<bool>() ? periodic_collision_mesh : collision_mesh, args["contact"]["dhat"],
242 avg_mass, args["contact"]["use_convergent_formulation"],
243 args["solver"]["contact"]["barrier_stiffness"],
244 args["solver"]["contact"]["CCD"]["broad_phase"],
245 args["solver"]["contact"]["CCD"]["tolerance"],
246 args["solver"]["contact"]["CCD"]["max_iterations"],
248 // Homogenization
250 // Periodic contact
251 args["contact"]["periodic"], periodic_collision_mesh_to_basis,
252 // Friction form
253 args["contact"]["friction_coefficient"],
254 args["contact"]["epsv"],
255 args["solver"]["contact"]["friction_iterations"],
256 // Rayleigh damping form
257 args["solver"]["rayleigh_damping"]);
258
259 for (const auto &form : forms)
260 form->set_output_dir(output_dir);
261
262 if (solve_data.contact_form != nullptr)
263 solve_data.contact_form->save_ccd_debug_meshes = args["output"]["advanced"]["save_ccd_debug_meshes"];
264
265 // --------------------------------------------------------------------
266 // Initialize nonlinear problems
267
268 const int ndof = n_bases * mesh->dimension();
269 solve_data.nl_problem = std::make_shared<NLProblem>(
272 solve_data.nl_problem->init(sol);
273 solve_data.nl_problem->update_quantities(t, sol);
274 // --------------------------------------------------------------------
275
276 stats.solver_info = json::array();
277 }
278
279 void State::solve_tensor_nonlinear(Eigen::MatrixXd &sol, const int t, const bool init_lagging)
280 {
281 assert(solve_data.nl_problem != nullptr);
282 NLProblem &nl_problem = *(solve_data.nl_problem);
283
284 assert(sol.size() == rhs.size());
285
286 if (nl_problem.uses_lagging())
287 {
288 if (init_lagging)
289 {
290 POLYFEM_SCOPED_TIMER("Initializing lagging");
291 nl_problem.init_lagging(sol); // TODO: this should be u_prev projected
292 }
293 logger().info("Lagging iteration 1:");
294 }
295
296 // ---------------------------------------------------------------------
297
298 // Save the subsolve sequence for debugging
299 int subsolve_count = 0;
300 save_subsolve(subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
301
302 // ---------------------------------------------------------------------
303
304 std::shared_ptr<polysolve::nonlinear::Solver> nl_solver = make_nl_solver(true);
305
306 ALSolver al_solver(
308 args["solver"]["augmented_lagrangian"]["initial_weight"],
309 args["solver"]["augmented_lagrangian"]["scaling"],
310 args["solver"]["augmented_lagrangian"]["max_weight"],
311 args["solver"]["augmented_lagrangian"]["eta"],
312 [&](const Eigen::VectorXd &x) {
313 this->solve_data.update_barrier_stiffness(sol);
314 });
315
316 al_solver.post_subsolve = [&](const double al_weight) {
317 stats.solver_info.push_back(
318 {{"type", al_weight > 0 ? "al" : "rc"},
319 {"t", t}, // TODO: null if static?
320 {"info", nl_solver->info()}});
321 if (al_weight > 0)
322 stats.solver_info.back()["weight"] = al_weight;
323 save_subsolve(++subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
324 };
325
326 Eigen::MatrixXd prev_sol = sol;
327 al_solver.solve_al(nl_solver, nl_problem, sol);
328
329 nl_solver = make_nl_solver(false);
330 al_solver.solve_reduced(nl_solver, nl_problem, sol);
331
332 // ---------------------------------------------------------------------
333
334 // TODO: Make this more general
335 const double lagging_tol = args["solver"]["contact"].value("friction_convergence_tol", 1e-2) * units.characteristic_length();
336
338 {
339 // Lagging loop (start at 1 because we already did an iteration above)
340 bool lagging_converged = !nl_problem.uses_lagging();
341 for (int lag_i = 1; !lagging_converged; lag_i++)
342 {
343 Eigen::VectorXd tmp_sol = nl_problem.full_to_reduced(sol);
344
345 // Update the lagging before checking for convergence
346 nl_problem.update_lagging(tmp_sol, lag_i);
347
348 // Check if lagging converged
349 Eigen::VectorXd grad;
350 nl_problem.gradient(tmp_sol, grad);
351 const double delta_x_norm = (prev_sol - sol).lpNorm<Eigen::Infinity>();
352 logger().debug("Lagging convergence grad_norm={:g} tol={:g} (||Δx||={:g})", grad.norm(), lagging_tol, delta_x_norm);
353 if (grad.norm() <= lagging_tol)
354 {
355 logger().info(
356 "Lagging converged in {:d} iteration(s) (grad_norm={:g} tol={:g})",
357 lag_i, grad.norm(), lagging_tol);
358 lagging_converged = true;
359 break;
360 }
361
362 if (delta_x_norm <= 1e-12)
363 {
364 logger().warn(
365 "Lagging produced tiny update between iterations {:d} and {:d} (grad_norm={:g} grad_tol={:g} ||Δx||={:g} Δx_tol={:g}); stopping early",
366 lag_i - 1, lag_i, grad.norm(), lagging_tol, delta_x_norm, 1e-6);
367 lagging_converged = false;
368 break;
369 }
370
371 // Check for convergence first before checking if we can continue
372 if (lag_i >= nl_problem.max_lagging_iterations())
373 {
374 logger().warn(
375 "Lagging failed to converge with {:d} iteration(s) (grad_norm={:g} tol={:g})",
376 lag_i, grad.norm(), lagging_tol);
377 lagging_converged = false;
378 break;
379 }
380
381 // Solve the problem with the updated lagging
382 logger().info("Lagging iteration {:d}:", lag_i + 1);
383 nl_problem.init(sol);
385 nl_solver->minimize(nl_problem, tmp_sol);
386 prev_sol = sol;
387 sol = nl_problem.reduced_to_full(tmp_sol);
388
389 // Save the subsolve sequence for debugging and info
390 stats.solver_info.push_back(
391 {{"type", "rc"},
392 {"t", t}, // TODO: null if static?
393 {"lag_i", lag_i},
394 {"info", nl_solver->info()}});
395 save_subsolve(++subsolve_count, t, sol, Eigen::MatrixXd()); // no pressure
396 }
397 }
398 }
399} // namespace polyfem
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
Eigen::MatrixXd initial_vel_update
Definition State.hpp:697
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:705
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:518
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:515
std::string output_dir
Directory for output files.
Definition State.hpp:573
solver::CacheLevel optimization_enabled
Definition State.hpp:647
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
bool remesh(const double time, const double dt, Eigen::MatrixXd &sol)
Remesh the FE space and update solution(s).
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:585
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:520
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:653
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:551
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:697
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:2933
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:2959
std::function< void(const double)> post_subsolve
Definition ALSolver.hpp:34
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