/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/model/ad.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/model/ad.hh Source File
DiFfRG
ad.hh
Go to the documentation of this file.
1#pragma once
2
3// DiFfRG
6
7// external libraries
8#include <autodiff/forward/dual.hpp>
9#include <autodiff/forward/real.hpp>
10#include <deal.II/base/point.h>
11#include <deal.II/base/tensor.h>
12#include <deal.II/lac/full_matrix.h>
13#include <deal.II/lac/vector.h>
14#include <tbb/tbb.h>
15
16namespace DiFfRG
17{
18 namespace def
19 {
20 using namespace dealii;
21 using std::get;
22
23 namespace internal
24 {
25 template <typename AD_type> struct AD_tools;
26
27 template <> struct AD_tools<autodiff::dual> {
28 template <uint n, typename Vector> static std::array<autodiff::dual, n> vector_to_AD(const Vector &v)
29 {
30 std::array<autodiff::dual, n> x;
31 for (uint i = 0; i < n; ++i) {
32 x[i] = v[i];
33 }
34 return x;
35 }
36
37 template <uint n, int dim, int r, typename NT>
38 static std::array<dealii::Tensor<r, dim, autodiff::dual>, n>
39 ten_to_AD(const std::vector<dealii::Tensor<r, dim, NT>> &v)
40 {
41 static_assert(r >= 1 && r <= 2, "Only rank 1 and 2 tensors are supported.");
42 std::array<dealii::Tensor<r, dim, autodiff::dual>, n> x;
43 for (uint i = 0; i < n; ++i) {
44 if constexpr (r == 1)
45 for (uint d = 0; d < dim; ++d)
46 x[i][d] = v[i][d];
47 else if constexpr (r == 2)
48 for (uint d1 = 0; d1 < dim; ++d1)
49 for (uint d2 = 0; d2 < dim; ++d2) {
50 x[i][d1][d2] = v[i][d1][d2];
51 }
52 }
53 return x;
54 }
55 };
56
57 template <> struct AD_tools<autodiff::real> {
58 template <uint n, typename Vector> static std::array<autodiff::real, n> vector_to_AD(const Vector &v)
59 {
60 std::array<autodiff::real, n> x;
61 for (uint i = 0; i < n; ++i)
62 x[i] = v[i];
63 return x;
64 }
65
66 template <uint n, int dim, int r, typename NT>
67 static std::array<dealii::Tensor<r, dim, autodiff::real>, n>
68 ten_to_AD(const std::vector<dealii::Tensor<r, dim, NT>> &v)
69 {
70 static_assert(r >= 1 && r <= 2, "Only rank 1 and 2 tensors are supported.");
71 std::array<dealii::Tensor<r, dim, autodiff::real>, n> x;
72 for (uint i = 0; i < n; ++i) {
73 if constexpr (r == 1) {
74 for (uint d = 0; d < dim; ++d)
75 x[i][d] = v[i][d];
76 } else if constexpr (r == 2) {
77 for (uint d1 = 0; d1 < dim; ++d1)
78 for (uint d2 = 0; d2 < dim; ++d2) {
79 x[i][d1][d2] = v[i][d1][d2];
80 }
81 }
82 }
83 return x;
84 }
85 };
86 } // namespace internal
87
88 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_flux
89 {
90 Model &asImp() { return static_cast<Model &>(*this); }
91 const Model &asImp() const { return static_cast<const Model &>(*this); }
93
94 public:
95 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
96 void jacobian_flux_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &jF, const Point<dim> &p,
97 const Vector &sol) const
98 {
99 using Components = typename Model::Components;
100 static_assert(n_to == Components::count_fe_functions(),
101 "jacobian_flux_grad: n_to must equal count_fe_functions()");
102 static_assert(n_from == Components::count_fe_functions(),
103 "jacobian_flux_grad: n_from must equal count_fe_functions()");
104
105 const auto &u = get<tup_idx>(sol);
106 auto du = AD_tools::template ten_to_AD<n_from>(u);
107 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
109 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
110 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
111 for (uint d = 0; d < dim; ++d) {
112 res = {};
113 seed(du[j][d]);
114 asImp().flux(res, p, Vector::as(ad_sol));
115 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
116 for (uint dd = 0; dd < dim; ++dd) {
117 jF(i, j)[dd][d] = grad(res[i][dd]);
118 }
119 }
120 unseed(du[j][d]);
121 }
122 }
123 }
124
125 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
126 void jacobian_flux_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &jF, const Point<dim> &p,
127 const Vector &sol) const
128 {
129 using Components = typename Model::Components;
130 static_assert(n_to == Components::count_fe_functions(),
131 "jacobian_flux_hess: n_to must equal count_fe_functions()");
132 static_assert(n_from == Components::count_fe_functions(),
133 "jacobian_flux_hess: n_from must equal count_fe_functions()");
134
135 const auto &u = get<tup_idx>(sol);
136 auto du = AD_tools::template ten_to_AD<n_from>(u);
137 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
139 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
140 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
141 for (uint d1 = 0; d1 < dim; ++d1)
142 for (uint d2 = 0; d2 < dim; ++d2) {
143 res = {};
144 seed(du[j][d1][d2]);
145 asImp().flux(res, p, Vector::as(ad_sol));
146 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
147 for (uint d = 0; d < dim; ++d) {
148 jF(i, j)[d][d1][d2] = grad(res[i][d]);
149 }
150 }
151 unseed(du[j][d1][d2]);
152 }
153 }
154 }
155
156 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
157 void jacobian_flux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jF, const Point<dim> &p,
158 const Vector &sol) const
159 {
160 using Components = typename Model::Components;
161 static_assert(n_to == Components::count_fe_functions(),
162 "jacobian_flux_extr: n_to must equal count_fe_functions()");
163 static_assert(n_from == Components::count_extractors(),
164 "jacobian_flux_extr: n_from must equal count_extractors()");
165
166 const auto &e = get<tup_idx>(sol);
167 auto de = AD_tools::template vector_to_AD<n_from>(e);
168 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(de),
170 for (uint j = 0; j < Components::count_extractors(); ++j) {
171 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
172 seed(de[j]);
173 asImp().flux(res, p, Vector::as(ad_sol));
174 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
175 for (uint d = 0; d < dim; ++d) {
176 jF(i, j)[d] = grad(res[i][d]);
177 }
178 }
179 unseed(de[j]);
180 }
181 }
182
183 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
184 void jacobian_flux(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jF, const Point<dim> &p,
185 const Vector &sol) const
186 {
187 using Components = typename Model::Components;
188 static_assert(n_to == Components::count_fe_functions(to),
189 "jacobian_flux: n_to must equal count_fe_functions(to)");
190 static_assert(n_from == Components::count_fe_functions(from),
191 "jacobian_flux: n_from must equal count_fe_functions(from)");
192
193 if constexpr (to == 0) {
194 const auto &u = get<from>(sol);
195 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
196 auto ad_sol = std::tuple_cat(tuple_first<from>(sol), std::tie(du),
198 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
199 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
200 seed(du[j]);
201 asImp().flux(res, p, Vector::as(ad_sol));
202 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
203 for (uint d = 0; d < dim; ++d) {
204 jF(i, j)[d] = grad(res[i][d]);
205 }
206 }
207 unseed(du[j]);
208 }
209 } else {
210 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol);
211 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
212 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
213 // take derivative with respect to jth variable
214 seed(du[j]);
215 asImp().template ldg_flux<to>(res, p, du);
216 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
217 for (uint d = 0; d < dim; ++d) {
218 jF(i, j)[d] = grad(res[i][d]);
219 }
220 }
221 unseed(du[j]);
222 }
223 }
224 }
225 };
226
227 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_source
228 {
229 Model &asImp() { return static_cast<Model &>(*this); }
230 const Model &asImp() const { return static_cast<const Model &>(*this); }
232
233 public:
234 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
235 void jacobian_source_grad(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jS, const Point<dim> &p,
236 const Vector &sol) const
237 {
238 using Components = typename Model::Components;
239 static_assert(n_to == Components::count_fe_functions(),
240 "jacobian_source_grad: n_to must equal count_fe_functions()");
241 static_assert(n_from == Components::count_fe_functions(),
242 "jacobian_source_grad: n_from must equal count_fe_functions()");
243
244 const auto &u = get<tup_idx>(sol);
245 auto du = AD_tools::template ten_to_AD<n_from>(u);
246 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
248 std::array<AD_type, Components::count_fe_functions()> res{{}};
249 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
250 for (uint d = 0; d < dim; ++d) {
251 res = {};
252 seed(du[j][d]);
253 asImp().source(res, p, Vector::as(ad_sol));
254 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
255 jS(i, j)[d] = grad(res[i]);
256 }
257 unseed(du[j][d]);
258 }
259 }
260 }
261
262 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
263 void jacobian_source_hess(SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &jS, const Point<dim> &p,
264 const Vector &sol) const
265 {
266 using Components = typename Model::Components;
267 static_assert(n_to == Components::count_fe_functions(),
268 "jacobian_source_hess: n_to must equal count_fe_functions()");
269 static_assert(n_from == Components::count_fe_functions(),
270 "jacobian_source_hess: n_from must equal count_fe_functions()");
271
272 const auto &u = get<tup_idx>(sol);
273 auto du = AD_tools::template ten_to_AD<n_from>(u);
274 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
276 std::array<AD_type, Components::count_fe_functions()> res{{}};
277 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
278 for (uint d1 = 0; d1 < dim; ++d1) {
279 for (uint d2 = 0; d2 < dim; ++d2) {
280 res = {};
281 seed(du[j][d1][d2]);
282 asImp().source(res, p, Vector::as(ad_sol));
283 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
284 jS(i, j)[d1][d2] = grad(res[i]);
285 }
286 unseed(du[j][d1][d2]);
287 }
288 }
289 }
290 }
291
292 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
293 void jacobian_source_extr(SimpleMatrix<NT, n_to, n_from> &jS, const Point<dim> &p, const Vector &sol) const
294 {
295 using Components = typename Model::Components;
296 static_assert(n_to == Components::count_fe_functions(),
297 "jacobian_source_extr: n_to must equal count_fe_functions()");
298 static_assert(n_from == Components::count_extractors(),
299 "jacobian_source_extr: n_from must equal count_extractors()");
300
301 const auto &e = get<tup_idx>(sol);
302 auto de = AD_tools::template vector_to_AD<n_from>(e);
303 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(de),
305 for (uint j = 0; j < n_from; ++j) {
306 std::array<AD_type, Components::count_fe_functions()> res{{}};
307 seed(de[j]);
308 asImp().source(res, p, Vector::as(ad_sol));
309 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
310 jS(i, j) = grad(res[i]);
311 }
312 unseed(de[j]);
313 }
314 }
315
316 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
317 void jacobian_source(SimpleMatrix<NT, n_to, n_from> &jS, const Point<dim> &p, const Vector &sol) const
318 {
319 using Components = typename Model::Components;
320 static_assert(n_to == Components::count_fe_functions(to),
321 "jacobian_source: n_to must equal count_fe_functions(to)");
322 static_assert(n_from == Components::count_fe_functions(from),
323 "jacobian_source: n_from must equal count_fe_functions(from)");
324
325 if constexpr (to == 0) {
326 const auto &u = get<from>(sol);
327 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
328 auto ad_sol = std::tuple_cat(tuple_first<from>(sol), std::tie(du),
330 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
331 std::array<AD_type, Components::count_fe_functions(to)> res{{}};
332 seed(du[j]);
333 asImp().source(res, p, Vector::as(ad_sol));
334 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
335 jS(i, j) = grad(res[i]);
336 }
337 unseed(du[j]);
338 }
339 } else {
340 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol);
341 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
342 std::array<AD_type, Components::count_fe_functions(to)> res{{}};
343 // take derivative with respect to jth variable
344 seed(du[j]);
345 asImp().template ldg_source<to>(res, p, du);
346 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
347 jS(i, j) = grad(res[i]);
348 }
349 unseed(du[j]);
350 }
351 }
352 }
353 };
354
355 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_flux_source
356 {
357 Model &asImp() { return static_cast<Model &>(*this); }
358 const Model &asImp() const { return static_cast<const Model &>(*this); }
360
361 public:
362 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
363 void jacobian_flux_source(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jF, SimpleMatrix<NT, n_to, n_from> &jS,
364 const Point<dim> &p, const Vector &sol) const
365 {
366 using Components = typename Model::Components;
367 static_assert(n_to == Components::count_fe_functions(to),
368 "jacobian_flux_source: n_to must equal count_fe_functions(to)");
369 static_assert(n_from == Components::count_fe_functions(from),
370 "jacobian_flux_source: n_from must equal count_fe_functions(from)");
371
372 if constexpr (to == 0) {
373 const auto &u = get<from>(sol);
374 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
375 auto ad_sol = std::tuple_cat(tuple_first<from>(sol), std::tie(du),
377 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
378 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res_flux{{}};
379 std::array<AD_type, Components::count_fe_functions(to)> res_source{{}};
380 seed(du[j]);
381 asImp().flux(res_flux, p, Vector::as(ad_sol));
382 asImp().source(res_source, p, Vector::as(ad_sol));
383 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
384 for (uint d = 0; d < dim; ++d)
385 jF(i, j)[d] = grad(res_flux[i][d]);
386 jS(i, j) = grad(res_source[i]);
387 }
388 unseed(du[j]);
389 }
390 } else {
391 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol);
392 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
393 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res_flux{{}};
394 std::array<AD_type, Components::count_fe_functions(to)> res_source{{}};
395 seed(du[j]);
396 asImp().template ldg_flux<to>(res_flux, p, du);
397 asImp().template ldg_source<to>(res_source, p, du);
398 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
399 for (uint d = 0; d < dim; ++d)
400 jF(i, j)[d] = grad(res_flux[i][d]);
401 jS(i, j) = grad(res_source[i]);
402 }
403 unseed(du[j]);
404 }
405 }
406 }
407
408 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
409 void jacobian_flux_source_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &jF,
410 SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jS, const Point<dim> &p,
411 const Vector &sol) const
412 {
413 using Components = typename Model::Components;
414 static_assert(n_to == Components::count_fe_functions(),
415 "jacobian_flux_source_grad: n_to must equal count_fe_functions()");
416 static_assert(n_from == Components::count_fe_functions(),
417 "jacobian_flux_source_grad: n_from must equal count_fe_functions()");
418
419 const auto &u = get<tup_idx>(sol);
420 auto du = AD_tools::template ten_to_AD<n_from>(u);
421 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
423 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res_flux{{}};
424 std::array<AD_type, Components::count_fe_functions()> res_source{{}};
425 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
426 for (uint d = 0; d < dim; ++d) {
427 res_flux = {};
428 res_source = {};
429 seed(du[j][d]);
430 asImp().flux(res_flux, p, Vector::as(ad_sol));
431 asImp().source(res_source, p, Vector::as(ad_sol));
432 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
433 for (uint dd = 0; dd < dim; ++dd)
434 jF(i, j)[dd][d] = grad(res_flux[i][dd]);
435 jS(i, j)[d] = grad(res_source[i]);
436 }
437 unseed(du[j][d]);
438 }
439 }
440 }
441
442 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
443 void jacobian_flux_source_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &jF,
444 SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &jS, const Point<dim> &p,
445 const Vector &sol) const
446 {
447 using Components = typename Model::Components;
448 static_assert(n_to == Components::count_fe_functions(),
449 "jacobian_flux_source_hess: n_to must equal count_fe_functions()");
450 static_assert(n_from == Components::count_fe_functions(),
451 "jacobian_flux_source_hess: n_from must equal count_fe_functions()");
452
453 const auto &u = get<tup_idx>(sol);
454 auto du = AD_tools::template ten_to_AD<n_from>(u);
455 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
457 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res_flux{{}};
458 std::array<AD_type, Components::count_fe_functions()> res_source{{}};
459 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
460 for (uint d1 = 0; d1 < dim; ++d1) {
461 for (uint d2 = 0; d2 < dim; ++d2) {
462 res_flux = {};
463 res_source = {};
464 seed(du[j][d1][d2]);
465 asImp().flux(res_flux, p, Vector::as(ad_sol));
466 asImp().source(res_source, p, Vector::as(ad_sol));
467 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
468 for (uint d = 0; d < dim; ++d)
469 jF(i, j)[d][d1][d2] = grad(res_flux[i][d]);
470 jS(i, j)[d1][d2] = grad(res_source[i]);
471 }
472 unseed(du[j][d1][d2]);
473 }
474 }
475 }
476 }
477
478 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
479 void jacobian_flux_source_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jF,
480 SimpleMatrix<NT, n_to, n_from> &jS, const Point<dim> &p,
481 const Vector &sol) const
482 {
483 using Components = typename Model::Components;
484 static_assert(n_to == Components::count_fe_functions(),
485 "jacobian_flux_source_extr: n_to must equal count_fe_functions()");
486 static_assert(n_from == Components::count_extractors(),
487 "jacobian_flux_source_extr: n_from must equal count_extractors()");
488
489 const auto &e = get<tup_idx>(sol);
490 auto de = AD_tools::template vector_to_AD<n_from>(e);
491 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(de),
493 for (uint j = 0; j < Components::count_extractors(); ++j) {
494 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res_flux{{}};
495 std::array<AD_type, Components::count_fe_functions()> res_source{{}};
496 seed(de[j]);
497 asImp().flux(res_flux, p, Vector::as(ad_sol));
498 asImp().source(res_source, p, Vector::as(ad_sol));
499 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
500 for (uint d = 0; d < dim; ++d)
501 jF(i, j)[d] = grad(res_flux[i][d]);
502 jS(i, j) = grad(res_source[i]);
503 }
504 unseed(de[j]);
505 }
506 }
507 };
508
509 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_numflux
510 {
511 Model &asImp() { return static_cast<Model &>(*this); }
512 const Model &asImp() const { return static_cast<const Model &>(*this); }
514
515 public:
516 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
517 void jacobian_numflux_grad(std::array<SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from>, 2> &jNF,
518 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector_s &sol_s,
519 const Vector_n &sol_n) const
520 {
521 using Components = typename Model::Components;
522 static_assert(n_to == Components::count_fe_functions(),
523 "jacobian_numflux_grad: n_to must equal count_fe_functions()");
524 static_assert(n_from == Components::count_fe_functions(),
525 "jacobian_numflux_grad: n_from must equal count_fe_functions()");
526
527 const auto &u_s = get<tup_idx>(sol_s);
528 const auto &u_n = get<tup_idx>(sol_n);
529
530 auto du_s = AD_tools::template ten_to_AD<n_from>(u_s);
531 auto ad_sol_s = std::tuple_cat(tuple_first<tup_idx>(sol_s), std::tie(du_s),
533 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
534 for (uint d = 0; d < dim; ++d) {
535 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
536 seed(du_s[j][d]);
537 asImp().numflux(res, normal, p, Vector_s::as(ad_sol_s), sol_n);
538 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
539 for (uint dd = 0; dd < dim; ++dd) {
540 jNF[0](i, j)[dd][d] = grad(res[i][dd]);
541 }
542 }
543 unseed(du_s[j][d]);
544 }
545 }
546 auto du_n = AD_tools::template ten_to_AD<n_from>(u_n);
547 auto ad_sol_n = std::tuple_cat(tuple_first<tup_idx>(sol_n), std::tie(du_n),
549 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
550 for (uint d = 0; d < dim; ++d) {
551 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
552 seed(du_n[j][d]);
553 asImp().numflux(res, normal, p, sol_s, Vector_n::as(ad_sol_n));
554 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
555 for (uint dd = 0; dd < dim; ++dd) {
556 jNF[1](i, j)[dd][d] = grad(res[i][dd]);
557 }
558 }
559 unseed(du_n[j][d]);
560 }
561 }
562 }
563
564 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
565 void jacobian_numflux_hess(std::array<SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from>, 2> &jNF,
566 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector_s &sol_s,
567 const Vector_n &sol_n) const
568 {
569 using Components = typename Model::Components;
570 static_assert(n_to == Components::count_fe_functions(),
571 "jacobian_numflux_hess: n_to must equal count_fe_functions()");
572 static_assert(n_from == Components::count_fe_functions(),
573 "jacobian_numflux_hess: n_from must equal count_fe_functions()");
574
575 const auto &u_s = get<tup_idx>(sol_s);
576 const auto &u_n = get<tup_idx>(sol_n);
577
578 auto du_s = AD_tools::template ten_to_AD<n_from>(u_s);
579 auto ad_sol_s = std::tuple_cat(tuple_first<tup_idx>(sol_s), std::tie(du_s),
581 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
582 for (uint d1 = 0; d1 < dim; ++d1)
583 for (uint d2 = 0; d2 < dim; ++d2) {
584 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
585 seed(du_s[j][d1][d2]);
586 asImp().numflux(res, normal, p, Vector_s::as(ad_sol_s), sol_n);
587 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
588 for (uint d = 0; d < dim; ++d) {
589 jNF[0](i, j)[d][d1][d2] = grad(res[i][d]);
590 }
591 }
592 unseed(du_s[j][d1][d2]);
593 }
594 }
595 auto du_n = AD_tools::template ten_to_AD<n_from>(u_n);
596 auto ad_sol_n = std::tuple_cat(tuple_first<tup_idx>(sol_n), std::tie(du_n),
598 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
599 for (uint d1 = 0; d1 < dim; ++d1)
600 for (uint d2 = 0; d2 < dim; ++d2) {
601 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
602 seed(du_n[j][d1][d2]);
603 asImp().numflux(res, normal, p, sol_s, Vector_n::as(ad_sol_n));
604 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
605 for (uint d = 0; d < dim; ++d) {
606 jNF[1](i, j)[d][d1][d2] = grad(res[i][d]);
607 }
608 }
609 unseed(du_n[j][d1][d2]);
610 }
611 }
612 }
613
614 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
615 void jacobian_numflux_extr(std::array<SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from>, 2> &jNF,
616 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector_s &sol_s,
617 const Vector_n &sol_n) const
618 {
619 using Components = typename Model::Components;
620 static_assert(n_to == Components::count_fe_functions(),
621 "jacobian_numflux_extr: n_to must equal count_fe_functions()");
622 static_assert(n_from == Components::count_extractors(),
623 "jacobian_numflux_extr: n_from must equal count_extractors()");
624
625 const auto &e = get<tup_idx>(sol_s);
626 auto de = AD_tools::template vector_to_AD<n_from>(e);
627 auto ad_sol_s = std::tuple_cat(tuple_first<tup_idx>(sol_s), std::tie(de),
629 auto ad_sol_n = std::tuple_cat(tuple_first<tup_idx>(sol_n), std::tie(de),
631 for (uint j = 0; j < n_from; ++j) {
632 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res_s{{}};
633 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res_n{{}};
634 seed(de[j]);
635 asImp().numflux(res_s, normal, p, Vector_s::as(ad_sol_s), sol_n);
636 asImp().numflux(res_n, normal, p, sol_s, Vector_n::as(ad_sol_n));
637 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
638 for (uint d = 0; d < dim; ++d) {
639 jNF[0](i, j)[d] = grad(res_s[i][d]);
640 jNF[1](i, j)[d] = grad(res_n[i][d]);
641 }
642 }
643 unseed(de[j]);
644 }
645 }
646
647 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
648 void jacobian_numflux(std::array<SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from>, 2> &jNF,
649 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector_s &sol_s,
650 const Vector_n &sol_n) const
651 {
652 using Components = typename Model::Components;
653 static_assert(n_to == Components::count_fe_functions(to),
654 "jacobian_numflux: n_to must equal count_fe_functions(to)");
655 static_assert(n_from == Components::count_fe_functions(from),
656 "jacobian_numflux: n_from must equal count_fe_functions(from)");
657
658 if constexpr (to == 0) {
659 const auto &u_s = get<from>(sol_s);
660 const auto &u_n = get<from>(sol_n);
661
662 auto du_s = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u_s);
663 auto ad_sol_s = std::tuple_cat(tuple_first<from>(sol_s), std::tie(du_s),
665 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
666 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
667 seed(du_s[j]);
668 asImp().numflux(res, normal, p, Vector_s::as(ad_sol_s), sol_n);
669 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
670 for (uint d = 0; d < dim; ++d) {
671 jNF[0](i, j)[d] = grad(res[i][d]);
672 }
673 }
674 unseed(du_s[j]);
675 }
676 auto du_n = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u_n);
677 auto ad_sol_n = std::tuple_cat(tuple_first<from>(sol_n), std::tie(du_n),
679 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
680 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
681 seed(du_n[j]);
682 asImp().numflux(res, normal, p, sol_s, Vector_n::as(ad_sol_n));
683 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
684 for (uint d = 0; d < dim; ++d) {
685 jNF[1](i, j)[d] = grad(res[i][d]);
686 }
687 }
688 unseed(du_n[j]);
689 }
690 } else {
691 auto du_s = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol_s);
692 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
693 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
694 // take derivative with respect to jth variable
695 seed(du_s[j]);
696 asImp().template ldg_numflux<to>(res, normal, p, du_s, sol_n);
697 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
698 for (uint d = 0; d < dim; ++d) {
699 jNF[0](i, j)[d] = grad(res[i][d]);
700 }
701 }
702 unseed(du_s[j]);
703 }
704 auto du_n = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol_n);
705 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
706 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
707 // take derivative with respect to jth variable
708 seed(du_n[j]);
709 asImp().template ldg_numflux<to>(res, normal, p, sol_s, du_n);
710 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
711 for (uint d = 0; d < dim; ++d) {
712 jNF[1](i, j)[d] = grad(res[i][d]);
713 }
714 }
715 unseed(du_n[j]);
716 }
717 }
718 }
719 };
720
721 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_boundary_numflux
722 {
723 Model &asImp() { return static_cast<Model &>(*this); }
724 const Model &asImp() const { return static_cast<const Model &>(*this); }
726
727 public:
728 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
729 void jacobian_boundary_numflux_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &jBNF,
730 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector &sol) const
731 {
732 using Components = typename Model::Components;
733 static_assert(n_to == Components::count_fe_functions(),
734 "jacobian_boundary_numflux_grad: n_to must equal count_fe_functions()");
735 static_assert(n_from == Components::count_fe_functions(),
736 "jacobian_boundary_numflux_grad: n_from must equal count_fe_functions()");
737
738 const auto &u = get<tup_idx>(sol);
739 auto du = AD_tools::template ten_to_AD<n_from>(u);
740 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
742 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
743 for (uint d = 0; d < dim; ++d) {
744 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
745 seed(du[j][d]);
746 asImp().boundary_numflux(res, normal, p, Vector::as(ad_sol));
747 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
748 for (uint dd = 0; dd < dim; ++dd) {
749 jBNF(i, j)[dd][d] = grad(res[i][dd]);
750 }
751 }
752 unseed(du[j][d]);
753 }
754 }
755 }
756
757 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
758 void jacobian_boundary_numflux_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &jBNF,
759 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector &sol) const
760 {
761 using Components = typename Model::Components;
762 static_assert(n_to == Components::count_fe_functions(),
763 "jacobian_boundary_numflux_hess: n_to must equal count_fe_functions()");
764 static_assert(n_from == Components::count_fe_functions(),
765 "jacobian_boundary_numflux_hess: n_from must equal count_fe_functions()");
766
767 const auto &u = get<tup_idx>(sol);
768 auto du = AD_tools::template ten_to_AD<n_from>(u);
769 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(du),
771 for (uint j = 0; j < Components::count_fe_functions(); ++j) {
772 for (uint d1 = 0; d1 < dim; ++d1)
773 for (uint d2 = 0; d2 < dim; ++d2) {
774 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
775 seed(du[j][d1][d2]);
776 asImp().boundary_numflux(res, normal, p, Vector::as(ad_sol));
777 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
778 for (uint d = 0; d < dim; ++d) {
779 jBNF(i, j)[d][d1][d2] = grad(res[i][d]);
780 }
781 }
782 unseed(du[j][d1][d2]);
783 }
784 }
785 }
786
787 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
788 void jacobian_boundary_numflux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jBNF,
789 const Tensor<1, dim> &normal, const Point<dim> &p, const Vector &sol) const
790 {
791 using Components = typename Model::Components;
792 static_assert(n_to == Components::count_fe_functions(),
793 "jacobian_boundary_numflux_extr: n_to must equal count_fe_functions()");
794 static_assert(n_from == Components::count_extractors(),
795 "jacobian_boundary_numflux_extr: n_from must equal count_extractors()");
796
797 const auto &e = get<tup_idx>(sol);
798 auto de = AD_tools::template vector_to_AD<n_from>(e);
799 auto ad_sol = std::tuple_cat(tuple_first<tup_idx>(sol), std::tie(de),
801 for (uint j = 0; j < n_from; ++j) {
802 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
803 seed(de[j]);
804 asImp().boundary_numflux(res, normal, p, Vector::as(ad_sol));
805 for (uint i = 0; i < Components::count_fe_functions(); ++i) {
806 for (uint d = 0; d < dim; ++d) {
807 jBNF(i, j)[d] = grad(res[i][d]);
808 }
809 }
810 unseed(de[j]);
811 }
812 }
813
814 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
815 void jacobian_boundary_numflux(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jBNF, const Tensor<1, dim> &normal,
816 const Point<dim> &p, const Vector &sol) const
817 {
818 using Components = typename Model::Components;
819 static_assert(n_to == Components::count_fe_functions(to),
820 "jacobian_boundary_numflux: n_to must equal count_fe_functions(to)");
821 static_assert(n_from == Components::count_fe_functions(from),
822 "jacobian_boundary_numflux: n_from must equal count_fe_functions(from)");
823
824 if constexpr (to == 0) {
825 const auto &u = get<from>(sol);
826 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
827 auto ad_sol = std::tuple_cat(tuple_first<from>(sol), std::tie(du),
829 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
830 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
831 seed(du[j]);
832 asImp().boundary_numflux(res, normal, p, Vector::as(ad_sol));
833 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
834 for (uint d = 0; d < dim; ++d) {
835 jBNF(i, j)[d] = grad(res[i][d]);
836 }
837 }
838 unseed(du[j]);
839 }
840 } else {
841 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(sol);
842 for (uint j = 0; j < Components::count_fe_functions(from); ++j) {
843 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions(to)> res{{}};
844 // take derivative with respect to jth variable
845 seed(du[j]);
846 asImp().template ldg_boundary_numflux<to>(res, normal, p, du);
847 for (uint i = 0; i < Components::count_fe_functions(to); ++i) {
848 for (uint d = 0; d < dim; ++d) {
849 jBNF(i, j)[d] = grad(res[i][d]);
850 }
851 }
852 unseed(du[j]);
853 }
854 }
855 }
856 };
857
858 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_mass
859 {
860 Model &asImp() { return static_cast<Model &>(*this); }
861 const Model &asImp() const { return static_cast<const Model &>(*this); }
863
864 public:
865 template <uint dot, uint n_from, uint n_to, int dim, typename NT, typename Vector>
866 void jacobian_mass(SimpleMatrix<NT, n_to, n_from> &jM, const Point<dim> &p, const Vector &u,
867 const Vector &u_dot) const
868 {
869 using Components = typename Model::Components;
870 static_assert(n_from == Components::count_fe_functions() && n_to == Components::count_fe_functions(),
871 "jacobian_mass: n_from and n_to must both equal count_fe_functions()");
872
873 if constexpr (dot == 0) {
874 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(0)>(u);
875 for (uint j = 0; j < Components::count_fe_functions(0); ++j) {
876 std::array<AD_type, Components::count_fe_functions(0)> res{{}};
877 // take derivative with respect to jth variable
878 seed(du[j]);
879 asImp().mass(res, p, du, u_dot);
880 for (uint i = 0; i < Components::count_fe_functions(0); ++i) {
881 jM(i, j) = grad(res[i]);
882 }
883 unseed(du[j]);
884 }
885 } else {
886 auto du_dot = AD_tools::template vector_to_AD<Components::count_fe_functions(0)>(u_dot);
887 for (uint j = 0; j < Components::count_fe_functions(0); ++j) {
888 std::array<AD_type, Components::count_fe_functions(0)> res{{}};
889 // take derivative with respect to jth variable
890 seed(du_dot[j]);
891 asImp().mass(res, p, u, du_dot);
892 for (uint i = 0; i < Components::count_fe_functions(0); ++i) {
893 jM(i, j) = grad(res[i]);
894 }
895 unseed(du_dot[j]);
896 }
897 }
898 }
899 };
900
901 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_variables
902 {
903 Model &asImp() { return static_cast<Model &>(*this); }
904 const Model &asImp() const { return static_cast<const Model &>(*this); }
906
907 public:
908 template <uint to, typename NT, typename Solution>
909 void jacobian_variables(FullMatrix<NT> &jac, const Solution &sol) const
910 {
911 const auto &variables = get<0>(sol);
912 const auto &extractors = get<1>(sol);
913 if constexpr (to == 0) {
914 AssertThrow(jac.m() == variables.size() && jac.n() == variables.size(),
915 ExcMessage("Assure that the jacobian has the right dimension!"));
916 } else if constexpr (to == 1) {
917 AssertThrow(jac.m() == variables.size() && jac.n() == extractors.size(),
918 ExcMessage("Assure that the jacobian has the right dimension!"));
919 }
920
921 if constexpr (to == 0) {
922 auto du = AD_tools::template vector_to_AD<Model::Components::count_variables()>(variables);
923 auto ad_sol = std::tuple_cat(tuple_first<to>(sol), std::tie(du),
925 for (uint j = 0; j < Model::Components::count_variables(); ++j) {
926 std::array<AD_type, Model::Components::count_variables()> res{{}};
927 seed(du[j]);
928 asImp().dt_variables(res, Solution::as(ad_sol));
929 for (uint i = 0; i < Model::Components::count_variables(); ++i) {
930 jac(i, j) = grad(res[i]);
931 }
932 unseed(du[j]);
933 }
934 } else if constexpr (to == 1) {
935 auto du = AD_tools::template vector_to_AD<Model::Components::count_extractors()>(extractors);
936 auto ad_sol = std::tuple_cat(tuple_first<to>(sol), std::tie(du),
938 for (uint j = 0; j < Model::Components::count_extractors(); ++j) {
939 std::array<AD_type, Model::Components::count_variables()> res{{}};
940 seed(du[j]);
941 asImp().dt_variables(res, Solution::as(ad_sol));
942 for (uint i = 0; i < Model::Components::count_extractors(); ++i) {
943 jac(i, j) = grad(res[i]);
944 }
945 unseed(du[j]);
946 }
947 }
948 }
949 };
950
951 template <typename Model, typename AD_type = autodiff::real> class ADjacobian_extractors
952 {
953 Model &asImp() { return static_cast<Model &>(*this); }
954 const Model &asImp() const { return static_cast<const Model &>(*this); }
956
957 public:
958 template <uint to, int dim, typename NT, typename Solution>
959 void jacobian_extractors(FullMatrix<NT> &jac, const Point<dim> &x, const Solution &sol) const
960 {
961 static_assert(std::is_same_v<NT, double>, "Only double is supported for now!");
962 const auto &fe_functions = get<0>(sol);
963 const auto &fe_derivatives = get<1>(sol);
964 const auto &fe_hessians = get<2>(sol);
965
966 if constexpr (to == 0) {
967 AssertThrow(jac.m() == Model::Components::count_extractors() && jac.n() == fe_functions.size(),
968 ExcMessage("Assure that the jacobian has the right dimension!"));
969 } else if constexpr (to == 1) {
970 AssertThrow(jac.m() == Model::Components::count_extractors() && jac.n() == fe_derivatives.size() * dim,
971 ExcMessage("Assure that the jacobian has the right dimension!"));
972 } else if constexpr (to == 2) {
973 AssertThrow(jac.m() == Model::Components::count_extractors() && jac.n() == fe_derivatives.size() * dim * dim,
974 ExcMessage("Assure that the jacobian has the right dimension!"));
975 }
976
977 if constexpr (to == 0) {
978 auto du = AD_tools::template vector_to_AD<Model::Components::count_fe_functions()>(fe_functions);
979 auto ad_sol = std::tuple_cat(tuple_first<to>(sol), std::tie(du),
981 for (uint j = 0; j < Model::Components::count_fe_functions(); ++j) {
982 std::array<AD_type, Model::Components::count_extractors()> res{{}};
983 seed(du[j]);
984 asImp().extract(res, x, Solution::as(ad_sol));
985 for (uint i = 0; i < Model::Components::count_extractors(); ++i) {
986 jac(i, j) = grad(res[i]);
987 }
988 unseed(du[j]);
989 }
990 } else if constexpr (to == 1) {
991 auto du = AD_tools::template ten_to_AD<Model::Components::count_fe_functions()>(fe_derivatives);
992 auto ad_sol = std::tuple_cat(tuple_first<to>(sol), std::tie(du),
994 for (uint j = 0; j < Model::Components::count_fe_functions(); ++j) {
995 for (uint d1 = 0; d1 < dim; ++d1) {
996 std::array<AD_type, Model::Components::count_extractors()> res{{}};
997 seed(du[j][d1]);
998 asImp().extract(res, x, Solution::as(ad_sol));
999 for (uint i = 0; i < Model::Components::count_extractors(); ++i) {
1000 jac(i, j * dim + d1) = grad(res[i]);
1001 }
1002 unseed(du[j][d1]);
1003 }
1004 }
1005 } else if constexpr (to == 2) {
1006 auto du = AD_tools::template ten_to_AD<Model::Components::count_fe_functions()>(fe_hessians);
1007 auto ad_sol = std::tuple_cat(tuple_first<to>(sol), std::tie(du),
1009 for (uint j = 0; j < Model::Components::count_fe_functions(); ++j) {
1010 for (uint d1 = 0; d1 < dim; ++d1)
1011 for (uint d2 = 0; d2 < dim; ++d2) {
1012 std::array<AD_type, Model::Components::count_extractors()> res{{}};
1013 seed(du[j][d1][d2]);
1014 asImp().extract(res, x, Solution::as(ad_sol));
1015 for (uint i = 0; i < Model::Components::count_extractors(); ++i) {
1016 jac(i, j * dim * dim + d1 * dim + d2) = grad(res[i]);
1017 }
1018 unseed(du[j][d1][d2]);
1019 }
1020 }
1021 }
1022 }
1023 };
1024
1025 template <typename Model>
1026 class AD_real : public ADjacobian_flux<Model, autodiff::real>,
1027 public ADjacobian_source<Model, autodiff::real>,
1028 public ADjacobian_flux_source<Model, autodiff::real>,
1029 public ADjacobian_numflux<Model, autodiff::real>,
1030 public ADjacobian_boundary_numflux<Model, autodiff::real>,
1031 public ADjacobian_mass<Model, autodiff::real>,
1032 public ADjacobian_variables<Model, autodiff::real>,
1033 public ADjacobian_extractors<Model, autodiff::real>
1034 {
1035 };
1036
1037 template <typename Model>
1038 class AD_dual : public ADjacobian_flux<Model, autodiff::dual>,
1039 public ADjacobian_source<Model, autodiff::dual>,
1040 public ADjacobian_flux_source<Model, autodiff::dual>,
1041 public ADjacobian_numflux<Model, autodiff::dual>,
1042 public ADjacobian_boundary_numflux<Model, autodiff::dual>,
1043 public ADjacobian_mass<Model, autodiff::dual>,
1044 public ADjacobian_variables<Model, autodiff::dual>,
1045 public ADjacobian_extractors<Model, autodiff::dual>
1046 {
1047 };
1048
1049 template <typename Model> using AD = AD_real<Model>;
1050
1051 template <typename Model>
1052 class FE_AD : public ADjacobian_flux<Model, autodiff::real>,
1053 public ADjacobian_source<Model, autodiff::real>,
1054 public ADjacobian_flux_source<Model, autodiff::real>,
1055 public ADjacobian_numflux<Model, autodiff::real>,
1056 public ADjacobian_boundary_numflux<Model, autodiff::real>,
1057 public ADjacobian_mass<Model, autodiff::real>
1058 {
1059 public:
1060 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1061 void jacobian_flux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Point<dim> &,
1062 const Vector &) const
1063 {
1064 }
1065 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1066 void jacobian_source_extr(SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &) const
1067 {
1068 }
1069 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1070 void jacobian_flux_source_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &,
1071 SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &) const
1072 {
1073 }
1074 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
1075 void jacobian_numflux_extr(std::array<SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from>, 2> &,
1076 const Tensor<1, dim> &, const Point<dim> &, const Vector_s &, const Vector_n &) const
1077 {
1078 }
1079 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1080 void jacobian_boundary_numflux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Tensor<1, dim> &,
1081 const Point<dim> &, const Vector &) const
1082 {
1083 }
1084 template <uint to, int dim, typename NT, typename Solution>
1085 void jacobian_extractors(FullMatrix<NT> &, const Point<dim> &, const Solution &) const
1086 {
1087 }
1088 template <uint to, typename NT, typename Solution>
1089 void jacobian_variables(FullMatrix<NT> &, const Solution &) const
1090 {
1091 }
1092 };
1093
1095 {
1096 public:
1097 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1098 void jacobian_flux_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &, const Point<dim> &,
1099 const Vector &) const
1100 {
1101 }
1102 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1103 void jacobian_flux_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &, const Point<dim> &,
1104 const Vector &) const
1105 {
1106 }
1107 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1108 void jacobian_flux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Point<dim> &,
1109 const Vector &) const
1110 {
1111 }
1112 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1113 void jacobian_flux(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Point<dim> &, const Vector &) const
1114 {
1115 }
1116 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1117 void jacobian_source_grad(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Point<dim> &,
1118 const Vector &) const
1119 {
1120 }
1121 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1122 void jacobian_source_hess(SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &, const Point<dim> &,
1123 const Vector &) const
1124 {
1125 }
1126 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1127 void jacobian_source_extr(SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &) const
1128 {
1129 }
1130 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1131 void jacobian_source(SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &) const
1132 {
1133 }
1134 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1135 void jacobian_flux_source(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, SimpleMatrix<NT, n_to, n_from> &,
1136 const Point<dim> &, const Vector &) const
1137 {
1138 }
1139 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1140 void jacobian_flux_source_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &,
1141 SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Point<dim> &,
1142 const Vector &) const
1143 {
1144 }
1145 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1146 void jacobian_flux_source_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &,
1147 SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &, const Point<dim> &,
1148 const Vector &) const
1149 {
1150 }
1151 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1152 void jacobian_flux_source_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &,
1153 SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &) const
1154 {
1155 }
1156 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
1157 void jacobian_numflux_grad(std::array<SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from>, 2> &,
1158 const Tensor<1, dim> &, const Point<dim> &, const Vector_s &, const Vector_n &) const
1159 {
1160 }
1161 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
1162 void jacobian_numflux_hess(std::array<SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from>, 2> &,
1163 const Tensor<1, dim> &, const Point<dim> &, const Vector_s &, const Vector_n &) const
1164 {
1165 }
1166 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
1167 void jacobian_numflux_extr(std::array<SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from>, 2> &,
1168 const Tensor<1, dim> &, const Point<dim> &, const Vector_s &, const Vector_n &) const
1169 {
1170 }
1171 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector_s, typename Vector_n>
1172 void jacobian_numflux(std::array<SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from>, 2> &, const Tensor<1, dim> &,
1173 const Point<dim> &, const Vector_s &, const Vector_n &) const
1174 {
1175 }
1176 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1177 void jacobian_boundary_numflux_grad(SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NT>>, n_to, n_from> &,
1178 const Tensor<1, dim> &, const Point<dim> &, const Vector &) const
1179 {
1180 }
1181 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1182 void jacobian_boundary_numflux_hess(SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NT>>, n_to, n_from> &,
1183 const Tensor<1, dim> &, const Point<dim> &, const Vector &) const
1184 {
1185 }
1186 template <uint tup_idx, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1187 void jacobian_boundary_numflux_extr(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Tensor<1, dim> &,
1188 const Point<dim> &, const Vector &) const
1189 {
1190 }
1191 template <uint from, uint to, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1192 void jacobian_boundary_numflux(SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &, const Tensor<1, dim> &,
1193 const Point<dim> &, const Vector &) const
1194 {
1195 }
1196 template <uint dot, uint n_from, uint n_to, int dim, typename NT, typename Vector>
1197 void jacobian_mass(SimpleMatrix<NT, n_to, n_from> &, const Point<dim> &, const Vector &, const Vector &) const
1198 {
1199 }
1200
1201 template <uint to, int dim, typename NT, typename Solution>
1202 void jacobian_extractors(FullMatrix<NT> &, const Point<dim> &, const Solution &) const
1203 {
1204 }
1205 template <uint to, typename NT, typename Solution>
1206 void jacobian_variables(FullMatrix<NT> &, const Solution &) const
1207 {
1208 }
1209 };
1210
1211 } // namespace def
1212} // namespace DiFfRG
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
Definition ad.hh:1046
Definition ad.hh:1034
void jacobian_boundary_numflux_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &jBNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:729
void jacobian_boundary_numflux_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &jBNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:758
const Model & asImp() const
Definition ad.hh:724
void jacobian_boundary_numflux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jBNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:788
Model & asImp()
Definition ad.hh:723
void jacobian_boundary_numflux(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jBNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:815
Model & asImp()
Definition ad.hh:953
void jacobian_extractors(FullMatrix< NT > &jac, const Point< dim > &x, const Solution &sol) const
Definition ad.hh:959
const Model & asImp() const
Definition ad.hh:954
Model & asImp()
Definition ad.hh:357
void jacobian_flux_source_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &jF, SimpleMatrix< Tensor< 2, dim, NT >, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:443
void jacobian_flux_source(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jF, SimpleMatrix< NT, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:363
void jacobian_flux_source_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jF, SimpleMatrix< NT, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:479
void jacobian_flux_source_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &jF, SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:409
const Model & asImp() const
Definition ad.hh:358
Definition ad.hh:89
void jacobian_flux_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &jF, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:96
void jacobian_flux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jF, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:157
void jacobian_flux(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jF, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:184
void jacobian_flux_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &jF, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:126
const Model & asImp() const
Definition ad.hh:91
Model & asImp()
Definition ad.hh:90
Definition ad.hh:859
void jacobian_mass(SimpleMatrix< NT, n_to, n_from > &jM, const Point< dim > &p, const Vector &u, const Vector &u_dot) const
Definition ad.hh:866
Model & asImp()
Definition ad.hh:860
const Model & asImp() const
Definition ad.hh:861
const Model & asImp() const
Definition ad.hh:512
void jacobian_numflux_extr(std::array< SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from >, 2 > &jNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector_s &sol_s, const Vector_n &sol_n) const
Definition ad.hh:615
Model & asImp()
Definition ad.hh:511
void jacobian_numflux(std::array< SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from >, 2 > &jNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector_s &sol_s, const Vector_n &sol_n) const
Definition ad.hh:648
void jacobian_numflux_grad(std::array< SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from >, 2 > &jNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector_s &sol_s, const Vector_n &sol_n) const
Definition ad.hh:517
void jacobian_numflux_hess(std::array< SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from >, 2 > &jNF, const Tensor< 1, dim > &normal, const Point< dim > &p, const Vector_s &sol_s, const Vector_n &sol_n) const
Definition ad.hh:565
Definition ad.hh:228
void jacobian_source_grad(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:235
void jacobian_source(SimpleMatrix< NT, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:317
Model & asImp()
Definition ad.hh:229
void jacobian_source_extr(SimpleMatrix< NT, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:293
const Model & asImp() const
Definition ad.hh:230
void jacobian_source_hess(SimpleMatrix< Tensor< 2, dim, NT >, n_to, n_from > &jS, const Point< dim > &p, const Vector &sol) const
Definition ad.hh:263
Model & asImp()
Definition ad.hh:903
void jacobian_variables(FullMatrix< NT > &jac, const Solution &sol) const
Definition ad.hh:909
const Model & asImp() const
Definition ad.hh:904
Definition ad.hh:1058
void jacobian_boundary_numflux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1080
void jacobian_extractors(FullMatrix< NT > &, const Point< dim > &, const Solution &) const
Definition ad.hh:1085
void jacobian_flux_source_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1070
void jacobian_flux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1061
void jacobian_variables(FullMatrix< NT > &, const Solution &) const
Definition ad.hh:1089
void jacobian_source_extr(SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1066
void jacobian_numflux_extr(std::array< SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from >, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector_s &, const Vector_n &) const
Definition ad.hh:1075
Definition ad.hh:1095
void jacobian_boundary_numflux_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1182
void jacobian_boundary_numflux(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1192
void jacobian_flux_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1103
void jacobian_mass(SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &, const Vector &) const
Definition ad.hh:1197
void jacobian_numflux_hess(std::array< SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from >, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector_s &, const Vector_n &) const
Definition ad.hh:1162
void jacobian_flux_source_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1152
void jacobian_source(SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1131
void jacobian_flux_source_hess(SimpleMatrix< Tensor< 1, dim, Tensor< 2, dim, NT > >, n_to, n_from > &, SimpleMatrix< Tensor< 2, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1146
void jacobian_source_grad(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1117
void jacobian_variables(FullMatrix< NT > &, const Solution &) const
Definition ad.hh:1206
void jacobian_source_hess(SimpleMatrix< Tensor< 2, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1122
void jacobian_boundary_numflux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1187
void jacobian_flux_extr(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1108
void jacobian_extractors(FullMatrix< NT > &, const Point< dim > &, const Solution &) const
Definition ad.hh:1202
void jacobian_numflux_grad(std::array< SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from >, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector_s &, const Vector_n &) const
Definition ad.hh:1157
void jacobian_flux_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1098
void jacobian_source_extr(SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1127
void jacobian_numflux_extr(std::array< SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from >, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector_s &, const Vector_n &) const
Definition ad.hh:1167
void jacobian_flux(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1113
void jacobian_boundary_numflux_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1177
void jacobian_flux_source(SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, SimpleMatrix< NT, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1135
void jacobian_numflux(std::array< SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from >, 2 > &, const Tensor< 1, dim > &, const Point< dim > &, const Vector_s &, const Vector_n &) const
Definition ad.hh:1172
void jacobian_flux_source_grad(SimpleMatrix< Tensor< 1, dim, Tensor< 1, dim, NT > >, n_to, n_from > &, SimpleMatrix< Tensor< 1, dim, NT >, n_to, n_from > &, const Point< dim > &, const Vector &) const
Definition ad.hh:1140
Definition complex_math.hh:10
constexpr auto tuple_first(const std::tuple< Head, Tail... > &t)
Definition tuples.hh:247
NT dot(const A1 &a1, const A2 &a2)
A dot product which takes the dot product between a1 and a2, assuming each has n entries which can be...
Definition math.hh:211
constexpr auto & get(named_tuple< tuple_type, strSet > &ob)
get a reference to the element with the given name
Definition tuples.hh:111
unsigned int uint
Definition utils.hh:24
KOKKOS_FORCEINLINE_FUNCTION auto real(const autodiff::Real< N, T > &a)
Definition complex_math.hh:96
constexpr auto tuple_last(const std::tuple< Head, Tail... > &t)
Definition tuples.hh:233
Definition complex_math.hh:19
static std::array< dealii::Tensor< r, dim, autodiff::dual >, n > ten_to_AD(const std::vector< dealii::Tensor< r, dim, NT > > &v)
Definition ad.hh:39
static std::array< autodiff::dual, n > vector_to_AD(const Vector &v)
Definition ad.hh:28
static std::array< autodiff::real, n > vector_to_AD(const Vector &v)
Definition ad.hh:58
static std::array< dealii::Tensor< r, dim, autodiff::real >, n > ten_to_AD(const std::vector< dealii::Tensor< r, dim, NT > > &v)
Definition ad.hh:68