PolyFEM
Loading...
Searching...
No Matches
BodyForm.cpp
Go to the documentation of this file.
1#include "BodyForm.hpp"
2
6
10
12
14
16
17namespace polyfem::solver
18{
19 namespace
20 {
21 class LocalThreadVecStorage
22 {
23 public:
24 Eigen::MatrixXd vec;
25 assembler::ElementAssemblyValues vals, gvals;
27
28 LocalThreadVecStorage(const int size)
29 {
30 vec.resize(size, 1);
31 vec.setZero();
32 }
33 };
34 } // namespace
35 BodyForm::BodyForm(const int ndof,
36 const int n_pressure_bases,
37 const std::vector<int> &boundary_nodes,
38 const std::vector<mesh::LocalBoundary> &local_boundary,
39 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
40 const int n_boundary_samples,
41 const Eigen::MatrixXd &rhs,
42 const assembler::RhsAssembler &rhs_assembler,
43 const assembler::Density &density,
44 const bool is_formulation_mixed,
45 const bool is_time_dependent)
46 : ndof_(ndof),
47 n_pressure_bases_(n_pressure_bases),
48 boundary_nodes_(boundary_nodes),
49 local_boundary_(local_boundary),
50 local_neumann_boundary_(local_neumann_boundary),
51 n_boundary_samples_(n_boundary_samples),
52 rhs_(rhs),
53 rhs_assembler_(rhs_assembler),
54 density_(density),
55 is_formulation_mixed_(is_formulation_mixed)
56 {
57 t_ = 0;
58 if (!is_time_dependent)
59 update_current_rhs(Eigen::VectorXd());
60 }
61
66
67 void BodyForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
68 {
69 // REMEMBER -!!!!!
70 gradv = -current_rhs_;
71 }
72
73 void BodyForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
74 {
75 hessian.resize(x.size(), x.size());
76 }
77
78 void BodyForm::update_quantities(const double t, const Eigen::VectorXd &x)
79 {
80 this->t_ = t;
81 this->x_prev_ = x;
83 }
84
85 void BodyForm::update_current_rhs(const Eigen::VectorXd &x)
86 {
91
93 {
94 current_rhs_.conservativeResize(
96 current_rhs_.bottomRows(n_pressure_bases_).setZero();
97 }
98
99 // Apply Neumann boundary conditions
101 std::vector<mesh::LocalBoundary>(), std::vector<int>(),
103 current_rhs_, x, t_);
104 }
105
106 void BodyForm::force_shape_derivative(const int n_verts, const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
107 {
108 const auto &bases = rhs_assembler_.bases();
109 const auto &gbases = rhs_assembler_.gbases();
110 const int dim = rhs_assembler_.mesh().dimension();
111 const int actual_dim = rhs_assembler_.problem().is_scalar() ? 1 : dim;
112
113 const int n_elements = int(bases.size());
114 term.setZero(n_verts * dim, 1);
115
116 auto storage = utils::create_thread_storage(LocalThreadVecStorage(term.size()));
117
118 utils::maybe_parallel_for(n_elements, [&](int start, int end, int thread_id) {
119 LocalThreadVecStorage &local_storage = utils::get_local_thread_storage(storage, thread_id);
120
121 for (int e = start; e < end; ++e)
122 {
123 assembler::ElementAssemblyValues &vals = local_storage.vals;
124 rhs_assembler_.ass_vals_cache().compute(e, rhs_assembler_.mesh().is_volume(), bases[e], gbases[e], vals);
125 assembler::ElementAssemblyValues &gvals = local_storage.gvals;
126 gvals.compute(e, rhs_assembler_.mesh().is_volume(), vals.quadrature.points, gbases[e], gbases[e]);
127
128 const quadrature::Quadrature &quadrature = vals.quadrature;
129 local_storage.da = vals.det.array() * quadrature.weights.array();
130
131 Eigen::MatrixXd p, grad_p;
132 io::Evaluator::interpolate_at_local_vals(e, dim, actual_dim, vals, adjoint, p, grad_p);
133
134 Eigen::MatrixXd rhs_function;
135 rhs_assembler_.problem().rhs(rhs_assembler_.assembler(), vals.val, t, rhs_function);
136 rhs_function *= -1;
137 for (int q = 0; q < vals.val.rows(); q++)
138 {
139 const double rho = density_(quadrature.points.row(q), vals.val.row(q), t, e);
140 rhs_function.row(q) *= rho;
141 }
142
143 for (int q = 0; q < local_storage.da.size(); ++q)
144 {
145 const double value = p.row(q).dot(rhs_function.row(q)) * local_storage.da(q);
146 for (auto &v : gvals.basis_values)
147 {
148 local_storage.vec.block(v.global[0].index * dim, 0, dim, 1) += value * v.grad_t_m.row(q).transpose();
149 }
150 }
151 }
152 });
153
154 // Zero entries in p that correspond to nodes lying between dirichlet and neumann surfaces
155 // DO NOT PARALLELIZE since they both write to the same location
156 Eigen::MatrixXd adjoint_zeroed = adjoint;
157 adjoint_zeroed(boundary_nodes_, Eigen::all).setZero();
158
159 utils::maybe_parallel_for(local_neumann_boundary_.size(), [&](int start, int end, int thread_id) {
160 LocalThreadVecStorage &local_storage = utils::get_local_thread_storage(storage, thread_id);
161
162 Eigen::MatrixXd uv, points, normals;
163 Eigen::VectorXd weights;
164 Eigen::VectorXi global_ids;
165 assembler::ElementAssemblyValues vals, gvals;
166
167 for (int lb_id = start; lb_id < end; ++lb_id)
168 {
169 const auto &lb = local_neumann_boundary_[lb_id];
170 const int e = lb.element_id();
171
172 for (int i = 0; i < lb.size(); i++)
173 {
174 const int global_primitive_id = lb.global_primitive_id(i);
175
176 utils::BoundarySampler::boundary_quadrature(lb, n_boundary_samples_, rhs_assembler_.mesh(), i, false, uv, points, normals, weights);
177 global_ids.setConstant(points.rows(), global_primitive_id);
178
179 Eigen::MatrixXd reference_normals = normals;
180
181 vals.compute(e, rhs_assembler_.mesh().is_volume(), points, bases[e], gbases[e]);
182 gvals.compute(e, rhs_assembler_.mesh().is_volume(), points, gbases[e], gbases[e]);
183
184 std::vector<std::vector<Eigen::MatrixXd>> grad_normal;
185 for (int n = 0; n < gvals.jac_it.size(); ++n)
186 {
187 Eigen::MatrixXd trafo = gvals.jac_it[n].inverse();
188
189 Eigen::MatrixXd deform_mat(dim, dim);
190 deform_mat.setZero();
191 Eigen::MatrixXd jac_mat(dim, gvals.basis_values.size());
192 int b_idx = 0;
193 for (const auto &b : gvals.basis_values)
194 {
195 jac_mat.col(b_idx++) = b.grad.row(n);
196
197 for (const auto &g : b.global)
198 for (int d = 0; d < dim; ++d)
199 deform_mat.row(d) += x(g.index * dim + d) * b.grad.row(n);
200 }
201
202 trafo += deform_mat;
203 trafo = trafo.inverse();
204
205 Eigen::VectorXd displaced_normal = normals.row(n) * trafo;
206 normals.row(n) = displaced_normal / displaced_normal.norm();
207
208 std::vector<Eigen::MatrixXd> grad;
209 {
210 Eigen::MatrixXd vec = -(jac_mat.transpose() * trafo * reference_normals.row(n).transpose());
211 // Gradient of the displaced normal computation
212 for (int k = 0; k < dim; ++k)
213 {
214 Eigen::MatrixXd grad_i(jac_mat.rows(), jac_mat.cols());
215 grad_i.setZero();
216 for (int m = 0; m < jac_mat.rows(); ++m)
217 for (int l = 0; l < jac_mat.cols(); ++l)
218 grad_i(m, l) = -(reference_normals.row(n) * trafo)(m) * (jac_mat.transpose() * trafo)(l, k);
219 grad.push_back(grad_i);
220 }
221 }
222
223 {
224 Eigen::MatrixXd normalization_chain_rule = (normals.row(n).transpose() * normals.row(n));
225 normalization_chain_rule = Eigen::MatrixXd::Identity(dim, dim) - normalization_chain_rule;
226 normalization_chain_rule /= displaced_normal.norm();
227
228 Eigen::VectorXd vec(dim);
229 b_idx = 0;
230 for (const auto &b : gvals.basis_values)
231 {
232 for (int d = 0; d < dim; ++d)
233 {
234 for (int k = 0; k < dim; ++k)
235 vec(k) = grad[k](d, b_idx);
236 vec = normalization_chain_rule * vec;
237 for (int k = 0; k < dim; ++k)
238 grad[k](d, b_idx) = vec(k);
239 }
240 b_idx++;
241 }
242 }
243
244 grad_normal.push_back(grad);
245 }
246
247 Eigen::MatrixXd neumann_val;
248 rhs_assembler_.problem().neumann_bc(rhs_assembler_.mesh(), global_ids, uv, vals.val, normals, t, neumann_val);
249
250 Eigen::MatrixXd p, grad_p;
251 io::Evaluator::interpolate_at_local_vals(e, dim, actual_dim, vals, adjoint_zeroed, p, grad_p);
252
253 const Eigen::VectorXi geom_nodes = gbases[e].local_nodes_for_primitive(global_primitive_id, rhs_assembler_.mesh());
254
255 if (geom_nodes.size() != dim)
256 log_and_throw_error("Only linear geometry is supported in shape derivative!");
257
258 Eigen::VectorXd value = (p.array() * neumann_val.array()).rowwise().sum() * weights.array();
259
260 Eigen::VectorXd pressure_bc;
261 pressure_bc = (neumann_val.array() * normals.array()).rowwise().sum();
262
263 for (long n = 0; n < geom_nodes.size(); ++n)
264 {
265 const assembler::AssemblyValues &v = gvals.basis_values[geom_nodes(n)];
266
267 Eigen::MatrixXd velocity_div_mat;
268 {
269 if (rhs_assembler_.mesh().is_volume())
270 {
271 Eigen::MatrixXd V(3, 3);
272 V.row(0) = gbases[e].bases[geom_nodes(0)].global()[0].node;
273 V.row(1) = gbases[e].bases[geom_nodes(1)].global()[0].node;
274 V.row(2) = gbases[e].bases[geom_nodes(2)].global()[0].node;
275
276 velocity_div_mat = AdjointTools::face_velocity_divergence(V);
277 }
278 else
279 {
280 Eigen::MatrixXd V(2, 2);
281 V.row(0) = gbases[e].bases[geom_nodes(0)].global()[0].node;
282 V.row(1) = gbases[e].bases[geom_nodes(1)].global()[0].node;
283
284 velocity_div_mat = AdjointTools::edge_velocity_divergence(V);
285 }
286 }
287
288 // integrate j * div(gbases) over the whole boundary
289 for (int d = 0; d < dim; d++)
290 {
291 assert(v.global.size() == 1);
292 const int g_index = v.global[0].index * dim + d;
293
294 Eigen::MatrixXd gvals(v.val.size(), dim);
295 gvals.setZero();
296 gvals.col(d) = v.val;
297
298 local_storage.vec(g_index) += value.sum() * velocity_div_mat(n, d);
299 const bool is_pressure = rhs_assembler_.problem().is_boundary_pressure(rhs_assembler_.mesh().get_boundary_id(global_primitive_id));
300 // if (is_pressure)
301 // {
302 // for (int q = 0; q < vals.jac_it.size(); ++q)
303 // {
304 // Eigen::MatrixXd grad_x_f(dim, dim);
305 // for (int ki = 0; ki < dim; ++ki)
306 // for (int kj = 0; kj < dim; ++kj)
307 // grad_x_f(ki, kj) = grad_normal[q][ki](kj, geom_nodes(n));
308 // local_storage.vec(g_index) += pressure_bc(q) * (grad_x_f * gvals.row(q).transpose()).dot(p.row(q)) * weights(q);
309 // }
310 // }
311 }
312 }
313 }
314 }
315 });
316
317 for (const LocalThreadVecStorage &local_storage : storage)
318 term += local_storage.vec;
319
320 term *= -1;
321
322 if (actual_dim == dim)
323 {
324 StiffnessMatrix hess;
325 hessian_wrt_u_prev(x, t, hess);
326 if (hess.cols() == term.size())
327 term += hess.transpose() * adjoint_zeroed;
328 // TODO: fix me for P2 basis
329 }
330 }
331
332 void BodyForm::hessian_wrt_u_prev(const Eigen::VectorXd &u_prev, const double t, StiffnessMatrix &hessian) const
333 {
335 hessian *= -1;
336 }
337} // namespace polyfem::solver
Eigen::MatrixXd vec
Definition Assembler.cpp:72
QuadratureVector da
Definition Assembler.cpp:23
ElementAssemblyValues vals
Definition Assembler.cpp:22
assembler::ElementAssemblyValues gvals
Definition BodyForm.cpp:25
Quadrature quadrature
int x
void compute(const int el_index, const bool is_volume, const basis::ElementBases &basis, const basis::ElementBases &gbasis, ElementAssemblyValues &vals) const
retrieves cached basis evaluation and geometric for the given element if it doesn't exist,...
stores per element basis values at given quadrature points and geometric mapping
virtual bool is_scalar() const =0
virtual void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const =0
void set_bc(const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, Eigen::MatrixXd &rhs, const Eigen::MatrixXd &displacement=Eigen::MatrixXd(), const double t=1) const
const AssemblyValsCache & ass_vals_cache() const
const mesh::Mesh & mesh() const
const Problem & problem() const
void compute_energy_grad(const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const Density &density, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Eigen::MatrixXd &final_rhs, const double t, Eigen::MatrixXd &rhs) const
double compute_energy(const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Density &density, const int resolution, const double t) const
void compute_energy_hess(const std::vector< int > &bounday_nodes, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Eigen::MatrixXd &displacement, const double t, const bool project_to_psd, StiffnessMatrix &hess) const
const Assembler & assembler() const
const std::vector< basis::ElementBases > & gbases() const
const std::vector< basis::ElementBases > & bases() const
static void interpolate_at_local_vals(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const int el_index, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, Eigen::MatrixXd &result_grad)
interpolate solution and gradient at element (calls interpolate_at_local_vals with sol)
virtual bool is_volume() const =0
checks if mesh is volume
int dimension() const
utily for dimension
Definition Mesh.hpp:151
const assembler::RhsAssembler & rhs_assembler_
Reference to the RHS assembler.
Definition BodyForm.hpp:86
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
Definition BodyForm.cpp:73
bool is_formulation_mixed_
True if the formulation is mixed.
Definition BodyForm.hpp:89
const std::vector< int > & boundary_nodes_
Definition BodyForm.hpp:80
double t_
Current time.
Definition BodyForm.hpp:74
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time dependent quantities.
Definition BodyForm.cpp:78
void force_shape_derivative(const int n_verts, const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
Compute the derivative of the force wrt vertex positions, then multiply the resulting matrix with adj...
Definition BodyForm.cpp:106
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
Definition BodyForm.cpp:67
const std::vector< mesh::LocalBoundary > & local_neumann_boundary_
Definition BodyForm.hpp:82
const int ndof_
Number of degrees of freedom.
Definition BodyForm.hpp:75
const assembler::Density & density_
Definition BodyForm.hpp:87
BodyForm(const int ndof, 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 assembler::RhsAssembler &rhs_assembler, const assembler::Density &density, const bool is_formulation_mixed, const bool is_time_dependent)
Construct a new Body Form object.
Definition BodyForm.cpp:35
Eigen::MatrixXd current_rhs_
Cached RHS for the current time.
Definition BodyForm.hpp:91
const std::vector< mesh::LocalBoundary > & local_boundary_
Definition BodyForm.hpp:81
void update_current_rhs(const Eigen::VectorXd &x)
Update current_rhs.
Definition BodyForm.cpp:85
const Eigen::MatrixXd & rhs_
static RHS for the current time
Definition BodyForm.hpp:85
Eigen::MatrixXd x_prev_
Cached previous solution.
Definition BodyForm.hpp:78
void hessian_wrt_u_prev(const Eigen::VectorXd &u_prev, const double t, StiffnessMatrix &hessian) const
Definition BodyForm.cpp:332
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the body force form.
Definition BodyForm.cpp:62
virtual double value(const Eigen::VectorXd &x) const
Compute the value of the form multiplied with the weigth.
Definition Form.hpp:26
auto & get_local_thread_storage(Storages &storage, int thread_id)
auto create_thread_storage(const LocalStorage &initial_local_storage)
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, MAX_QUAD_POINTS, 1 > QuadratureVector
Definition Types.hpp:17
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22