PolyFEM
Loading...
Searching...
No Matches
DiffCache.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
4#include <polyfem/State.hpp>
5
7
8#include <ipc/ipc.hpp>
9#include <ipc/collisions/normal/normal_collisions.hpp>
10#include <ipc/collisions/tangential/tangential_collisions.hpp>
11
12#include <Eigen/Core>
13
14#include <cassert>
15#include <vector>
16
17namespace polyfem
18{
21 {
22 public:
25
26 void cache_adjoints(const Eigen::MatrixXd &adjoint_mat);
27
36 void cache_transient(
37 int step,
38 State &state,
39 const Eigen::MatrixXd &sol,
40 const Eigen::MatrixXd *disp_grad,
41 const Eigen::MatrixXd *pressure);
42
43 const Eigen::MatrixXd &adjoint_mat() const { return adjoint_mat_; }
44
46
47 inline int size() const { return cur_size_; }
48 inline int bdf_order(int step) const
49 {
50 assert(step < size());
51 if (step < 0)
52 step += bdf_order_.size();
53 return bdf_order_(step);
54 }
55
56 Eigen::MatrixXd disp_grad(int step = 0) const
57 {
58 assert(step < size());
59 if (step < 0)
60 step += disp_grad_.size();
61 return disp_grad_[step];
62 }
63
64 Eigen::VectorXd u(int step) const
65 {
66 assert(step < size());
67 if (step < 0)
68 step += u_.cols();
69 return u_.col(step);
70 }
71 Eigen::VectorXd v(int step) const
72 {
73 assert(step < size());
74 if (step < 0)
75 step += v_.cols();
76 return v_.col(step);
77 }
78 Eigen::VectorXd acc(int step) const
79 {
80 assert(step < size());
81 if (step < 0)
82 step += acc_.cols();
83 return acc_.col(step);
84 }
85
86 const StiffnessMatrix &gradu_h(int step) const
87 {
88 assert(step < size());
89 if (step < 0)
90 step += gradu_h_.size();
91 return gradu_h_[step];
92 }
93 // const StiffnessMatrix &gradu_h_prev(const int step) const { assert(step < size()); return gradu_h_prev_[step]; }
94
95 const ipc::NormalCollisions &collision_set(int step) const
96 {
97 assert(step < size());
98 if (step < 0)
99 step += collision_set_.size();
100 return collision_set_[step];
101 }
102 const ipc::SmoothCollisions &smooth_collision_set(int step) const
103 {
104 assert(step < size());
105 if (step < 0)
106 step += smooth_collision_set_.size();
107 return smooth_collision_set_[step];
108 }
109 const ipc::TangentialCollisions &friction_collision_set(int step) const
110 {
111 assert(step < size());
112 if (step < 0)
113 step += friction_collision_set_.size();
114 return friction_collision_set_[step];
115 }
116 const ipc::NormalCollisions &normal_adhesion_collision_set(int step) const
117 {
118 assert(step < size());
119 if (step < 0)
120 step += normal_adhesion_collision_set_.size();
122 }
123 const ipc::TangentialCollisions &tangential_adhesion_collision_set(int step) const
124 {
125 assert(step < size());
126 if (step < 0)
129 }
130
131 private:
133 int cur_size_ = 0;
134
135 // Mapping from positions of FE basis nodes to positions of geometry nodes.
137
138 std::vector<Eigen::MatrixXd> disp_grad_; // macro linear displacement in homogenization
139 Eigen::MatrixXd u_; // PDE solution
140 Eigen::MatrixXd v_; // velocity in transient elastic simulations
141 Eigen::MatrixXd acc_; // acceleration in transient elastic simulations
142
143 Eigen::VectorXi bdf_order_; // BDF orders used at each time step in forward simulation
144
145 std::vector<StiffnessMatrix> gradu_h_; // gradient of force at time T wrt. u at time T
146 // std::vector<StiffnessMatrix> gradu_h_prev_; // gradient of force at time T wrt. u at time (T-1) in transient simulations
147
148 std::vector<ipc::NormalCollisions> collision_set_;
149 std::vector<ipc::SmoothCollisions> smooth_collision_set_;
150 std::vector<ipc::TangentialCollisions> friction_collision_set_;
151
152 std::vector<ipc::NormalCollisions> normal_adhesion_collision_set_;
153 std::vector<ipc::TangentialCollisions> tangential_adhesion_collision_set_;
154
155 Eigen::MatrixXd adjoint_mat_;
156
157 void init(const int dimension, const int ndof, const int n_time_steps = 0);
158
160 const Eigen::MatrixXd &u,
162 const ipc::NormalCollisions &collision_set,
163 const ipc::SmoothCollisions &smooth_collision_set,
164 const ipc::TangentialCollisions &friction_constraint_set,
165 const ipc::NormalCollisions &normal_adhesion_set,
166 const ipc::TangentialCollisions &tangential_adhesion_set,
167 const Eigen::MatrixXd &disp_grad);
168
170 const int cur_step,
171 const int cur_bdf_order,
172 const Eigen::MatrixXd &u,
173 const Eigen::MatrixXd &v,
174 const Eigen::MatrixXd &acc,
176 // const StiffnessMatrix &gradu_h_prev,
177 const ipc::NormalCollisions &collision_set,
178 const ipc::SmoothCollisions &smooth_collision_set,
179 const ipc::TangentialCollisions &friction_collision_set);
180
182 const int cur_step,
183 const Eigen::MatrixXd &u,
185 const ipc::NormalCollisions &collision_set,
186 const ipc::SmoothCollisions &smooth_collision_set,
187 const ipc::NormalCollisions &normal_adhesion_set,
188 const Eigen::MatrixXd &disp_grad);
189 };
190} // namespace polyfem
Storage for additional data required by differntial code.
Definition DiffCache.hpp:21
void cache_transient(int step, State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd *disp_grad, const Eigen::MatrixXd *pressure)
Cache time-dependent adjoint optimization data.
int bdf_order(int step) const
Definition DiffCache.hpp:48
int size() const
Definition DiffCache.hpp:47
Eigen::MatrixXd v_
std::vector< StiffnessMatrix > gradu_h_
std::vector< ipc::NormalCollisions > normal_adhesion_collision_set_
Eigen::MatrixXd acc_
std::vector< ipc::TangentialCollisions > friction_collision_set_
const ipc::NormalCollisions & collision_set(int step) const
Definition DiffCache.hpp:95
Eigen::MatrixXd disp_grad(int step=0) const
Definition DiffCache.hpp:56
Eigen::MatrixXd adjoint_mat_
Eigen::VectorXd v(int step) const
Definition DiffCache.hpp:71
StiffnessMatrix basis_nodes_to_gbasis_nodes_
std::vector< Eigen::MatrixXd > disp_grad_
const ipc::TangentialCollisions & friction_collision_set(int step) const
const Eigen::MatrixXd & adjoint_mat() const
Definition DiffCache.hpp:43
void cache_quantities_static(const Eigen::MatrixXd &u, const StiffnessMatrix &gradu_h, const ipc::NormalCollisions &collision_set, const ipc::SmoothCollisions &smooth_collision_set, const ipc::TangentialCollisions &friction_constraint_set, const ipc::NormalCollisions &normal_adhesion_set, const ipc::TangentialCollisions &tangential_adhesion_set, const Eigen::MatrixXd &disp_grad)
Eigen::VectorXd acc(int step) const
Definition DiffCache.hpp:78
std::vector< ipc::SmoothCollisions > smooth_collision_set_
std::vector< ipc::TangentialCollisions > tangential_adhesion_collision_set_
InitialConditionOverride initial_condition_override
Initial-condition override storage for initial condition optimization.
Definition DiffCache.hpp:24
const StiffnessMatrix & basis_nodes_to_gbasis_nodes() const
Eigen::VectorXi bdf_order_
void cache_quantities_quasistatic(const int cur_step, const Eigen::MatrixXd &u, const StiffnessMatrix &gradu_h, const ipc::NormalCollisions &collision_set, const ipc::SmoothCollisions &smooth_collision_set, const ipc::NormalCollisions &normal_adhesion_set, const Eigen::MatrixXd &disp_grad)
std::vector< ipc::NormalCollisions > collision_set_
const ipc::SmoothCollisions & smooth_collision_set(int step) const
Eigen::VectorXd u(int step) const
Definition DiffCache.hpp:64
void cache_adjoints(const Eigen::MatrixXd &adjoint_mat)
const StiffnessMatrix & gradu_h(int step) const
Definition DiffCache.hpp:86
void init(const int dimension, const int ndof, const int n_time_steps=0)
const ipc::NormalCollisions & normal_adhesion_collision_set(int step) const
Eigen::MatrixXd u_
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::NormalCollisions &collision_set, const ipc::SmoothCollisions &smooth_collision_set, const ipc::TangentialCollisions &friction_collision_set)
const ipc::TangentialCollisions & tangential_adhesion_collision_set(int step) const
Runtime override for initial-condition histories.
Definition State.hpp:90
main class that contains the polyfem solver and all its state
Definition State.hpp:113
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24