12#include <polysolve/linear/FEMSolver.hpp>
28#include <ipc/potentials/friction_potential.hpp>
45 reduced_mat.resize(mat.rows(), mat.cols());
47 std::vector<bool> mask(mat.rows(),
false);
51 std::vector<Eigen::Triplet<double>> coeffs;
52 for (
int k = 0; k < mat.outerSize(); ++k)
54 for (StiffnessMatrix::InnerIterator it(mat, k); it; ++it)
58 if (it.row() == it.col())
59 coeffs.emplace_back(it.row(), it.col(), 1.0);
62 coeffs.emplace_back(it.row(), it.col(), it.value());
65 reduced_mat.setFromTriplets(coeffs.begin(), coeffs.end());
68 void compute_force_jacobian_prev(
const State &state,
const DiffCache &diff_cache,
const int force_step,
const int sol_step,
StiffnessMatrix &hessian_prev)
70 assert(force_step > 0);
71 assert(force_step > sol_step);
75 if (s.assembler->is_linear() && !s.is_contact_enabled())
81 const Eigen::MatrixXd u = diff_cache.
u(force_step);
82 const Eigen::MatrixXd u_prev = diff_cache.
u(sol_step);
84 const double dt = s.solve_data.time_integrator->dt();
87 if (s.problem->is_time_dependent())
89 if (s.solve_data.friction_form)
91 if (sol_step == force_step - 1)
93 Eigen::MatrixXd surface_solution_prev = s.collision_mesh.vertices(
utils::unflatten(u_prev, s.mesh->dimension()));
94 Eigen::MatrixXd surface_solution = s.collision_mesh.vertices(
utils::unflatten(u, s.mesh->dimension()));
97 const Eigen::MatrixXd surface_velocities = (surface_solution - surface_solution_prev) / dt;
98 const double dv_dut = -1 / dt;
103 s.solve_data.friction_form->friction_potential().force_jacobian(
106 s.collision_mesh.rest_positions(),
107 surface_solution_prev,
109 barrier_contact->barrier_potential(),
110 barrier_contact->barrier_stiffness(),
111 ipc::FrictionPotential::DiffWRT::LAGGED_DISPLACEMENTS)
112 + s.solve_data.friction_form->friction_potential().force_jacobian(
115 s.collision_mesh.rest_positions(),
116 surface_solution_prev,
118 barrier_contact->barrier_potential(),
119 barrier_contact->barrier_stiffness(),
120 ipc::FrictionPotential::DiffWRT::VELOCITIES)
157 hessian_prev = s.collision_mesh.to_full_dof(hessian_prev);
170 if (s.solve_data.tangential_adhesion_form)
173 if (sol_step == force_step - 1)
177 Eigen::MatrixXd surface_solution_prev = s.collision_mesh.vertices(
utils::unflatten(u_prev, s.mesh->dimension()));
178 Eigen::MatrixXd surface_solution = s.collision_mesh.vertices(
utils::unflatten(u, s.mesh->dimension()));
181 const Eigen::MatrixXd surface_velocities = (surface_solution - surface_solution_prev) / dt;
182 const double dv_dut = -1 / dt;
184 adhesion_hessian_prev =
185 s.solve_data.tangential_adhesion_form->tangential_adhesion_potential().force_jacobian(
188 s.collision_mesh.rest_positions(),
189 surface_solution_prev,
191 s.solve_data.normal_adhesion_form->normal_adhesion_potential(),
193 ipc::TangentialPotential::DiffWRT::LAGGED_DISPLACEMENTS)
194 + s.solve_data.tangential_adhesion_form->tangential_adhesion_potential().force_jacobian(
197 s.collision_mesh.rest_positions(),
198 surface_solution_prev,
200 s.solve_data.normal_adhesion_form->normal_adhesion_potential(),
202 ipc::TangentialPotential::DiffWRT::VELOCITIES)
205 adhesion_hessian_prev *= -1;
207 adhesion_hessian_prev = s.collision_mesh.to_full_dof(adhesion_hessian_prev);
209 hessian_prev += adhesion_hessian_prev;
213 if (s.damping_assembler->is_valid() && sol_step == force_step - 1)
217 s.damping_prev_assembler->assemble_hessian(s.mesh->is_volume(), s.n_bases,
false, s.bases, s.geom_bases(), s.ass_vals_cache, force_step * s.args[
"time"][
"dt"].get<
double>() + s.args[
"time"][
"t0"].get<
double>(), dt, u, u_prev, mat_cache, damping_hessian_prev);
219 hessian_prev += damping_hessian_prev;
222 if (sol_step == force_step - 1)
225 s.solve_data.body_form->hessian_wrt_u_prev(u_prev, force_step * dt, body_force_hessian);
226 hessian_prev += body_force_hessian;
232 Eigen::MatrixXd solve_static_adjoint(
const State &state,
const DiffCache &diff_cache,
const Eigen::MatrixXd &adjoint_rhs)
236 Eigen::MatrixXd
b = adjoint_rhs;
238 Eigen::MatrixXd adjoint;
239 if (s.static_linear_solver_cache)
241 b(s.boundary_nodes, Eigen::all).setZero();
244 const int full_size = A.rows();
245 const int problem_dim = s.problem->is_scalar() ? 1 : s.mesh->dimension();
246 int precond_num = problem_dim * s.n_bases;
248 b.conservativeResizeLike(Eigen::MatrixXd::Zero(A.rows(),
b.cols()));
250 std::vector<int> boundary_nodes_tmp;
251 if (s.has_periodic_bc())
253 boundary_nodes_tmp = s.periodic_bc->full_to_periodic(s.boundary_nodes);
254 precond_num = s.periodic_bc->full_to_periodic(A);
255 b = s.periodic_bc->full_to_periodic(b,
true);
258 boundary_nodes_tmp = s.boundary_nodes;
260 adjoint.setZero(s.ndof(), adjoint_rhs.cols());
261 for (
int i = 0; i <
b.cols(); i++)
263 Eigen::VectorXd
x,
tmp;
265 dirichlet_solve_prefactorized(*s.static_linear_solver_cache, A, tmp, boundary_nodes_tmp,
x);
267 if (s.has_periodic_bc())
268 adjoint.col(i) = s.periodic_bc->periodic_to_full(full_size,
x);
275 auto solver = polysolve::linear::Solver::create(s.args[
"solver"][
"adjoint_linear"],
adjoint_logger());
283 if (!s.is_homogenization())
285 adjoint.setZero(s.ndof(), adjoint_rhs.cols());
286 for (
int i = 0; i <
b.cols(); i++)
288 Eigen::VectorXd
tmp =
b.col(i);
289 tmp(s.boundary_nodes).setZero();
292 x.setZero(
tmp.size());
293 dirichlet_solve(*solver, A, tmp, s.boundary_nodes,
x, A.rows(),
"",
false,
false,
false);
296 adjoint(s.boundary_nodes, i) = -
b(s.boundary_nodes, i);
301 solver->analyze_pattern(A, A.rows());
302 solver->factorize(A);
304 adjoint.setZero(adjoint_rhs.rows(), adjoint_rhs.cols());
305 for (
int i = 0; i <
b.cols(); i++)
307 Eigen::MatrixXd
tmp =
b.col(i);
310 x.setZero(
tmp.size());
311 solver->solve(tmp,
x);
312 x.conservativeResize(adjoint.rows());
322 Eigen::MatrixXd solve_transient_adjoint(
const State &state,
const DiffCache &diff_cache,
const Eigen::MatrixXd &adjoint_rhs)
326 const double dt = s.
args[
"time"][
"dt"];
327 const int time_steps = s.args[
"time"][
"time_steps"];
330 if (s.args[
"time"][
"integrator"].is_string())
332 else if (s.args[
"time"][
"integrator"][
"type"] ==
"ImplicitEuler")
334 else if (s.args[
"time"][
"integrator"][
"type"] ==
"BDF")
335 bdf_order = s.args[
"time"][
"integrator"][
"steps"].get<
int>();
339 assert(adjoint_rhs.cols() == time_steps + 1);
341 const int cols_per_adjoint = time_steps + 1;
342 Eigen::MatrixXd adjoints;
343 adjoints.setZero(s.ndof(), cols_per_adjoint * 2);
347 replace_rows_by_identity(reduced_mass, s.mass, s.boundary_nodes);
349 Eigen::MatrixXd sum_alpha_p, sum_alpha_nu;
350 for (
int i = time_steps; i >= 0; --i)
353 sum_alpha_p.setZero(s.ndof(), 1);
354 sum_alpha_nu.setZero(s.ndof(), 1);
356 const int num = std::min(bdf_order, time_steps - i);
358 Eigen::VectorXd bdf_coeffs = Eigen::VectorXd::Zero(num);
359 for (
int j = 0; j < bdf_order && i + j < time_steps; ++j)
362 sum_alpha_p = adjoints.middleCols(i + 1, num) * bdf_coeffs;
363 sum_alpha_nu = adjoints.middleCols(cols_per_adjoint + i + 1, num) * bdf_coeffs;
366 Eigen::VectorXd rhs_ = -reduced_mass.transpose() * sum_alpha_nu - adjoint_rhs.col(i);
367 for (
int j = 1; j <= bdf_order; j++)
369 if (i + j > time_steps)
373 compute_force_jacobian_prev(state, diff_cache, i + j, i, gradu_h_prev);
375 tmp(s.boundary_nodes).setZero();
376 rhs_ += -gradu_h_prev.transpose() *
tmp;
383 rhs_ += (1. / beta_dt) * (diff_cache.
gradu_h(i) - reduced_mass).transpose() * sum_alpha_p;
387 Eigen::VectorXd b_ = rhs_;
388 b_(s.boundary_nodes).setZero();
390 auto solver = polysolve::linear::Solver::create(s.args[
"solver"][
"adjoint_linear"],
adjoint_logger());
393 dirichlet_solve(*solver, A, b_, s.boundary_nodes,
x, A.rows(),
"",
false,
false,
false);
394 adjoints.col(i + cols_per_adjoint) =
x;
398 Eigen::VectorXd
tmp = rhs_(s.boundary_nodes);
399 if (i + 1 < cols_per_adjoint)
400 tmp += (-2. / beta_dt) * adjoints(s.boundary_nodes, i + 1);
401 if (i + 2 < cols_per_adjoint)
402 tmp += (1. / beta_dt) * adjoints(s.boundary_nodes, i + 2);
404 tmp -= (diff_cache.
gradu_h(i).transpose() * adjoints.col(i + cols_per_adjoint))(s.boundary_nodes);
405 adjoints(s.boundary_nodes, i + cols_per_adjoint) =
tmp;
406 adjoints.col(i) = beta_dt * adjoints.col(i + cols_per_adjoint) - sum_alpha_p;
410 adjoints.col(i) = -reduced_mass.transpose() * sum_alpha_p;
411 adjoints.col(i + cols_per_adjoint) = rhs_;
417 Eigen::MatrixXd solve_adjoint(
const State &state,
const DiffCache &diff_cache,
const Eigen::MatrixXd &rhs)
419 if (state.
problem->is_time_dependent())
420 return solve_transient_adjoint(state, diff_cache, rhs);
422 return solve_static_adjoint(state, diff_cache, rhs);
444 if (s.problem->is_time_dependent())
464 for (
const auto &lb : s.total_local_boundary)
466 const int e = lb.element_id();
467 for (
int i = 0; i < lb.size(); ++i)
469 const int primitive_global_id = lb.global_primitive_id(i);
470 const int boundary_id = s.mesh->get_boundary_id(primitive_global_id);
471 const auto nodes = gbases[e].local_nodes_for_primitive(primitive_global_id, *s.mesh);
473 if (boundary_id == surface_selection)
475 for (
long n = 0; n < nodes.size(); ++n)
477 const int g_id = gbases[e].bases[nodes(n)].global()[0].index;
479 if (std::count(node_ids.begin(), node_ids.end(), g_id) == 0)
480 node_ids.push_back(g_id);
494 for (
const auto &lb : s.total_local_boundary)
496 const int e = lb.element_id();
497 for (
int i = 0; i < lb.size(); ++i)
499 const int primitive_global_id = lb.global_primitive_id(i);
500 const auto nodes = gbases[e].local_nodes_for_primitive(primitive_global_id, *s.mesh);
502 for (
long n = 0; n < nodes.size(); ++n)
504 const int g_id = gbases[e].bases[nodes(n)].global()[0].index;
506 if (std::count(node_ids.begin(), node_ids.end(), g_id) == 0)
507 node_ids.push_back(g_id);
520 for (
int e = 0; e < gbases.size(); e++)
522 const int body_id = s.mesh->get_body_id(e);
523 if (body_id == volume_selection)
524 for (
const auto &gbs : gbases[e].bases)
525 for (
const auto &g : gbs.global())
526 node_ids.push_back(g.index);
Storage for additional data required by differntial code.
int bdf_order(int step) const
const ipc::TangentialCollisions & friction_collision_set(int step) const
const Eigen::MatrixXd & adjoint_mat() const
Eigen::VectorXd u(int step) const
void cache_adjoints(const Eigen::MatrixXd &adjoint_mat)
const StiffnessMatrix & gradu_h(int step) const
const ipc::TangentialCollisions & tangential_adhesion_collision_set(int step) const
main class that contains the polyfem solver and all its state
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
json args
main input arguments containing all defaults
static double betas(const int i)
Retrieve the value of beta used for BDF with i steps.
static const std::vector< double > & alphas(const int i)
Retrieve the alphas used for BDF with i steps.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
void compute_volume_node_ids(const State &state, const int volume_selection, std::vector< int > &node_ids)
void compute_surface_node_ids(const State &state, const int surface_selection, std::vector< int > &node_ids)
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
void log_and_throw_adjoint_error(const std::string &msg)
Eigen::MatrixXd get_adjoint_mat(const State &state, const DiffCache &diff_cache, int type)
Get adjoint parameter nu or p.
void compute_total_surface_node_ids(const State &state, std::vector< int > &node_ids)
void solve_adjoint_cached(const State &state, DiffCache &diff_cache, const Eigen::MatrixXd &rhs)
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix