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 });
68 }
69
70 Eigen::VectorXd SpatialIntegralForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const
71 {
72 if (&state != &state_)
73 return Eigen::VectorXd::Zero(state.ndof());
74
75 assert(time_step < state_.diff_cached.size());
76
77 Eigen::VectorXd rhs;
79
80 return rhs * weight();
81 }
82
84 {
86
87 const std::string formulation = state_.formulation();
88
89 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) {
90 val.setZero(grad_u.rows(), 1);
91 const int dim = u.cols();
92 Eigen::MatrixXd grad_u_q;
93 for (int q = 0; q < grad_u.rows(); q++)
94 {
95 grad_u_q = utils::unflatten(grad_u.row(q), u.cols());
96 if (formulation == "LinearElasticity")
97 {
98 grad_u_q = (grad_u_q + grad_u_q.transpose()).eval() / 2.;
99 val(q) = mu(q) * (grad_u_q.transpose() * grad_u_q).trace() + lambda(q) / 2. * grad_u_q.trace() * grad_u_q.trace();
100 }
101 else if (formulation == "NeoHookean")
102 {
103 Eigen::MatrixXd def_grad = grad_u_q + Eigen::MatrixXd::Identity(dim, dim);
104 double log_det_j = log(def_grad.determinant());
105 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;
106 }
107 else
108 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
109 }
110 });
111
112 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) {
113 val.setZero(grad_u.rows(), grad_u.cols());
114 Eigen::MatrixXd grad_u_q, def_grad, FmT, stress;
115 for (int q = 0; q < grad_u.rows(); q++)
116 {
117 grad_u_q = utils::unflatten(grad_u.row(q), u.cols());
118 if (formulation == "LinearElasticity")
119 {
120 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());
121 }
122 else if (formulation == "NeoHookean")
123 {
124 def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
125 FmT = def_grad.inverse().transpose();
126 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
127 }
128 else
129 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
130 val.row(q) = utils::flatten(stress);
131 }
132 });
133
134 return j;
135 }
136
137 void ElasticEnergyForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
138 {
141 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
142 return Eigen::VectorXd::Zero(0).eval();
143 });
144 }
145
147 {
149
150 const std::string formulation = state_.formulation();
151 const int power = in_power_;
152
153 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) {
154 val.setZero(grad_u.rows(), 1);
155 const double dt = state.problem->is_time_dependent() ? state.args["time"]["dt"].get<double>() : 0;
156
157 Eigen::MatrixXd grad_u_q, stress, grad_unused;
158 for (int q = 0; q < grad_u.rows(); q++)
159 {
160 if (formulation == "Laplacian")
161 stress = grad_u.row(q);
162 else
163 {
164 vector2matrix(grad_u.row(q), grad_u_q);
165 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);
166 }
167 val(q) = pow(stress.squaredNorm(), power / 2.);
168 }
169 });
170
171 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) {
172 val.setZero(grad_u.rows(), grad_u.cols());
173 const double dt = state.problem->is_time_dependent() ? state.args["time"]["dt"].get<double>() : 0;
174 const int dim = state.mesh->dimension();
175
176 if (formulation == "Laplacian")
177 {
178 for (int q = 0; q < grad_u.rows(); q++)
179 {
180 const double coef = power * pow(grad_u.row(q).squaredNorm(), power / 2. - 1.);
181 val.row(q) = coef * grad_u.row(q);
182 }
183 }
184 else
185 {
186 Eigen::MatrixXd grad_u_q, stress, stress_dstress;
187 for (int q = 0; q < grad_u.rows(); q++)
188 {
189 vector2matrix(grad_u.row(q), grad_u_q);
190 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);
191
192 const double coef = power * pow(stress.squaredNorm(), power / 2. - 1.);
193 val.row(q) = coef * utils::flatten(stress_dstress);
194 }
195 }
196 });
197
198 return j;
199 }
200
201 void StressNormForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
202 {
205 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
206 return Eigen::VectorXd::Zero(0).eval();
207 });
208 }
209
211 {
213
214 if (state_.formulation() != "LinearElasticity")
215 log_and_throw_adjoint_error("[{}] Only Linear Elasticity formulation is supported!", name());
216
217 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) {
218 val.setZero(grad_u.rows(), 1);
219
220 Eigen::MatrixXd grad_u_q, stress;
221 for (int q = 0; q < grad_u.rows(); q++)
222 {
223 vector2matrix(grad_u.row(q), grad_u_q);
224 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());
225 val(q) = (stress.array() * grad_u_q.array()).sum();
226 }
227 });
228
229 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) {
230 val.setZero(grad_u.rows(), grad_u.cols());
231 const int dim = u.cols();
232
233 Eigen::MatrixXd grad_u_q, stress;
234 for (int q = 0; q < grad_u.rows(); q++)
235 {
236 vector2matrix(grad_u.row(q), grad_u_q);
237 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());
238 val.row(q) = 2 * utils::flatten(stress);
239 }
240 });
241
242 return j;
243 }
244
245 void ComplianceForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
246 {
247 const double dt = state_.problem->is_time_dependent() ? state_.args["time"]["dt"].get<double>() : 0;
248 const double t = state_.problem->is_time_dependent() ? dt * time_step + state_.args["time"]["t0"].get<double>() : 0;
249
252 const auto &bases = state_.bases;
253 Eigen::VectorXd term = Eigen::VectorXd::Zero(bases.size() * 2);
254 const int dim = state_.mesh->dimension();
255
256 for (int e = 0; e < bases.size(); e++)
257 {
259 state_.ass_vals_cache.compute(e, state_.mesh->is_volume(), bases[e], state_.geom_bases()[e], vals);
260
262 Eigen::VectorXd da = vals.det.array() * quadrature.weights.array();
263
264 Eigen::MatrixXd u, grad_u;
265 io::Evaluator::interpolate_at_local_vals(e, dim, dim, vals, state_.diff_cached.u(time_step), u, grad_u);
266
267 Eigen::MatrixXd grad_u_q;
268 for (int q = 0; q < quadrature.weights.size(); q++)
269 {
270 double lambda, mu;
271 lambda = state_.assembler->parameters().at("lambda")(quadrature.points.row(q), vals.val.row(q), t, e);
272 mu = state_.assembler->parameters().at("mu")(quadrature.points.row(q), vals.val.row(q), t, e);
273
274 vector2matrix(grad_u.row(q), grad_u_q);
275
276 Eigen::MatrixXd f_prime_dmu, f_prime_dlambda;
277 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);
278
279 term(e + bases.size()) += dot(f_prime_dmu, grad_u_q) * da(q);
280 term(e) += dot(f_prime_dlambda, grad_u_q) * da(q);
281 }
282 }
283 return term;
284 });
285 }
286
288 {
290 const int dim = dim_;
291
292 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) {
293 val = u.col(dim) + pts.col(dim);
294 });
295
296 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) {
297 val.setZero(u.rows(), u.cols());
298 val.col(dim).setOnes();
299 });
300
301 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) {
302 val.setZero(pts.rows(), pts.cols());
303 val.col(dim).setOnes();
304 });
305
306 return j;
307 }
308
310 {
312 const int dim = this->dim_;
313
314 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) {
315 Eigen::MatrixXd acc, grad_acc;
316 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);
317
318 val = acc.col(dim);
319 });
320
321 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) {
322 val.setZero(u.rows(), u.cols());
323 log_and_throw_adjoint_error("[{}] Gradient not implemented!", name());
324 });
325
326 return j;
327 }
328
330 {
332
333 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) {
334 Eigen::MatrixXd v, grad_v;
335 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);
336
337 val.setZero(u.rows(), 1);
338 for (int q = 0; q < v.rows(); q++)
339 {
340 const double rho = state_.mass_matrix_assembler->density()(local_pts.row(q), pts.row(q), params.t, params.elem);
341 val(q) = 0.5 * rho * v.row(q).squaredNorm();
342 }
343 });
344
345 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) {
346 val.setZero(u.rows(), u.cols());
347
348 log_and_throw_adjoint_error("[{}] Gradient not implemented!", name());
349 });
350
351 return j;
352 }
353
354 void StressForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
355 {
358 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. material!", name());
359 return Eigen::VectorXd::Zero(0).eval();
360 });
361 }
362
364 {
366
367 std::string formulation = state_.formulation();
368 auto dimensions = dimensions_;
369
370 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) {
371 val.setZero(grad_u.rows(), 1);
372
373 Eigen::MatrixXd grad_u_q, stress;
374 for (int q = 0; q < grad_u.rows(); q++)
375 {
376 vector2matrix(grad_u.row(q), grad_u_q);
377 if (formulation == "LinearElasticity")
378 {
379 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());
380 }
381 else if (formulation == "NeoHookean")
382 {
383 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
384 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
385 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
386 }
387 else
388 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
389 val(q) = stress(dimensions[0], dimensions[1]);
390 }
391 });
392
393 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) {
394 val.setZero(grad_u.rows(), grad_u.cols());
395
396 const int dim = state_.mesh->dimension();
397 Eigen::MatrixXd grad_u_q, stiffness, stress;
398 for (int q = 0; q < grad_u.rows(); q++)
399 {
400 stiffness.setZero(1, dim * dim * dim * dim);
401 vector2matrix(grad_u.row(q), grad_u_q);
402
403 if (formulation == "LinearElasticity")
404 {
405 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());
406 for (int i = 0, idx = 0; i < dim; i++)
407 for (int j = 0; j < dim; j++)
408 for (int k = 0; k < dim; k++)
409 for (int l = 0; l < dim; l++)
410 {
411 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);
412 }
413 }
414 else if (formulation == "NeoHookean")
415 {
416 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_q.rows(), grad_u_q.cols()) + grad_u_q;
417 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
418 stress = mu(q) * (def_grad - FmT) + lambda(q) * std::log(def_grad.determinant()) * FmT;
419 Eigen::VectorXd FmT_vec = utils::flatten(FmT);
420 double J = def_grad.determinant();
421 double tmp1 = mu(q) - lambda(q) * std::log(J);
422 for (int i = 0, idx = 0; i < dim; i++)
423 for (int j = 0; j < dim; j++)
424 for (int k = 0; k < dim; k++)
425 for (int l = 0; l < dim; l++)
426 {
427 stiffness(idx++) = mu(q) * delta(i, k) * delta(j, l) + tmp1 * FmT(i, l) * FmT(k, j);
428 }
429 stiffness += lambda(q) * utils::flatten(FmT_vec * FmT_vec.transpose()).transpose();
430 }
431 else
432 log_and_throw_adjoint_error("[{}] Unknown formulation {}!", name(), formulation);
433
434 val.row(q) = stiffness.block(0, (dimensions[0] * dim + dimensions[1]) * dim * dim, 1, dim * dim);
435 }
436 });
437
438 return j;
439 }
440
442 {
444
445 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) {
446 val.setOnes(grad_u.rows(), 1);
447 });
448
449 return j;
450 }
451} // 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:387
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::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:669
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
int ndof() const
Definition State.hpp:673
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::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:126
IntegrableFunctional get_integral_functional() const override
IntegrableFunctional get_integral_functional() const override
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
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:15
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Parameters for the functional evaluation.