PolyFEM
Loading...
Searching...
No Matches
SpatialIntegralForms.cpp
Go to the documentation of this file.
4
5#include <polyfem/State.hpp>
7
10
12
13using namespace polyfem::utils;
14
15namespace polyfem::solver
16{
17 namespace
18 {
19 bool delta(int i, int j)
20 {
21 return (i == j) ? true : false;
22 }
23
24 double dot(const Eigen::MatrixXd &A, const Eigen::MatrixXd &B) { return (A.array() * B.array()).sum(); }
25
26 class LocalThreadScalarStorage
27 {
28 public:
29 double val;
32
33 LocalThreadScalarStorage()
34 {
35 val = 0;
36 }
37 };
38
39 class LocalThreadVecStorage
40 {
41 public:
42 Eigen::MatrixXd vec;
45
46 LocalThreadVecStorage(const int size)
47 {
48 vec.resize(size, 1);
49 vec.setZero();
50 }
51 };
52 } // namespace
53
54 double SpatialIntegralForm::value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const
55 {
56 assert(time_step < state_.diff_cached.size());
58 }
59
60 void SpatialIntegralForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
61 {
62 assert(time_step < state_.diff_cached.size());
64 Eigen::VectorXd term;
65 AdjointTools::compute_shape_derivative_functional_term(this->state_, this->state_.diff_cached.u(time_step), this->get_integral_functional(), this->ids_, this->spatial_integral_type_, term, time_step);
66 return term;
67 });
69 Eigen::VectorXd term;
70 AdjointTools::compute_shape_derivative_functional_term(this->state_, this->state_.diff_cached.u(time_step), this->get_integral_functional(), this->ids_, this->spatial_integral_type_, term, time_step);
71 term *= this->weight();
72
73 const Eigen::VectorXd adjoint_rhs = this->compute_adjoint_rhs_step(time_step, x, state_);
74 const NLHomoProblem &homo_problem = *std::dynamic_pointer_cast<NLHomoProblem>(state_.solve_data.nl_problem);
75 const Eigen::VectorXd full_shape_deriv = homo_problem.reduced_to_full_shape_derivative(state_.diff_cached.disp_grad(), adjoint_rhs);
76 term += utils::flatten(utils::unflatten(full_shape_deriv, state_.mesh->dimension())(state_.primitive_to_node(), Eigen::all));
77
78 return term;
79 });
80 }
81
82 Eigen::VectorXd SpatialIntegralForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const
83 {
84 if (&state != &state_)
85 return Eigen::VectorXd::Zero(state.ndof());
86
87 assert(time_step < state_.diff_cached.size());
88
89 Eigen::VectorXd rhs;
91
92 return rhs * weight();
93 }
94
96 {
98
99 const std::string formulation = state_.formulation();
100
101 j.set_j([formulation, this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
102 val.setZero(grad_u.rows(), 1);
103 const int dim = u.cols();
104 Eigen::MatrixXd grad_u_q;
105 for (int q = 0; q < grad_u.rows(); q++)
106 {
107 grad_u_q = utils::unflatten(grad_u.row(q), u.cols());
108 if (formulation == "LinearElasticity")
109 {
110 grad_u_q = (grad_u_q + grad_u_q.transpose()).eval() / 2.;
111 val(q) = mu(q) * (grad_u_q.transpose() * grad_u_q).trace() + lambda(q) / 2. * grad_u_q.trace() * grad_u_q.trace();
112 }
113 else if (formulation == "NeoHookean")
114 {
115 Eigen::MatrixXd def_grad = grad_u_q + Eigen::MatrixXd::Identity(dim, dim);
116 double log_det_j = log(def_grad.determinant());
117 val(q) = mu(q) / 2 * ((def_grad.transpose() * def_grad).trace() - dim - 2 * log_det_j) + lambda(q) / 2 * log_det_j * log_det_j;
118 }
119 else
120 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
121 }
122 });
123
124 j.set_dj_dgradu([formulation, this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
125 val.setZero(grad_u.rows(), grad_u.cols());
126 Eigen::MatrixXd grad_u_q, def_grad, FmT, stress;
127 for (int q = 0; q < grad_u.rows(); q++)
128 {
129 grad_u_q = utils::unflatten(grad_u.row(q), u.cols());
130 if (formulation == "LinearElasticity")
131 {
132 stress = mu(q) * (grad_u_q + grad_u_q.transpose()) + lambda(q) * grad_u_q.trace() * Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols());
133 }
134 else if (formulation == "NeoHookean")
135 {
136 def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
137 FmT = def_grad.inverse().transpose();
138 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
139 }
140 else
141 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
142 val.row(q) = utils::flatten(stress);
143 }
144 });
145
146 return j;
147 }
148
149 void ElasticEnergyForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
150 {
153 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
154 return Eigen::VectorXd::Zero(0).eval();
155 });
156 }
157
159 {
161
162 const std::string formulation = state_.formulation();
163 const int power = in_power_;
164
165 j.set_j([formulation, power, &state = std::as_const(state_)](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
166 val.setZero(grad_u.rows(), 1);
167 const double dt = state.problem->is_time_dependent() ? state.args["time"]["dt"].get<double>() : 0;
168
169 Eigen::MatrixXd grad_u_q, stress, grad_unused;
170 for (int q = 0; q < grad_u.rows(); q++)
171 {
172 if (formulation == "Laplacian")
173 stress = grad_u.row(q);
174 else
175 {
176 vector2matrix(grad_u.row(q), grad_u_q);
177 state.assembler->compute_stress_grad_multiply_mat(OptAssemblerData(params.t, dt, params.elem, local_pts.row(q), pts.row(q), grad_u_q), Eigen::MatrixXd::Zero(grad_u_q.rows(), grad_u_q.cols()), stress, grad_unused);
178 }
179 val(q) = pow(stress.squaredNorm(), power / 2.);
180 }
181 });
182
183 j.set_dj_dgradu([formulation, power, &state = std::as_const(state_)](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
184 val.setZero(grad_u.rows(), grad_u.cols());
185 const double dt = state.problem->is_time_dependent() ? state.args["time"]["dt"].get<double>() : 0;
186 const int dim = state.mesh->dimension();
187
188 if (formulation == "Laplacian")
189 {
190 for (int q = 0; q < grad_u.rows(); q++)
191 {
192 const double coef = power * pow(grad_u.row(q).squaredNorm(), power / 2. - 1.);
193 val.row(q) = coef * grad_u.row(q);
194 }
195 }
196 else
197 {
198 Eigen::MatrixXd grad_u_q, stress, stress_dstress;
199 for (int q = 0; q < grad_u.rows(); q++)
200 {
201 vector2matrix(grad_u.row(q), grad_u_q);
202 state.assembler->compute_stress_grad_multiply_stress(OptAssemblerData(params.t, dt, params.elem, local_pts.row(q), pts.row(q), grad_u_q), stress, stress_dstress);
203
204 const double coef = power * pow(stress.squaredNorm(), power / 2. - 1.);
205 val.row(q) = coef * utils::flatten(stress_dstress);
206 }
207 }
208 });
209
210 return j;
211 }
212
213 void StressNormForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
214 {
217 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
218 return Eigen::VectorXd::Zero(0).eval();
219 });
220 }
221
223 {
225
226 if (state_.formulation() != "LinearElasticity")
227 log_and_throw_adjoint_error("[{}] Only Linear Elasticity formulation is supported!", name());
228
229 j.set_j([this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
230 val.setZero(grad_u.rows(), 1);
231
232 Eigen::MatrixXd grad_u_q, stress;
233 for (int q = 0; q < grad_u.rows(); q++)
234 {
235 vector2matrix(grad_u.row(q), grad_u_q);
236 stress = mu(q) * (grad_u_q + grad_u_q.transpose()) + lambda(q) * grad_u_q.trace() * Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols());
237 val(q) = (stress.array() * grad_u_q.array()).sum();
238 }
239 });
240
241 j.set_dj_dgradu([](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
242 val.setZero(grad_u.rows(), grad_u.cols());
243 const int dim = u.cols();
244
245 Eigen::MatrixXd grad_u_q, stress;
246 for (int q = 0; q < grad_u.rows(); q++)
247 {
248 vector2matrix(grad_u.row(q), grad_u_q);
249 stress = mu(q) * (grad_u_q + grad_u_q.transpose()) + lambda(q) * grad_u_q.trace() * Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols());
250 val.row(q) = 2 * utils::flatten(stress);
251 }
252 });
253
254 return j;
255 }
256
257 void ComplianceForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
258 {
259 const double dt = state_.problem->is_time_dependent() ? state_.args["time"]["dt"].get<double>() : 0;
260 const double t = state_.problem->is_time_dependent() ? dt * time_step + state_.args["time"]["t0"].get<double>() : 0;
261
264 const auto &bases = state_.bases;
265 Eigen::VectorXd term = Eigen::VectorXd::Zero(bases.size() * 2);
266 const int dim = state_.mesh->dimension();
267
268 for (int e = 0; e < bases.size(); e++)
269 {
271 state_.ass_vals_cache.compute(e, state_.mesh->is_volume(), bases[e], state_.geom_bases()[e], vals);
272
274 Eigen::VectorXd da = vals.det.array() * quadrature.weights.array();
275
276 Eigen::MatrixXd u, grad_u;
277 io::Evaluator::interpolate_at_local_vals(e, dim, dim, vals, state_.diff_cached.u(time_step), u, grad_u);
278
279 Eigen::MatrixXd grad_u_q;
280 for (int q = 0; q < quadrature.weights.size(); q++)
281 {
282 double lambda, mu;
283 lambda = state_.assembler->parameters().at("lambda")(quadrature.points.row(q), vals.val.row(q), t, e);
284 mu = state_.assembler->parameters().at("mu")(quadrature.points.row(q), vals.val.row(q), t, e);
285
286 vector2matrix(grad_u.row(q), grad_u_q);
287
288 Eigen::MatrixXd f_prime_dmu, f_prime_dlambda;
289 state_.assembler->compute_dstress_dmu_dlambda(OptAssemblerData(t, dt, e, quadrature.points.row(q), vals.val.row(q), grad_u_q), f_prime_dmu, f_prime_dlambda);
290
291 term(e + bases.size()) += dot(f_prime_dmu, grad_u_q) * da(q);
292 term(e) += dot(f_prime_dlambda, grad_u_q) * da(q);
293 }
294 }
295 return term;
296 });
297 }
298
300 {
302 const int dim = dim_;
303
304 j.set_j([dim](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
305 val = u.col(dim) + pts.col(dim);
306 });
307
308 j.set_dj_du([dim](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
309 val.setZero(u.rows(), u.cols());
310 val.col(dim).setOnes();
311 });
312
313 j.set_dj_dx([dim](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
314 val.setZero(pts.rows(), pts.cols());
315 val.col(dim).setOnes();
316 });
317
318 return j;
319 }
320
322 {
324 const int dim = this->dim_;
325
326 j.set_j([dim, &state = std::as_const(this->state_)](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
327 Eigen::MatrixXd acc, grad_acc;
328 io::Evaluator::interpolate_at_local_vals(*(state.mesh), state.problem->is_scalar(), state.bases, state.geom_bases(), params.elem, local_pts, state.diff_cached.acc(params.step), acc, grad_acc);
329
330 val = acc.col(dim);
331 });
332
333 j.set_dj_du([dim, this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
334 val.setZero(u.rows(), u.cols());
335 log_and_throw_adjoint_error("[{}] Gradient not implemented!", name());
336 });
337
338 return j;
339 }
340
342 {
344
345 j.set_j([this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
346 Eigen::MatrixXd v, grad_v;
347 io::Evaluator::interpolate_at_local_vals(*(state_.mesh), state_.problem->is_scalar(), state_.bases, state_.geom_bases(), params.elem, local_pts, state_.diff_cached.v(params.step), v, grad_v);
348
349 val.setZero(u.rows(), 1);
350 for (int q = 0; q < v.rows(); q++)
351 {
352 const double rho = state_.mass_matrix_assembler->density()(local_pts.row(q), pts.row(q), params.t, params.elem);
353 val(q) = 0.5 * rho * v.row(q).squaredNorm();
354 }
355 });
356
357 j.set_dj_du([this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
358 val.setZero(u.rows(), u.cols());
359
360 log_and_throw_adjoint_error("[{}] Gradient not implemented!", name());
361 });
362
363 return j;
364 }
365
366 void StressForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
367 {
370 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
371 return Eigen::VectorXd::Zero(0).eval();
372 });
373 }
374
376 {
378
379 std::string formulation = state_.formulation();
380 auto dimensions = dimensions_;
381
382 j.set_j([formulation, dimensions, this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
383 val.setZero(grad_u.rows(), 1);
384
385 Eigen::MatrixXd grad_u_q, stress;
386 for (int q = 0; q < grad_u.rows(); q++)
387 {
388 vector2matrix(grad_u.row(q), grad_u_q);
389 if (formulation == "LinearElasticity")
390 {
391 stress = mu(q) * (grad_u_q + grad_u_q.transpose()) + lambda(q) * grad_u_q.trace() * Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols());
392 }
393 else if (formulation == "NeoHookean")
394 {
395 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
396 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
397 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
398 }
399 else
400 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
401 val(q) = stress(dimensions[0], dimensions[1]);
402 }
403 });
404
405 j.set_dj_dgradu([formulation, dimensions, this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
406 val.setZero(grad_u.rows(), grad_u.cols());
407
408 const int dim = state_.mesh->dimension();
409 Eigen::MatrixXd grad_u_q, stiffness, stress;
410 for (int q = 0; q < grad_u.rows(); q++)
411 {
412 stiffness.setZero(1, dim * dim * dim * dim);
413 vector2matrix(grad_u.row(q), grad_u_q);
414
415 if (formulation == "LinearElasticity")
416 {
417 stress = mu(q) * (grad_u_q + grad_u_q.transpose()) + lambda(q) * grad_u_q.trace() * Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols());
418 for (int i = 0, idx = 0; i < dim; i++)
419 for (int j = 0; j < dim; j++)
420 for (int k = 0; k < dim; k++)
421 for (int l = 0; l < dim; l++)
422 {
423 stiffness(idx++) = mu(q) * delta(i, k) * delta(j, l) + mu(q) * delta(i, l) * delta(j, k) + lambda(q) * delta(i, j) * delta(k, l);
424 }
425 }
426 else if (formulation == "NeoHookean")
427 {
428 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
429 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
430 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
431 Eigen::VectorXd FmT_vec = utils::flatten(FmT);
432 double J = def_grad.determinant();
433 double tmp1 = mu(q) - lambda(q) * std::log(J);
434 for (int i = 0, idx = 0; i < dim; i++)
435 for (int j = 0; j < dim; j++)
436 for (int k = 0; k < dim; k++)
437 for (int l = 0; l < dim; l++)
438 {
439 stiffness(idx++) = mu(q) * delta(i, k) * delta(j, l) + tmp1 * FmT(i, l) * FmT(k, j);
440 }
441 stiffness += lambda(q) * utils::flatten(FmT_vec * FmT_vec.transpose()).transpose();
442 }
443 else
444 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
445
446 val.row(q) = stiffness.block(0, (dimensions[0] * dim + dimensions[1]) * dim * dim, 1, dim * dim);
447 }
448 });
449
450 return j;
451 }
452
454 {
456
457 j.set_j([](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
458 val.setOnes(grad_u.rows(), 1);
459 });
460
461 return j;
462 }
463} // namespace polyfem::solver
Eigen::MatrixXd vec
Definition Assembler.cpp:72
double val
Definition Assembler.cpp:86
QuadratureVector da
Definition Assembler.cpp:23
ElementAssemblyValues vals
Definition Assembler.cpp:22
Quadrature quadrature
int x
void set_dj_du(const functionalType &dj_du)
void set_dj_dgradu(const functionalType &dj_dgradu)
void set_j(const functionalType &j)
void set_dj_dx(const functionalType &dj_dx)
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:328
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
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::Assembler > assembler
assemblers
Definition State.hpp:155
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
std::vector< int > primitive_to_node() const
Definition State.cpp:228
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
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
int ndof() const
Definition State.hpp:653
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
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
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)
IntegrableFunctional get_integral_functional() const override
const VariableToSimulationGroup variable_to_simulations_
IntegrableFunctional get_integral_functional() const override
std::string name() const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Eigen::MatrixXd disp_grad(int step=0) const
Eigen::VectorXd v(int step) const
Eigen::VectorXd u(int step) const
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
IntegrableFunctional get_integral_functional() const override
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
IntegrableFunctional get_integral_functional() const override
TVector reduced_to_full_shape_derivative(const Eigen::MatrixXd &disp_grad, const TVector &adjoint_full) const
IntegrableFunctional get_integral_functional() const override
std::shared_ptr< solver::NLProblem > nl_problem
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
virtual IntegrableFunctional get_integral_functional() const =0
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
IntegrableFunctional get_integral_functional() const override
std::string name() const override
std::string name() const override
IntegrableFunctional get_integral_functional() const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
virtual Eigen::VectorXd apply_parametrization_jacobian(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, const std::function< Eigen::VectorXd()> &grad) const
Maps the partial gradient wrt.
IntegrableFunctional get_integral_functional() const override
double integrate_objective(const State &state, const IntegrableFunctional &j, const Eigen::MatrixXd &solution, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, const int cur_step=0)
void compute_shape_derivative_functional_term(const State &state, const Eigen::MatrixXd &solution, const IntegrableFunctional &j, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, Eigen::VectorXd &term, const int cur_time_step)
void dJ_du_step(const State &state, const IntegrableFunctional &j, const Eigen::MatrixXd &solution, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, const int cur_step, Eigen::VectorXd &term)
void vector2matrix(const Eigen::VectorXd &vec, Eigen::MatrixXd &mat)
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::Matrix< double, Eigen::Dynamic, 1, 0, MAX_QUAD_POINTS, 1 > QuadratureVector
Definition Types.hpp:17
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Parameters for the functional evaluation.