PolyFEM
Loading...
Searching...
No Matches
StateDiff.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
3
6
8#include <polysolve/linear/FEMSolver.hpp>
12
15
23
24#include <ipc/ipc.hpp>
25#include <ipc/barrier/barrier.hpp>
26#include <ipc/utils/local_to_global.hpp>
27#include <ipc/potentials/friction_potential.hpp>
28
29#include <Eigen/Dense>
30#include <unsupported/Eigen/SparseExtra>
31#include <deque>
32#include <map>
33#include <algorithm>
34
35#include <fstream>
36
37using namespace polyfem::basis;
38
39namespace polyfem
40{
41 namespace
42 {
43 void replace_rows_by_identity(StiffnessMatrix &reduced_mat, const StiffnessMatrix &mat, const std::vector<int> &rows)
44 {
45 reduced_mat.resize(mat.rows(), mat.cols());
46
47 std::vector<bool> mask(mat.rows(), false);
48 for (int i : rows)
49 mask[i] = true;
50
51 std::vector<Eigen::Triplet<double>> coeffs;
52 for (int k = 0; k < mat.outerSize(); ++k)
53 {
54 for (StiffnessMatrix::InnerIterator it(mat, k); it; ++it)
55 {
56 if (mask[it.row()])
57 {
58 if (it.row() == it.col())
59 coeffs.emplace_back(it.row(), it.col(), 1.0);
60 }
61 else
62 coeffs.emplace_back(it.row(), it.col(), it.value());
63 }
64 }
65 reduced_mat.setFromTriplets(coeffs.begin(), coeffs.end());
66 }
67 } // namespace
68
69 void State::get_vertices(Eigen::MatrixXd &vertices) const
70 {
71 vertices.setZero(mesh->n_vertices(), mesh->dimension());
72
73 for (int v = 0; v < mesh->n_vertices(); v++)
74 vertices.row(v) = mesh->point(v);
75 }
76
77 void State::get_elements(Eigen::MatrixXi &elements) const
78 {
79 assert(mesh->is_simplicial());
80
81 auto node_to_primitive_map = node_to_primitive();
82
83 const auto &gbases = geom_bases();
84 int dim = mesh->dimension();
85 elements.setZero(gbases.size(), dim + 1);
86 for (int e = 0; e < gbases.size(); e++)
87 {
88 int i = 0;
89 for (const auto &gbs : gbases[e].bases)
90 elements(e, i++) = node_to_primitive_map[gbs.global()[0].index];
91 }
92 }
93
94 void State::set_mesh_vertex(int v_id, const Eigen::VectorXd &vertex)
95 {
96 assert(vertex.size() == mesh->dimension());
97 mesh->set_point(v_id, vertex);
98 }
99
100 void State::cache_transient_adjoint_quantities(const int current_step, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad)
101 {
102 StiffnessMatrix gradu_h(sol.size(), sol.size());
103 if (current_step == 0)
104 diff_cached.init(mesh->dimension(), ndof(), problem->is_time_dependent() ? args["time"]["time_steps"].get<int>() : 0);
105
106 ipc::Collisions cur_collision_set;
107 ipc::FrictionCollisions cur_friction_set;
108
110 {
111 if (!problem->is_time_dependent() || current_step > 0)
112 compute_force_jacobian(sol, disp_grad, gradu_h);
113
114 cur_collision_set = solve_data.contact_form ? solve_data.contact_form->collision_set() : ipc::Collisions();
115 cur_friction_set = solve_data.friction_form ? solve_data.friction_form->friction_collision_set() : ipc::FrictionCollisions();
116 }
117 else
118 {
119 cur_collision_set = ipc::Collisions();
120 cur_friction_set = ipc::FrictionCollisions();
121 }
122
123 if (problem->is_time_dependent())
124 {
125 if (args["time"]["quasistatic"].get<bool>())
126 {
127 diff_cached.cache_quantities_quasistatic(current_step, sol, gradu_h, cur_collision_set, disp_grad);
128 }
129 else
130 {
131 Eigen::MatrixXd vel, acc;
132 if (current_step == 0)
133 {
134 if (dynamic_cast<time_integrator::BDF *>(solve_data.time_integrator.get()))
135 {
136 const auto bdf_integrator = dynamic_cast<time_integrator::BDF *>(solve_data.time_integrator.get());
137 vel = bdf_integrator->weighted_sum_v_prevs();
138 }
139 else if (dynamic_cast<time_integrator::ImplicitEuler *>(solve_data.time_integrator.get()))
140 {
141 const auto euler_integrator = dynamic_cast<time_integrator::ImplicitEuler *>(solve_data.time_integrator.get());
142 vel = euler_integrator->v_prev();
143 }
144 else
145 log_and_throw_error("Differentiable code doesn't support this time integrator!");
146
147 acc.setZero(ndof(), 1);
148 }
149 else
150 {
151 vel = solve_data.time_integrator->compute_velocity(sol);
152 acc = solve_data.time_integrator->compute_acceleration(vel);
153 }
154
155 diff_cached.cache_quantities_transient(current_step, solve_data.time_integrator->steps(), sol, vel, acc, gradu_h, cur_collision_set, cur_friction_set);
156 }
157 }
158 else
159 {
160 diff_cached.cache_quantities_static(sol, gradu_h, cur_collision_set, cur_friction_set, disp_grad);
161 }
162 }
163
164 void State::compute_force_jacobian(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad, StiffnessMatrix &hessian)
165 {
166 if (problem->is_time_dependent())
167 {
168 if (assembler->is_linear() && !is_contact_enabled())
169 log_and_throw_adjoint_error("Differentiable transient linear solve is not supported!");
170
171 StiffnessMatrix tmp_hess;
172 solve_data.nl_problem->set_project_to_psd(false);
173 solve_data.nl_problem->FullNLProblem::solution_changed(sol);
174 solve_data.nl_problem->FullNLProblem::hessian(sol, tmp_hess);
175 hessian.setZero();
176 replace_rows_by_identity(hessian, tmp_hess, boundary_nodes);
177 }
178 else // static formulation
179 {
180 if (assembler->is_linear() && !is_contact_enabled() && !is_homogenization())
181 {
182 hessian.setZero();
183 StiffnessMatrix stiffness;
184 build_stiffness_mat(stiffness);
185 replace_rows_by_identity(hessian, stiffness, boundary_nodes);
186 }
187 else
188 {
189 solve_data.nl_problem->set_project_to_psd(false);
190 if (is_homogenization())
191 {
192 Eigen::VectorXd reduced;
193 std::shared_ptr<solver::NLHomoProblem> homo_problem = std::dynamic_pointer_cast<solver::NLHomoProblem>(solve_data.nl_problem);
194 reduced = homo_problem->full_to_reduced(sol, disp_grad);
195 solve_data.nl_problem->solution_changed(reduced);
196 solve_data.nl_problem->hessian(reduced, hessian);
197 }
198 else
199 {
200 StiffnessMatrix tmp_hess;
201 solve_data.nl_problem->FullNLProblem::solution_changed(sol);
202 solve_data.nl_problem->FullNLProblem::hessian(sol, tmp_hess);
203 hessian.setZero();
204 replace_rows_by_identity(hessian, tmp_hess, boundary_nodes);
205 }
206 }
207 }
208 }
209
210 void State::compute_force_jacobian_prev(const int force_step, const int sol_step, StiffnessMatrix &hessian_prev) const
211 {
212 assert(force_step > 0);
213 assert(force_step > sol_step);
214 if (assembler->is_linear() && !is_contact_enabled())
215 {
216 hessian_prev = StiffnessMatrix(ndof(), ndof());
217 }
218 else
219 {
220 const Eigen::MatrixXd u = diff_cached.u(force_step);
221 const Eigen::MatrixXd u_prev = diff_cached.u(sol_step);
222 const double beta = time_integrator::BDF::betas(diff_cached.bdf_order(force_step) - 1);
223 const double dt = solve_data.time_integrator->dt();
224
225 hessian_prev = StiffnessMatrix(u.size(), u.size());
226 if (problem->is_time_dependent())
227 {
229 {
230 if (sol_step == force_step - 1)
231 {
232 Eigen::MatrixXd surface_solution_prev = collision_mesh.vertices(utils::unflatten(u_prev, mesh->dimension()));
233 Eigen::MatrixXd surface_solution = collision_mesh.vertices(utils::unflatten(u, mesh->dimension()));
234
235 // TODO: use the time integration to compute the velocity
236 const Eigen::MatrixXd surface_velocities = (surface_solution - surface_solution_prev) / dt;
237 const double dv_dut = -1 / dt;
238
239 hessian_prev =
240 solve_data.friction_form->friction_potential().force_jacobian(
243 collision_mesh.rest_positions(),
244 /*lagged_displacements=*/surface_solution_prev,
245 surface_velocities,
246 solve_data.contact_form->barrier_potential(),
247 solve_data.contact_form->barrier_stiffness(),
248 ipc::FrictionPotential::DiffWRT::LAGGED_DISPLACEMENTS)
249 + solve_data.friction_form->friction_potential().force_jacobian(
252 collision_mesh.rest_positions(),
253 /*lagged_displacements=*/surface_solution_prev,
254 surface_velocities,
255 solve_data.contact_form->barrier_potential(),
256 solve_data.contact_form->barrier_stiffness(),
257 ipc::FrictionPotential::DiffWRT::VELOCITIES)
258 * dv_dut;
259
260 hessian_prev *= -1;
261
262 // {
263 // Eigen::MatrixXd X = collision_mesh.rest_positions();
264 // Eigen::VectorXd x = utils::flatten(surface_solution_prev);
265 // const double barrier_stiffness = solve_data.contact_form->barrier_stiffness();
266 // const double dhat = solve_data.contact_form->dhat();
267 // const double mu = solve_data.friction_form->mu();
268 // const double epsv = solve_data.friction_form->epsv();
269
270 // Eigen::MatrixXd fgrad;
271 // fd::finite_jacobian(
272 // x, [&](const Eigen::VectorXd &y) -> Eigen::VectorXd
273 // {
274 // Eigen::MatrixXd fd_Ut = utils::unflatten(y, surface_solution_prev.cols());
275
276 // ipc::FrictionCollisions fd_friction_constraints;
277 // ipc::Collisions fd_constraints;
278 // fd_constraints.set_use_convergent_formulation(solve_data.contact_form->use_convergent_formulation());
279 // fd_constraints.set_are_shape_derivatives_enabled(true);
280 // fd_constraints.build(collision_mesh, X + fd_Ut, dhat);
281
282 // fd_friction_constraints.build(
283 // collision_mesh, X + fd_Ut, fd_constraints, dhat, barrier_stiffness,
284 // mu);
285
286 // return fd_friction_constraints.compute_potential_gradient(collision_mesh, (surface_solution - fd_Ut) / dt, epsv);
287
288 // }, fgrad, fd::AccuracyOrder::SECOND, 1e-8);
289
290 // std::cout << "force Ut derivative error " << (fgrad - hessian_prev).norm() << " " << hessian_prev.norm() << "\n";
291 // }
292
293 hessian_prev = collision_mesh.to_full_dof(hessian_prev); // / (beta * dt) / (beta * dt);
294 }
295 else
296 {
297 // const double alpha = time_integrator::BDF::alphas(std::min(diff_cached.bdf_order(force_step), force_step) - 1)[force_step - sol_step - 1];
298 // Eigen::MatrixXd velocity = collision_mesh.map_displacements(utils::unflatten(diff_cached.v(force_step), collision_mesh.dim()));
299 // hessian_prev = diff_cached.friction_collision_set(force_step).compute_potential_hessian( //
300 // collision_mesh, velocity, solve_data.friction_form->epsv(), false) * (-alpha / beta / dt);
301
302 // hessian_prev = collision_mesh.to_full_dof(hessian_prev);
303 }
304 }
305
306 if (damping_assembler->is_valid() && sol_step == force_step - 1) // velocity in damping uses BDF1
307 {
308 utils::SparseMatrixCache mat_cache;
309 StiffnessMatrix damping_hessian_prev(u.size(), u.size());
310 damping_prev_assembler->assemble_hessian(mesh->is_volume(), n_bases, false, bases, geom_bases(), ass_vals_cache, force_step * args["time"]["dt"].get<double>() + args["time"]["t0"].get<double>(), dt, u, u_prev, mat_cache, damping_hessian_prev);
311
312 hessian_prev += damping_hessian_prev;
313 }
314
315 if (sol_step == force_step - 1)
316 {
317 StiffnessMatrix body_force_hessian(u.size(), u.size());
318 solve_data.body_form->hessian_wrt_u_prev(u_prev, force_step * dt, body_force_hessian);
319 hessian_prev += body_force_hessian;
320 }
321 }
322 }
323 }
324
325 void State::solve_adjoint_cached(const Eigen::MatrixXd &rhs)
326 {
328 }
329
330 Eigen::MatrixXd State::solve_adjoint(const Eigen::MatrixXd &rhs) const
331 {
332 if (problem->is_time_dependent())
334 else
336 }
337
338 Eigen::MatrixXd State::solve_static_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
339 {
340 Eigen::MatrixXd b = adjoint_rhs;
341
342 Eigen::MatrixXd adjoint;
344 {
345 b(boundary_nodes, Eigen::all).setZero();
346
348 const int full_size = A.rows();
349 const int problem_dim = problem->is_scalar() ? 1 : mesh->dimension();
350 int precond_num = problem_dim * n_bases;
351
352 b.conservativeResizeLike(Eigen::MatrixXd::Zero(A.rows(), b.cols()));
353
354 std::vector<int> boundary_nodes_tmp;
355 if (has_periodic_bc())
356 {
357 boundary_nodes_tmp = periodic_bc->full_to_periodic(boundary_nodes);
358 precond_num = periodic_bc->full_to_periodic(A);
359 b = periodic_bc->full_to_periodic(b, true);
360 }
361 else
362 boundary_nodes_tmp = boundary_nodes;
363
364 adjoint.setZero(ndof(), adjoint_rhs.cols());
365 for (int i = 0; i < b.cols(); i++)
366 {
367 Eigen::VectorXd x, tmp;
368 tmp = b.col(i);
369 dirichlet_solve_prefactorized(*lin_solver_cached, A, tmp, boundary_nodes_tmp, x);
370
371 if (has_periodic_bc())
372 adjoint.col(i) = periodic_bc->periodic_to_full(full_size, x);
373 else
374 adjoint.col(i) = x;
375 }
376 }
377 else
378 {
379 auto solver = polysolve::linear::Solver::create(args["solver"]["adjoint_linear"], adjoint_logger());
380
381 StiffnessMatrix A = diff_cached.gradu_h(0); // This should be transposed, but A is symmetric in hyper-elastic and diffusion problems
382
383 /*
384 For non-periodic problems, the adjoint solution p's size is the full size in NLProblem
385 For periodic problems, the adjoint solution p's size is the reduced size in NLProblem
386 */
387 if (!is_homogenization())
388 {
389 adjoint.setZero(ndof(), adjoint_rhs.cols());
390 for (int i = 0; i < b.cols(); i++)
391 {
392 Eigen::VectorXd tmp = b.col(i);
393 tmp(boundary_nodes).setZero();
394
395 Eigen::VectorXd x;
396 x.setZero(tmp.size());
397 dirichlet_solve(*solver, A, tmp, boundary_nodes, x, A.rows(), "", false, false, false);
398
399 adjoint.col(i) = x;
400 }
401 // NLProblem sets dirichlet values to forward BC values, but we want zero in adjoint
402 adjoint(boundary_nodes, Eigen::all).setZero();
403 }
404 else
405 {
406 solver->analyze_pattern(A, A.rows());
407 solver->factorize(A);
408
409 adjoint.setZero(adjoint_rhs.rows(), adjoint_rhs.cols());
410 for (int i = 0; i < b.cols(); i++)
411 {
412 Eigen::MatrixXd tmp = b.col(i);
413
414 Eigen::VectorXd x;
415 x.setZero(tmp.size());
416 solver->solve(tmp, x);
417 x.conservativeResize(adjoint.rows());
418
419 adjoint.col(i) = x;
420 }
421 }
422 }
423
424 return adjoint;
425 }
426
427 Eigen::MatrixXd State::solve_transient_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
428 {
429 const double dt = args["time"]["dt"];
430 const int time_steps = args["time"]["time_steps"];
431
432 int bdf_order = 1;
433 if (args["time"]["integrator"].is_string())
434 bdf_order = 1;
435 else if (args["time"]["integrator"]["type"] == "ImplicitEuler")
436 bdf_order = 1;
437 else if (args["time"]["integrator"]["type"] == "BDF")
438 bdf_order = args["time"]["integrator"]["steps"].get<int>();
439 else
440 log_and_throw_adjoint_error("Integrator type not supported for differentiability.");
441
442 assert(adjoint_rhs.cols() == time_steps + 1);
443
444 const int cols_per_adjoint = time_steps + 1;
445 Eigen::MatrixXd adjoints;
446 adjoints.setZero(ndof(), cols_per_adjoint * 2);
447
448 // set dirichlet rows of mass to identity
449 StiffnessMatrix reduced_mass;
450 replace_rows_by_identity(reduced_mass, mass, boundary_nodes);
451
452 Eigen::MatrixXd sum_alpha_p, sum_alpha_nu;
453 for (int i = time_steps; i >= 0; --i)
454 {
455 {
456 sum_alpha_p.setZero(ndof(), 1);
457 sum_alpha_nu.setZero(ndof(), 1);
458
459 const int num = std::min(bdf_order, time_steps - i);
460
461 Eigen::VectorXd bdf_coeffs = Eigen::VectorXd::Zero(num);
462 for (int j = 0; j < bdf_order && i + j < time_steps; ++j)
463 bdf_coeffs(j) = -time_integrator::BDF::alphas(std::min(bdf_order - 1, i + j))[j];
464
465 sum_alpha_p = adjoints.middleCols(i + 1, num) * bdf_coeffs;
466 sum_alpha_nu = adjoints.middleCols(cols_per_adjoint + i + 1, num) * bdf_coeffs;
467 }
468
469 Eigen::VectorXd rhs_ = -reduced_mass.transpose() * sum_alpha_nu - adjoint_rhs.col(i);
470 for (int j = 1; j <= bdf_order; j++)
471 {
472 if (i + j > time_steps)
473 break;
474
475 StiffnessMatrix gradu_h_prev;
476 compute_force_jacobian_prev(i + j, i, gradu_h_prev);
477 Eigen::VectorXd tmp = adjoints.col(i + j) * (time_integrator::BDF::betas(diff_cached.bdf_order(i + j) - 1) * dt);
478 tmp(boundary_nodes).setZero();
479 rhs_ += -gradu_h_prev.transpose() * tmp;
480 }
481
482 if (i > 0)
483 {
484 double beta_dt = time_integrator::BDF::betas(diff_cached.bdf_order(i) - 1) * dt;
485
486 rhs_ += (1. / beta_dt) * (diff_cached.gradu_h(i) - reduced_mass).transpose() * sum_alpha_p;
487
488 {
489 StiffnessMatrix A = diff_cached.gradu_h(i).transpose();
490 Eigen::VectorXd b_ = rhs_;
491 b_(boundary_nodes).setZero();
492
493 auto solver = polysolve::linear::Solver::create(args["solver"]["adjoint_linear"], adjoint_logger());
494
495 Eigen::VectorXd x;
496 dirichlet_solve(*solver, A, b_, boundary_nodes, x, A.rows(), "", false, false, false);
497 adjoints.col(i + cols_per_adjoint) = x;
498 }
499
500 // TODO: generalize to BDFn
501 Eigen::VectorXd tmp = rhs_(boundary_nodes);
502 if (i + 1 < cols_per_adjoint)
503 tmp += (-2. / beta_dt) * adjoints(boundary_nodes, i + 1);
504 if (i + 2 < cols_per_adjoint)
505 tmp += (1. / beta_dt) * adjoints(boundary_nodes, i + 2);
506
507 tmp -= (diff_cached.gradu_h(i).transpose() * adjoints.col(i + cols_per_adjoint))(boundary_nodes);
508 adjoints(boundary_nodes, i + cols_per_adjoint) = tmp;
509 adjoints.col(i) = beta_dt * adjoints.col(i + cols_per_adjoint) - sum_alpha_p;
510 }
511 else
512 {
513 adjoints.col(i) = -reduced_mass.transpose() * sum_alpha_p;
514 adjoints.col(i + cols_per_adjoint) = rhs_; // adjoint_nu[0] actually stores adjoint_mu[0]
515 }
516 }
517 return adjoints;
518 }
519
520 void State::compute_surface_node_ids(const int surface_selection, std::vector<int> &node_ids) const
521 {
522 node_ids = {};
523
524 const auto &gbases = geom_bases();
525 for (const auto &lb : total_local_boundary)
526 {
527 const int e = lb.element_id();
528 for (int i = 0; i < lb.size(); ++i)
529 {
530 const int primitive_global_id = lb.global_primitive_id(i);
531 const int boundary_id = mesh->get_boundary_id(primitive_global_id);
532 const auto nodes = gbases[e].local_nodes_for_primitive(primitive_global_id, *mesh);
533
534 if (boundary_id == surface_selection)
535 {
536 for (long n = 0; n < nodes.size(); ++n)
537 {
538 const int g_id = gbases[e].bases[nodes(n)].global()[0].index;
539
540 if (std::count(node_ids.begin(), node_ids.end(), g_id) == 0)
541 node_ids.push_back(g_id);
542 }
543 }
544 }
545 }
546 }
547
548 void State::compute_total_surface_node_ids(std::vector<int> &node_ids) const
549 {
550 node_ids = {};
551
552 const auto &gbases = geom_bases();
553 for (const auto &lb : total_local_boundary)
554 {
555 const int e = lb.element_id();
556 for (int i = 0; i < lb.size(); ++i)
557 {
558 const int primitive_global_id = lb.global_primitive_id(i);
559 const auto nodes = gbases[e].local_nodes_for_primitive(primitive_global_id, *mesh);
560
561 for (long n = 0; n < nodes.size(); ++n)
562 {
563 const int g_id = gbases[e].bases[nodes(n)].global()[0].index;
564
565 if (std::count(node_ids.begin(), node_ids.end(), g_id) == 0)
566 node_ids.push_back(g_id);
567 }
568 }
569 }
570 }
571
572 void State::compute_volume_node_ids(const int volume_selection, std::vector<int> &node_ids) const
573 {
574 node_ids = {};
575
576 const auto &gbases = geom_bases();
577 for (int e = 0; e < gbases.size(); e++)
578 {
579 const int body_id = mesh->get_body_id(e);
580 if (body_id == volume_selection)
581 for (const auto &gbs : gbases[e].bases)
582 for (const auto &g : gbs.global())
583 node_ids.push_back(g.index);
584 }
585 }
586
587} // namespace polyfem
int x
void get_vertices(Eigen::MatrixXd &vertices) const
Definition StateDiff.cpp:69
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)
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
void set_mesh_vertex(int v_id, const Eigen::VectorXd &vertex)
Definition StateDiff.cpp:94
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
std::shared_ptr< assembler::ViscousDamping > damping_assembler
Definition State.hpp:164
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:515
solver::CacheLevel optimization_enabled
Definition State.hpp:647
bool is_homogenization() const
Definition State.hpp:711
Eigen::MatrixXd solve_static_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
void compute_force_jacobian(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad, StiffnessMatrix &hessian)
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
bool has_periodic_bc() const
Definition State.hpp:390
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
std::vector< int > node_to_primitive() const
Definition State.cpp:235
json args
main input arguments containing all defaults
Definition State.hpp:101
solver::DiffCache diff_cached
Definition State.hpp:649
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
void solve_adjoint_cached(const Eigen::MatrixXd &rhs)
Eigen::MatrixXd solve_adjoint(const Eigen::MatrixXd &rhs) const
void get_elements(Eigen::MatrixXi &elements) const
Definition StateDiff.cpp:77
void build_stiffness_mat(StiffnessMatrix &stiffness)
utility that builds the stiffness matrix and collects stats, used only for linear problems
void compute_volume_node_ids(const int volume_selection, std::vector< int > &node_ids) const
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
Definition State.hpp:431
int ndof() const
Definition State.hpp:653
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:427
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
std::unique_ptr< polysolve::linear::Solver > lin_solver_cached
Definition State.hpp:651
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
Eigen::MatrixXd solve_transient_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
void compute_force_jacobian_prev(const int force_step, const int sol_step, StiffnessMatrix &hessian_prev) const
void compute_total_surface_node_ids(std::vector< int > &node_ids) const
void compute_surface_node_ids(const int surface_selection, std::vector< int > &node_ids) const
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
const StiffnessMatrix & gradu_h(int step) const
int bdf_order(int step) const
void cache_quantities_quasistatic(const int cur_step, const Eigen::MatrixXd &u, const StiffnessMatrix &gradu_h, const ipc::Collisions &contact_set, const Eigen::MatrixXd &disp_grad)
Definition DiffCache.hpp:83
void init(const int dimension, const int ndof, const int n_time_steps=0)
Definition DiffCache.hpp:21
void cache_quantities_static(const Eigen::MatrixXd &u, const StiffnessMatrix &gradu_h, const ipc::Collisions &contact_set, const ipc::FrictionCollisions &friction_constraint_set, const Eigen::MatrixXd &disp_grad)
Definition DiffCache.hpp:40
void cache_adjoints(const Eigen::MatrixXd &adjoint_mat)
Definition DiffCache.hpp:98
void cache_quantities_transient(const int cur_step, const int cur_bdf_order, const Eigen::MatrixXd &u, const Eigen::MatrixXd &v, const Eigen::MatrixXd &acc, const StiffnessMatrix &gradu_h, const ipc::Collisions &collision_set, const ipc::FrictionCollisions &friction_collision_set)
Definition DiffCache.hpp:57
Eigen::VectorXd u(int step) const
const ipc::FrictionCollisions & friction_collision_set(int step) const
std::shared_ptr< solver::FrictionForm > friction_form
std::shared_ptr< solver::BodyForm > body_form
std::shared_ptr< solver::NLProblem > nl_problem
std::shared_ptr< solver::ContactForm > contact_form
std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator
Backward Differential Formulas.
Definition BDF.hpp:14
static double betas(const int i)
Retrieve the value of beta used for BDF with i steps.
Definition BDF.cpp:35
Eigen::VectorXd weighted_sum_v_prevs() const
Compute the weighted sum of the previous velocities.
Definition BDF.cpp:62
static const std::vector< double > & alphas(const int i)
Retrieve the alphas used for BDF with i steps.
Definition BDF.cpp:21
Implicit Euler time integrator of a second order ODE (equivently a system of coupled first order ODEs...
const Eigen::VectorXd & v_prev() const
Get the most recent previous velocity value.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:28
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22