PolyFEM
Loading...
Searching...
No Matches
PeriodicContactForceDerivative.cpp
Go to the documentation of this file.
2
3#include <Eigen/Core>
5#include <polyfem/State.hpp>
9#include <ipc/collisions/normal/normal_collisions.hpp>
10
11namespace polyfem::solver
12{
14 const PeriodicContactForm &form,
15 const State &state,
16 const PeriodicMeshToMesh &periodic_mesh_map,
17 const Eigen::VectorXd &periodic_mesh_representation,
18 const ipc::NormalCollisions &contact_set,
19 const Eigen::VectorXd &solution,
20 const Eigen::VectorXd &adjoint_sol,
21 Eigen::VectorXd &term)
22 {
23 const int dim = form.collision_mesh_.dim();
24 const Eigen::MatrixXd displaced_surface = form.compute_displaced_surface(form.single_to_tiled(solution));
25
26 Eigen::VectorXd tiled_term;
27
28 {
29 StiffnessMatrix dq_h = form.collision_mesh_.to_full_dof(form.barrier_potential().shape_derivative(contact_set, form.collision_mesh_, displaced_surface));
30 tiled_term = dq_h.transpose() * form.single_to_tiled(adjoint_sol);
31 }
32
33 {
34 Eigen::VectorXd force;
35 force = form.barrier_potential().gradient(contact_set, form.collision_mesh_, displaced_surface);
36 force = form.collision_mesh_.to_full_dof(force);
37 Eigen::MatrixXd adjoint_affine = utils::unflatten(adjoint_sol.tail(dim * dim), dim);
38 for (int k = 0; k < form.collision_mesh_.num_vertices(); k++)
39 {
40 const int k_full = form.collision_mesh_.to_full_vertex_id(k);
41 tiled_term.segment(k_full * dim, dim) += adjoint_affine.transpose() * force.segment(k_full * dim, dim);
42 }
43 }
44
45 {
46 StiffnessMatrix hessian_full;
47 form.BarrierContactForm::second_derivative_unweighted(form.single_to_tiled(solution), hessian_full);
48 Eigen::VectorXd tmp = form.single_to_tiled(adjoint_sol).transpose() * hessian_full;
49 Eigen::MatrixXd sol_affine = utils::unflatten(solution.tail(dim * dim), dim);
50 for (int k = 0; k < form.collision_mesh_.num_vertices(); k++)
51 {
52 const int k_full = form.collision_mesh_.to_full_vertex_id(k);
53 tiled_term.segment(k_full * dim, dim) += sol_affine.transpose() * tmp.segment(k_full * dim, dim);
54 }
55 }
56
57 // chain rule from tiled to periodic
58 Eigen::VectorXd unit_term;
59 Eigen::MatrixXd affine_term;
60 affine_term.setZero(dim, dim);
61 unit_term.setZero(tiled_term.size());
62 Eigen::MatrixXd affine = utils::unflatten(periodic_mesh_representation.tail(dim * dim), dim).transpose();
63 for (int k = 0; k < form.collision_mesh_.num_vertices(); k++)
64 {
65 const int k_full = form.collision_mesh_.to_full_vertex_id(k);
66 affine_term += tiled_term.segment(k_full * dim, dim) * form.collision_mesh_.rest_positions().row(k);
67 unit_term.segment(k_full * dim, dim) += affine.transpose() * tiled_term.segment(k_full * dim, dim);
68 }
69
70 Eigen::VectorXd single_term = state.basis_nodes_to_gbasis_nodes * form.proj.topRows(dim * form.n_single_dof_) * unit_term;
71 single_term = utils::flatten(utils::unflatten(single_term, dim)(state.primitive_to_node(), Eigen::all));
72
73 term.setZero(periodic_mesh_representation.size());
74 for (int i = 0; i < state.n_geom_bases; i++)
75 term.segment(periodic_mesh_map.full_to_periodic(i) * dim, dim).array() += single_term.segment(i * dim, dim).array();
76 term.tail(dim * dim) = Eigen::Map<Eigen::VectorXd>(affine_term.data(), dim * dim, 1);
77
78 term *= form.weight();
79 }
80} // namespace polyfem::solver
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::vector< int > primitive_to_node() const
Definition State.cpp:227
int n_geom_bases
number of geometric bases
Definition State.hpp:182
StiffnessMatrix basis_nodes_to_gbasis_nodes
Definition State.hpp:708
const ipc::BarrierPotential & barrier_potential() const
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
double weight() const override
Get the form's multiplicative constant weight.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
static void force_shape_derivative(const PeriodicContactForm &form, const State &state, const PeriodicMeshToMesh &periodic_mesh_map, const Eigen::VectorXd &periodic_mesh_representation, const ipc::NormalCollisions &contact_set, const Eigen::VectorXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term)
Form representing the contact potential and forces on a periodic mesh This form has a different input...
Eigen::VectorXd single_to_tiled(const Eigen::VectorXd &x) const
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24