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>
20 using namespace dealii;
28 template <u
int n,
typename Vector>
static std::array<autodiff::dual, n>
vector_to_AD(
const Vector &v)
30 std::array<autodiff::dual, n> x;
31 for (
uint i = 0; i < n; ++i) {
37 template <u
int 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)
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) {
45 for (
uint d = 0; d < dim; ++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];
58 template <u
int n,
typename Vector>
static std::array<autodiff::real, n>
vector_to_AD(
const Vector &v)
60 std::array<autodiff::real, n> x;
61 for (
uint i = 0; i < n; ++i)
66 template <u
int 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)
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)
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];
88 template <
typename Model,
typename AD_type = autodiff::real>
class ADjacobian_flux
90 Model &
asImp() {
return static_cast<Model &
>(*this); }
91 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
95 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
97 const Vector &sol)
const
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()");
106 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
125 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
127 const Vector &sol)
const
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()");
136 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
151 unseed(du[j][d1][d2]);
156 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
158 const Vector &sol)
const
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()");
167 auto de = AD_tools::template vector_to_AD<n_from>(e);
170 for (
uint j = 0; j < Components::count_extractors(); ++j) {
171 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
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]);
183 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
185 const Vector &sol)
const
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)");
193 if constexpr (to == 0) {
195 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
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{{}};
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]);
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{{}};
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]);
229 Model &
asImp() {
return static_cast<Model &
>(*this); }
230 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
234 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
236 const Vector &sol)
const
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()");
245 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
262 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
264 const Vector &sol)
const
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()");
273 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
286 unseed(du[j][d1][d2]);
292 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
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()");
302 auto de = AD_tools::template vector_to_AD<n_from>(e);
305 for (
uint j = 0; j < n_from; ++j) {
306 std::array<AD_type, Components::count_fe_functions()> res{{}};
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]);
316 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
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)");
325 if constexpr (to == 0) {
327 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
330 for (
uint j = 0; j < Components::count_fe_functions(from); ++j) {
331 std::array<AD_type, Components::count_fe_functions(to)> res{{}};
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]);
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{{}};
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]);
357 Model &
asImp() {
return static_cast<Model &
>(*this); }
358 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
362 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
364 const Point<dim> &p,
const Vector &sol)
const
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)");
372 if constexpr (to == 0) {
374 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
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{{}};
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]);
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{{}};
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]);
408 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
410 SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &jS,
const Point<dim> &p,
411 const Vector &sol)
const
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()");
420 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
442 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
444 SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &jS,
const Point<dim> &p,
445 const Vector &sol)
const
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()");
454 auto du = AD_tools::template ten_to_AD<n_from>(u);
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) {
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]);
472 unseed(du[j][d1][d2]);
478 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
481 const Vector &sol)
const
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()");
490 auto de = AD_tools::template vector_to_AD<n_from>(e);
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{{}};
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]);
511 Model &
asImp() {
return static_cast<Model &
>(*this); }
512 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
516 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
518 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector_s &sol_s,
519 const Vector_n &sol_n)
const
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()");
530 auto du_s = AD_tools::template ten_to_AD<n_from>(u_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{{}};
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]);
546 auto du_n = AD_tools::template ten_to_AD<n_from>(u_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{{}};
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]);
564 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
566 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector_s &sol_s,
567 const Vector_n &sol_n)
const
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()");
578 auto du_s = AD_tools::template ten_to_AD<n_from>(u_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]);
592 unseed(du_s[j][d1][d2]);
595 auto du_n = AD_tools::template ten_to_AD<n_from>(u_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]);
609 unseed(du_n[j][d1][d2]);
614 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
616 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector_s &sol_s,
617 const Vector_n &sol_n)
const
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()");
626 auto de = AD_tools::template vector_to_AD<n_from>(e);
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{{}};
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]);
647 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
649 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector_s &sol_s,
650 const Vector_n &sol_n)
const
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)");
658 if constexpr (to == 0) {
662 auto du_s = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u_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{{}};
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]);
676 auto du_n = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u_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{{}};
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]);
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{{}};
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]);
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{{}};
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]);
723 Model &
asImp() {
return static_cast<Model &
>(*this); }
724 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
728 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
730 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector &sol)
const
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()");
739 auto du = AD_tools::template ten_to_AD<n_from>(u);
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{{}};
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]);
757 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
759 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector &sol)
const
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()");
768 auto du = AD_tools::template ten_to_AD<n_from>(u);
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{{}};
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]);
782 unseed(du[j][d1][d2]);
787 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
789 const Tensor<1, dim> &normal,
const Point<dim> &p,
const Vector &sol)
const
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()");
798 auto de = AD_tools::template vector_to_AD<n_from>(e);
801 for (
uint j = 0; j < n_from; ++j) {
802 std::array<Tensor<1, dim, AD_type>, Components::count_fe_functions()> res{{}};
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]);
814 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
816 const Point<dim> &p,
const Vector &sol)
const
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)");
824 if constexpr (to == 0) {
826 auto du = AD_tools::template vector_to_AD<Components::count_fe_functions(from)>(u);
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{{}};
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]);
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{{}};
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]);
860 Model &
asImp() {
return static_cast<Model &
>(*this); }
861 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
865 template <u
int dot, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
867 const Vector &u_dot)
const
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()");
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{{}};
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]);
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{{}};
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]);
903 Model &
asImp() {
return static_cast<Model &
>(*this); }
904 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
908 template <u
int to,
typename NT,
typename Solution>
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!"));
921 if constexpr (to == 0) {
922 auto du = AD_tools::template vector_to_AD<Model::Components::count_variables()>(variables);
925 for (
uint j = 0; j < Model::Components::count_variables(); ++j) {
926 std::array<AD_type, Model::Components::count_variables()> res{{}};
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]);
934 }
else if constexpr (to == 1) {
935 auto du = AD_tools::template vector_to_AD<Model::Components::count_extractors()>(extractors);
938 for (
uint j = 0; j < Model::Components::count_extractors(); ++j) {
939 std::array<AD_type, Model::Components::count_variables()> res{{}};
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]);
953 Model &
asImp() {
return static_cast<Model &
>(*this); }
954 const Model &
asImp()
const {
return static_cast<const Model &
>(*this); }
958 template <u
int to,
int dim,
typename NT,
typename Solution>
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);
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!"));
977 if constexpr (to == 0) {
978 auto du = AD_tools::template vector_to_AD<Model::Components::count_fe_functions()>(fe_functions);
981 for (
uint j = 0; j < Model::Components::count_fe_functions(); ++j) {
982 std::array<AD_type, Model::Components::count_extractors()> res{{}};
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]);
990 }
else if constexpr (to == 1) {
991 auto du = AD_tools::template ten_to_AD<Model::Components::count_fe_functions()>(fe_derivatives);
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{{}};
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]);
1005 }
else if constexpr (to == 2) {
1006 auto du = AD_tools::template ten_to_AD<Model::Components::count_fe_functions()>(fe_hessians);
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]);
1018 unseed(du[j][d1][d2]);
1025 template <
typename Model>
1037 template <
typename Model>
1051 template <
typename Model>
1060 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1062 const Vector &)
const
1065 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1069 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1074 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
1076 const Tensor<1, dim> &,
const Point<dim> &,
const Vector_s &,
const Vector_n &)
const
1079 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1081 const Point<dim> &,
const Vector &)
const
1084 template <u
int to,
int dim,
typename NT,
typename Solution>
1088 template <u
int to,
typename NT,
typename Solution>
1097 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1099 const Vector &)
const
1102 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1104 const Vector &)
const
1107 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1109 const Vector &)
const
1112 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1116 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1118 const Vector &)
const
1121 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1123 const Vector &)
const
1126 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1130 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1134 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1136 const Point<dim> &,
const Vector &)
const
1139 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1141 SimpleMatrix<Tensor<1, dim, NT>, n_to, n_from> &,
const Point<dim> &,
1142 const Vector &)
const
1145 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1147 SimpleMatrix<Tensor<2, dim, NT>, n_to, n_from> &,
const Point<dim> &,
1148 const Vector &)
const
1151 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1156 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
1158 const Tensor<1, dim> &,
const Point<dim> &,
const Vector_s &,
const Vector_n &)
const
1161 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
1163 const Tensor<1, dim> &,
const Point<dim> &,
const Vector_s &,
const Vector_n &)
const
1166 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
1168 const Tensor<1, dim> &,
const Point<dim> &,
const Vector_s &,
const Vector_n &)
const
1171 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector_s,
typename Vector_n>
1173 const Point<dim> &,
const Vector_s &,
const Vector_n &)
const
1176 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1178 const Tensor<1, dim> &,
const Point<dim> &,
const Vector &)
const
1181 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1183 const Tensor<1, dim> &,
const Point<dim> &,
const Vector &)
const
1186 template <u
int tup_
idx, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1188 const Point<dim> &,
const Vector &)
const
1191 template <u
int from, u
int to, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1193 const Point<dim> &,
const Vector &)
const
1196 template <u
int dot, u
int n_from, u
int n_to,
int dim,
typename NT,
typename Vector>
1201 template <u
int to,
int dim,
typename NT,
typename Solution>
1205 template <u
int to,
typename NT,
typename Solution>
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
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: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
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
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
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
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
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