Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
tuple_tensor_dual_functions.hpp
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Serac Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
7 #pragma once
8 
9 #include "mfem.hpp"
10 
14 
15 namespace serac {
16 
18 template <typename T>
20  static constexpr bool value = false;
21 };
22 
24 template <typename T, int... n>
26  static constexpr bool value = true;
27 };
28 
37 template <typename S, typename T, int m, int... n,
38  typename = std::enable_if_t<std::is_arithmetic_v<S> || is_dual_number<S>::value>>
39 SERAC_HOST_DEVICE constexpr auto operator*(S scale, const tensor<T, m, n...>& A)
40 {
41  tensor<decltype(S{} * T{}), m, n...> C{};
42  for (int i = 0; i < m; i++) {
43  C[i] = scale * A[i];
44  }
45  return C;
46 }
47 
56 template <typename S, typename T, int m, int... n,
57  typename = std::enable_if_t<std::is_arithmetic_v<S> || is_dual_number<S>::value>>
58 SERAC_HOST_DEVICE constexpr auto operator*(const tensor<T, m, n...>& A, S scale)
59 {
60  tensor<decltype(T{} * S{}), m, n...> C{};
61  for (int i = 0; i < m; i++) {
62  C[i] = A[i] * scale;
63  }
64  return C;
65 }
66 
75 template <typename S, typename T, int m, int... n,
76  typename = std::enable_if_t<std::is_arithmetic_v<S> || is_dual_number<S>::value>>
77 SERAC_HOST_DEVICE constexpr auto operator/(S scale, const tensor<T, m, n...>& A)
78 {
79  tensor<decltype(S{} * T{}), n...> C{};
80  for (int i = 0; i < m; i++) {
81  C[i] = scale / A[i];
82  }
83  return C;
84 }
85 
94 template <typename S, typename T, int m, int... n,
95  typename = std::enable_if_t<std::is_arithmetic_v<S> || is_dual_number<S>::value>>
96 SERAC_HOST_DEVICE constexpr auto operator/(const tensor<T, m, n...>& A, S scale)
97 {
98  tensor<decltype(T{} * S{}), m, n...> C{};
99  for (int i = 0; i < m; i++) {
100  C[i] = A[i] / scale;
101  }
102  return C;
103 }
104 
106 template <int i, typename S, typename T>
107 struct one_hot_helper;
108 
109 template <int i, int... I, typename T>
110 struct one_hot_helper<i, std::integer_sequence<int, I...>, T> {
111  using type = tuple<std::conditional_t<i == I, T, zero>...>;
112 };
113 
114 template <int i, int n, typename T>
115 struct one_hot : public one_hot_helper<i, std::make_integer_sequence<int, n>, T> {
116 };
118 
125 template <int i, int n, typename T>
127 
129 template <int i, int N>
130 SERAC_HOST_DEVICE constexpr auto make_dual_helper(zero /*arg*/)
131 {
132  return zero{};
133 }
134 
142 template <int i, int N>
143 SERAC_HOST_DEVICE constexpr auto make_dual_helper(double arg)
144 {
145  using gradient_t = one_hot_t<i, N, double>;
146  dual<gradient_t> arg_dual{};
147  arg_dual.value = arg;
148  serac::get<i>(arg_dual.gradient) = 1.0;
149  return arg_dual;
150 }
151 
159 template <int i, int N, typename T, int... n>
161 {
162  using gradient_t = one_hot_t<i, N, tensor<T, n...>>;
163  tensor<dual<gradient_t>, n...> arg_dual{};
164  for_constexpr<n...>([&](auto... j) {
165  arg_dual(j...).value = arg(j...);
166  serac::get<i>(arg_dual(j...).gradient)(j...) = 1.0;
167  });
168  return arg_dual;
169 }
170 
188 template <typename T0, typename T1>
189 SERAC_HOST_DEVICE constexpr auto make_dual(const tuple<T0, T1>& args)
190 {
191  return tuple{make_dual_helper<0, 2>(get<0>(args)), make_dual_helper<1, 2>(get<1>(args))};
192 }
193 
195 template <typename T0, typename T1, typename T2>
196 SERAC_HOST_DEVICE constexpr auto make_dual(const tuple<T0, T1, T2>& args)
197 {
198  return tuple{make_dual_helper<0, 3>(get<0>(args)), make_dual_helper<1, 3>(get<1>(args)),
199  make_dual_helper<2, 3>(get<2>(args))};
200 }
201 
209 template <bool dualify, typename T>
211 {
212  if constexpr (dualify) {
213  return make_dual(x);
214  }
215  if constexpr (!dualify) {
216  return x;
217  }
218 }
219 
228 template <bool dualify, typename T, int n>
230 {
231  if constexpr (dualify) {
232  using return_type = decltype(make_dual(T{}));
233  tensor<return_type, n> output;
234  for (int i = 0; i < n; i++) {
235  output[i] = make_dual(x[i]);
236  }
237  return output;
238  }
239  if constexpr (!dualify) {
240  return x;
241  }
242 }
243 
245 template <int n, typename... T, int... i>
246 SERAC_HOST_DEVICE constexpr auto make_dual_helper(const serac::tuple<T...>& args, std::integer_sequence<int, i...>)
247 {
248  // Sam: it took me longer than I'd like to admit to find this issue, so here's an explanation
249  //
250  // note: we use serac::make_tuple(...) instead of serac::tuple{...} here because if
251  // the first argument passed in is of type `serac::tuple < serac::tuple < T ... > >`
252  // then doing something like
253  //
254  // serac::tuple{serac::get<i>(args)...};
255  //
256  // will be expand to something like
257  //
258  // serac::tuple{serac::tuple< T ... >{}};
259  //
260  // which invokes the copy ctor, returning a `serac::tuple< T ... >`
261  // instead of `serac::tuple< serac::tuple < T ... > >`
262  //
263  // but serac::make_tuple(serac::get<i>(args)...) will never accidentally trigger the copy ctor
264  return serac::make_tuple(promote_to_dual_when<i == n>(serac::get<i>(args))...);
265 }
266 
274 template <int n, typename... T>
275 constexpr auto make_dual_wrt(const serac::tuple<T...>& args)
276 {
277  return make_dual_helper<n>(args, std::make_integer_sequence<int, static_cast<int>(sizeof...(T))>{});
278 }
279 
289 template <typename T1, typename T2, int n>
291 {
292  tensor<decltype(get_value(tuple<T1, T2>{})), n> output{};
293  for (int i = 0; i < n; i++) {
294  output[i] = get_value(input[i]);
295  }
296  return output;
297 }
298 
304 template <typename... T>
305 SERAC_HOST_DEVICE auto get_value(const serac::tuple<T...>& tuple_of_values)
306 {
307  return serac::apply([](const auto&... each_value) { return serac::tuple{get_value(each_value)...}; },
308  tuple_of_values);
309 }
310 
315 template <typename... T>
317 {
318  return serac::apply([](auto... each_value) { return serac::tuple{each_value...}; }, arg.gradient);
319 }
320 
322 template <typename... T, int... n>
324 {
325  serac::tuple<outer_product_t<tensor<double, n...>, T>...> g{};
326  for_constexpr<n...>([&](auto... i) {
327  for_constexpr<sizeof...(T)>([&](auto j) { serac::get<j>(g)(i...) = serac::get<j>(arg(i...).gradient); });
328  });
329  return g;
330 }
331 
333 template <typename... T>
335 {
336  return serac::apply([](auto... each_value) { return serac::tuple{get_gradient(each_value)...}; }, tuple_of_values);
337 }
338 
344 template <int... n>
346 {
347  tensor<dual<tensor<double, n...>>, n...> A_dual{};
348  for_constexpr<n...>([&](auto... i) {
349  A_dual(i...).value = A(i...);
350  A_dual(i...).gradient(i...) = 1.0;
351  });
352  return A_dual;
353 }
354 
364 template <typename T, int n>
366 {
367  constexpr auto abs = [](double x) { return (x < 0) ? -x : x; };
368  constexpr auto swap = [](auto& x, auto& y) {
369  auto tmp = x;
370  x = y;
371  y = tmp;
372  };
373 
374  auto U = A;
375  // initialize L to Identity
376  auto L = tensor<T, n, n>{};
377  // This handles the case if T is a dual number
378  // TODO - BT: make a dense identity that is templated on type
379  for (int i = 0; i < n; i++) {
380  if constexpr (is_dual_number<T>::value) {
381  L[i][i].value = 1.0;
382  } else {
383  L[i][i] = 1.0;
384  }
385  }
386  tensor<int, n> P(make_tensor<n>([](auto i) { return i; }));
387 
388  for (int i = 0; i < n; i++) {
389  // Search for maximum in this column
390  double max_val = abs(get_value(U[i][i]));
391 
392  int max_row = i;
393  for (int j = i + 1; j < n; j++) {
394  auto U_ji = get_value(U[j][i]);
395  if (abs(U_ji) > max_val) {
396  max_val = abs(U_ji);
397  max_row = j;
398  }
399  }
400 
401  swap(P[max_row], P[i]);
402  swap(U[max_row], U[i]);
403  }
404 
405  for (int i = 0; i < n; i++) {
406  // zero entries below in this column in U
407  // and fill in L entries
408  for (int j = i + 1; j < n; j++) {
409  auto c = U[j][i] / U[i][i];
410  L[j][i] = c;
411  U[j] -= c * U[i];
412  U[j][i] = T{};
413  }
414  }
415 
416  return {P, L, U};
417 }
418 
425 template <typename S, typename T, int n, int... m>
427 {
428  // We want to avoid accumulating the derivative through the
429  // LU factorization, because it is computationally expensive.
430  // Instead, we perform the LU factorization on the values of
431  // A, and then two backsolves: one to compute the primal (x),
432  // and another to compute its derivative (dx).
433  // If A is not dual, the second solve is a no-op.
434 
435  // Strip off derivatives, if any, and compute only x (ie no derivative)
436  auto lu_factors = factorize_lu(get_value(A));
437  auto x = linear_solve(lu_factors, get_value(b));
438 
439  // Compute directional derivative of x.
440  // If both b and A are not dual, the zero type
441  // makes these no-ops.
442  auto r = get_gradient(b) - dot(get_gradient(A), x);
443  auto dx = linear_solve(lu_factors, r);
444 
445  if constexpr (is_zero<decltype(dx)>{}) {
446  return x;
447  } else {
448  return make_dual(x, dx);
449  }
450 }
451 
455 template <typename T, int n>
456 SERAC_HOST_DEVICE constexpr auto make_dual(const tensor<T, n>& x, const tensor<T, n>& dx)
457 {
458  return make_tensor<n>([&](int i) { return dual<T>{x[i], dx[i]}; });
459 }
460 
464 template <typename T, int m, int n>
465 SERAC_HOST_DEVICE constexpr auto make_dual(const tensor<T, m, n>& x, const tensor<T, m, n>& dx)
466 {
467  return make_tensor<m, n>([&](int i, int j) { return dual<T>{x[i][j], dx[i][j]}; });
468 }
469 
479 template <typename gradient_type, int n>
481 {
482  auto invA = inv(get_value(A));
483  return make_tensor<n, n>([&](int i, int j) {
484  auto value = invA[i][j];
485  gradient_type gradient{};
486  for (int k = 0; k < n; k++) {
487  for (int l = 0; l < n; l++) {
488  gradient -= invA[i][k] * A[k][l].gradient * invA[l][j];
489  }
490  }
491  return dual<gradient_type>{value, gradient};
492  });
493 }
494 
499 template <typename T, int... n>
500 SERAC_HOST_DEVICE auto get_value(const tensor<dual<T>, n...>& arg)
501 {
502  tensor<double, n...> value{};
503  for_constexpr<n...>([&](auto... i) { value(i...) = arg(i...).value; });
504  return value;
505 }
506 
511 template <int... n>
512 SERAC_HOST_DEVICE constexpr auto get_gradient(const tensor<dual<double>, n...>& arg)
513 {
514  tensor<double, n...> g{};
515  for_constexpr<n...>([&](auto... i) { g(i...) = arg(i...).gradient; });
516  return g;
517 }
518 
520 template <int... n, int... m>
522 {
523  tensor<double, n..., m...> g{};
524  for_constexpr<n...>([&](auto... i) { g(i...) = arg(i...).gradient; });
525  return g;
526 }
527 
531 struct SolverStatus {
532  bool converged;
533  unsigned int iterations;
534  double residual;
535 };
536 
541  double xtol;
542  double rtol;
543  unsigned int max_iter;
544 };
545 
547 const ScalarSolverOptions default_solver_options{.xtol = 1e-8, .rtol = 0, .max_iter = 25};
548 
576 template <typename function, typename... ParamTypes>
577 auto solve_scalar_equation(const function& f, double x0, double lower_bound, double upper_bound,
578  ScalarSolverOptions options, ParamTypes... params)
579 {
580  double x, df_dx;
581  double fl = f(lower_bound, get_value(params)...);
582  double fh = f(upper_bound, get_value(params)...);
583 
584  SLIC_ERROR_ROOT_IF(fl * fh > 0, "solve_scalar_equation: root not bracketed by input bounds.");
585 
586  unsigned int iterations = 0;
587  bool converged = false;
588 
589  // handle corner cases where one of the brackets is the root
590  if (fl == 0) {
591  x = lower_bound;
592  converged = true;
593  } else if (fh == 0) {
594  x = upper_bound;
595  converged = true;
596  }
597 
598  if (converged) {
599  df_dx = get_gradient(f(make_dual(x), get_value(params)...));
600 
601  } else {
602  // orient search so that f(xl) < 0
603  double xl = lower_bound;
604  double xh = upper_bound;
605  if (fl > 0) {
606  xl = upper_bound;
607  xh = lower_bound;
608  }
609 
610  // move initial guess if it is not between brackets
611  if (x0 < lower_bound || x0 > upper_bound) {
612  x0 = 0.5 * (lower_bound + upper_bound);
613  }
614 
615  x = x0;
616  double delta_x_old = std::abs(upper_bound - lower_bound);
617  double delta_x = delta_x_old;
618  auto R = f(make_dual(x), get_value(params)...);
619  auto fval = get_value(R);
620  df_dx = get_gradient(R);
621 
622  while (!converged) {
623  if (iterations == options.max_iter) {
624  SLIC_WARNING("solve_scalar_equation failed to converge in allotted iterations.");
625  break;
626  }
627 
628  // use bisection if Newton oversteps brackets or is not decreasing sufficiently
629  if ((x - xh) * df_dx - fval > 0 || (x - xl) * df_dx - fval < 0 ||
630  std::abs(2. * fval) > std::abs(delta_x_old * df_dx)) {
631  delta_x_old = delta_x;
632  delta_x = 0.5 * (xh - xl);
633  x = xl + delta_x;
634  converged = (x == xl);
635  } else { // use Newton step
636  delta_x_old = delta_x;
637  delta_x = fval / df_dx;
638  auto temp = x;
639  x -= delta_x;
640  converged = (x == temp);
641  }
642 
643  // function and jacobian evaluation
644  R = f(make_dual(x), get_value(params)...);
645  fval = get_value(R);
646  df_dx = get_gradient(R);
647 
648  // convergence check
649  converged = converged || (std::abs(delta_x) < options.xtol) || (std::abs(fval) < options.rtol);
650 
651  // maintain bracket on root
652  if (fval < 0) {
653  xl = x;
654  } else {
655  xh = x;
656  }
657 
658  ++iterations;
659  }
660  }
661 
662  // Accumulate derivatives so that the user can get derivatives
663  // with respect to parameters, subject to constraing that f(x, p) = 0 for all p
664  // Conceptually, we're doing the following:
665  // [fval, df_dp] = f(get_value(x), p)
666  // df = 0
667  // for p in params:
668  // df += inner(df_dp, dp)
669  // dx = -df / df_dx
670  constexpr bool contains_duals =
672  if constexpr (contains_duals) {
673  auto [fval, df] = f(x, params...);
674  auto dx = -df / df_dx;
675  SolverStatus status{.converged = converged, .iterations = iterations, .residual = fval};
676  return tuple{dual{x, dx}, status};
677  }
678  if constexpr (!contains_duals) {
679  auto fval = f(x, params...);
680  SolverStatus status{.converged = converged, .iterations = iterations, .residual = fval};
681  return tuple{x, status};
682  }
683 }
684 
697 template <typename function, int n>
698 auto find_root(const function& f, tensor<double, n> x0)
699 {
700  static_assert(std::is_same_v<decltype(f(x0)), tensor<double, n>>,
701  "error: f(x) must have the same number of equations as unknowns");
702 
703  double epsilon = 1.0e-8;
704  int max_iterations = 10;
705 
706  auto x = x0;
707 
708  for (int k = 0; k < max_iterations; k++) {
709  auto output = f(make_dual(x));
710  auto r = get_value(output);
711  if (norm(r) < epsilon) break;
712  auto J = get_gradient(output);
713  x -= linear_solve(J, r);
714  }
715 
716  return x;
717 };
718 
727 template <typename T, int size>
729 {
730  // put tensor values in an mfem::DenseMatrix
731  mfem::DenseMatrix matA(size, size);
732  for (int i = 0; i < size; i++) {
733  for (int j = 0; j < size; j++) {
734  if constexpr (is_dual_number<T>::value) {
735  matA(i, j) = A[i][j].value;
736  } else {
737  matA(i, j) = A[i][j];
738  }
739  }
740  }
741 
742  // compute eigendecomposition
743  mfem::DenseMatrixEigensystem eig_sys(matA);
744  eig_sys.Eval();
745 
746  serac::tensor<T, size> output;
747 
748  for (int k = 0; k < size; k++) {
749  // extract eigenvalues
750  output[k] = eig_sys.Eigenvalue(k);
751 
752  // and calculate their derivatives, when appropriate
753  if constexpr (is_dual_number<T>::value) {
754  tensor<double, size> phi = make_tensor<size>([&](int i) { return eig_sys.Eigenvector(k)[i]; });
755  auto dA = make_tensor<size, size>([&](int i, int j) { return A(i, j).gradient; });
756  output[k].gradient = dot(phi, dA, phi);
757  }
758  }
759 
760  return output;
761 }
762 
769 template <typename T>
770 int sgn(T val)
771 {
772  // Should we implement the derivative?
773  // It should be NaN when val = 0
774  return (T(0) < val) - (val < T(0));
775 }
776 
783 template <typename T>
785 {
786  auto swap = [](int& first, int& second) {
787  int tmp = first;
788  first = second;
789  second = tmp;
790  };
791  tensor<int, 3> order{0, 1, 2};
792  if (v[0] > v[1]) swap(order[0], order[1]);
793  if (v[order[1]] > v[order[2]]) swap(order[1], order[2]);
794  if (v[order[0]] > v[order[1]]) swap(order[0], order[1]);
795  return order;
796 }
797 
809 {
810  // We know of optimizations for this routine. When this becomes the
811  // bottleneck, we can revisit. See OptimiSM for details.
812 
813  tensor<double, 3> eta{};
814  tensor<double, 3, 3> Q = DenseIdentity<3>();
815 
816  auto A_dev = dev(A);
817  double J2 = 0.5 * inner(A_dev, A_dev);
818  double J3 = det(A_dev);
819 
820  if (J2 > 0.0) {
821  // angle used to find eigenvalues
822  double tmp = (0.5 * J3) * std::pow(3.0 / J2, 1.5);
823  double alpha = std::acos(fmin(fmax(tmp, -1.0), 1.0)) / 3.0;
824 
825  // consider the most distinct eigenvalue first
826  if (6.0 * alpha < M_PI) {
827  eta[0] = 2 * std::sqrt(J2 / 3.0) * std::cos(alpha);
828  } else {
829  eta[0] = 2 * std::sqrt(J2 / 3.0) * std::cos(alpha + 2.0 * M_PI / 3.0);
830  }
831 
832  // find the eigenvector for that eigenvalue
833  mat3 r;
834 
835  int imax = -1;
836  double norm_max = -1.0;
837 
838  for (int i = 0; i < 3; i++) {
839  for (int j = 0; j < 3; j++) {
840  r[i][j] = A_dev(j, i) - (i == j) * eta(0);
841  }
842 
843  double norm_r = norm(r[i]);
844  if (norm_max < norm_r) {
845  imax = i;
846  norm_max = norm_r;
847  }
848  }
849 
850  vec3 s0, s1, t1, t2, v0, v1, v2, w;
851 
852  s0 = normalize(r[imax]);
853  t1 = r[(imax + 1) % 3] - dot(r[(imax + 1) % 3], s0) * s0;
854  t2 = r[(imax + 2) % 3] - dot(r[(imax + 2) % 3], s0) * s0;
855  s1 = normalize((norm(t1) > norm(t2)) ? t1 : t2);
856 
857  // record the first eigenvector
858  v0 = cross(s0, s1);
859  for (int i = 0; i < 3; i++) {
860  Q[i][0] = v0[i];
861  }
862 
863  // get the other two eigenvalues by solving the
864  // remaining quadratic characteristic polynomial
865  auto A_dev_s0 = dot(A_dev, s0);
866  auto A_dev_s1 = dot(A_dev, s1);
867 
868  double A11 = dot(s0, A_dev_s0);
869  double A12 = dot(s0, A_dev_s1);
870  double A21 = A12;
871  double A22 = dot(s1, A_dev_s1);
872 
873  double delta = 0.5 * std::sqrt((A11 - A22) * (A11 - A22) + 4 * A12 * A21);
874 
875  eta(1) = 0.5 * (A11 + A22) - delta;
876  eta(2) = 0.5 * (A11 + A22) + delta;
877 
878  // if the remaining eigenvalues are exactly the same
879  // then just use the basis for the orthogonal complement
880  // found earlier
881  if (fabs(delta) <= 1.0e-15) {
882  for (int i = 0; i < 3; i++) {
883  Q[i][1] = s0(i);
884  Q[i][2] = s1(i);
885  }
886 
887  // otherwise compute the remaining eigenvectors
888  } else {
889  t1 = A_dev_s0 - eta(1) * s0;
890  t2 = A_dev_s1 - eta(1) * s1;
891 
892  w = normalize((norm(t1) > norm(t2)) ? t1 : t2);
893 
894  v1 = normalize(cross(w, v0));
895  for (int i = 0; i < 3; i++) Q[i][1] = v1(i);
896 
897  // define the last eigenvector as
898  // the direction perpendicular to the
899  // first two directions
900  v2 = normalize(cross(v0, v1));
901  for (int i = 0; i < 3; i++) Q[i][2] = v2(i);
902  }
903  }
904  // eta are actually eigenvalues of A_dev, so
905  // shift them to get eigenvalues of A
906  for (int i = 0; i < 3; i++) eta[i] += tr(A) / 3.0;
907 
908  // sort eigenvalues into ascending order
909  auto order = argsort(eta);
910  vec3 eigvals{{eta[order[0]], eta[order[1]], eta[order[2]]}};
911  // clang-format off
912  mat3 eigvecs{{{Q[0][order[0]], Q[0][order[1]], Q[0][order[2]]},
913  {Q[1][order[0]], Q[1][order[1]], Q[1][order[2]]},
914  {Q[2][order[0]], Q[2][order[1]], Q[2][order[2]]}}};
915  // clang-format on
916 
917  return {eigvals, eigvecs};
918 }
919 
920 /*
921 // Should we provide this fallback, or force the author to consider how to
922 // write a numerically stable version on a case-by-case basis?
923 // The convenience of this is somewhat undermined by the fact that it would
924 // only work for functions that already have a dual number overload.
925 template <typename Function>
926 double generic_eigenvalue_tangent(double lam1, double lam2, const Function& f)
927 {
928  if (lam1 == lam2) {
929  return f(make_dual(lam1));
930  } else {
931  return (f(lam1) - f(lam2))/(lam1 - lam2);
932  }
933 }
934 */
935 
968 template <typename T, typename Function, typename EigvalSecantFunction>
969 auto symmetric_mat3_function(tensor<T, 3, 3> A, const Function& f, const EigvalSecantFunction& g)
970 {
971  auto [lambda, Q] = eig_symm(get_value(A));
972  vec3 y;
973  for (int i = 0; i < 3; i++) {
974  y[i] = f(lambda[i]);
975  }
976  auto f_A = dot(Q, dot(diag(y), transpose(Q)));
977 
978  if constexpr (!is_dual_number<T>::value) {
979  return f_A;
980  } else {
981  return symmetric_mat3_function_with_derivative(A, f_A, lambda, Q, g);
982  }
983 }
984 
986 template <typename Gradient, typename Function>
988  tensor<double, 3, 3> f_A, vec3 lambda, mat3 Q,
989  const Function& g)
990 {
991  return make_tensor<3, 3>([&](int i, int j) {
992  auto value = f_A[i][j];
993  Gradient gradient{};
994  for (int k = 0; k < 3; k++) {
995  for (int l = 0; l < 3; l++) {
996  for (int a = 0; a < 3; a++) {
997  for (int b = 0; b < 3; b++) {
998  gradient += g(lambda[a], lambda[b]) * Q[k][a] * Q[l][b] * Q[i][a] * Q[j][b] * A[k][l].gradient;
999  }
1000  }
1001  }
1002  }
1003  return dual<Gradient>{value, gradient};
1004  });
1005 }
1006 
1013 template <typename T>
1015 {
1016  auto g = [](double lam1, double lam2) {
1017  if (lam1 == lam2) {
1018  return 1 / lam1;
1019  } else {
1020  double y = lam1 / lam2;
1021  return (std::log(y) / (y - 1.0)) / lam2;
1022  }
1023  };
1024  return symmetric_mat3_function(
1025  A, [](double x) { return std::log(x); }, g);
1026 }
1027 
1034 template <typename T>
1036 {
1037  auto g = [](double lam1, double lam2) {
1038  if (lam1 == lam2) {
1039  return std::exp(lam1);
1040  } else {
1041  double arg = lam1 - lam2;
1042  return std::exp(lam2) * std::expm1(arg) / arg;
1043  }
1044  };
1045  return symmetric_mat3_function(
1046  A, [](double x) { return std::exp(x); }, g);
1047 }
1048 
1055 template <typename T>
1057 {
1058  auto g = [](double lam1, double lam2) { return 1.0 / (std::sqrt(lam1) + std::sqrt(lam2)); };
1059  return symmetric_mat3_function(
1060  A, [](double x) { return std::sqrt(x); }, g);
1061 }
1062 
1063 } // namespace serac
#define SERAC_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
Definition: accelerator.hpp:38
This file contains the declaration of a dual number class.
constexpr SERAC_HOST_DEVICE void for_constexpr(const lambda &f)
multidimensional loop tool that evaluates the lambda body inside the innermost loop.
Accelerator functionality.
Definition: serac.cpp:36
auto solve_scalar_equation(const function &f, double x0, double lower_bound, double upper_bound, ScalarSolverOptions options, ParamTypes... params)
Solves a nonlinear scalar-valued equation and gives derivatives of solution to parameters.
constexpr SERAC_HOST_DEVICE auto inner(const dual< S > &A, const dual< T > &B)
Definition: dual.hpp:281
constexpr SERAC_HOST_DEVICE auto make_dual(double x)
promote a value to a dual number of the appropriate type
Definition: dual.hpp:441
SERAC_HOST_DEVICE tuple< vec3, mat3 > eig_symm(const mat3 &A)
SERAC_HOST_DEVICE tensor< int, 3 > argsort(const tensor< T, 3 > &v)
Find indices that would sort a 3-vector.
typename detail::outer_prod< T1, T2 >::type outer_product_t
a type function that returns the tensor type of an outer product of two tensors
Definition: tensor.hpp:1824
constexpr SERAC_HOST_DEVICE auto tr(const isotropic_tensor< T, m, m > &I)
calculate the trace of an isotropic tensor
constexpr SERAC_HOST_DEVICE auto dot(const isotropic_tensor< S, m, m > &I, const tensor< T, m, n... > &A)
dot product between an isotropic and (nonisotropic) tensor
SERAC_HOST_DEVICE auto apply(lambda f, tuple< T... > &args)
a way of passing an n-tuple to a function that expects n separate arguments
Definition: tuple.hpp:782
auto find_root(const function &f, tensor< double, n > x0)
Finds a root of a vector-valued nonlinear function.
constexpr SERAC_HOST_DEVICE auto operator*(const dual< gradient_type > &a, double b)
multiplication of a dual number and a non-dual number
Definition: dual.hpp:108
SERAC_HOST_DEVICE auto cos(dual< gradient_type > a)
implementation of cosine for dual numbers
Definition: dual.hpp:316
auto exp_symm(tensor< T, 3, 3 > A)
Exponential of a symmetric matrix.
SERAC_HOST_DEVICE auto promote_each_to_dual_when(const tensor< T, n > &x)
a function that optionally (decided at compile time) converts a list of values to their dual types
auto log_symm(tensor< T, 3, 3 > A)
Logarithm of a symmetric matrix.
constexpr SERAC_HOST_DEVICE auto make_dual_helper(zero)
This is an overloaded member function, provided for convenience. It differs from the above function o...
SERAC_HOST_DEVICE auto pow(dual< gradient_type > a, dual< gradient_type > b)
implementation of a (dual) raised to the b (dual) power
Definition: dual.hpp:405
constexpr SERAC_HOST_DEVICE tensor< T, n, n > diag(const tensor< T, n > &d)
Returns a square diagonal matrix by specifying the diagonal entries.
Definition: tensor.hpp:1226
constexpr auto make_dual_wrt(const serac::tuple< T... > &args)
take a tuple of values, and promote the nth one to a one-hot dual number of the appropriate type
SERAC_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
Definition: dual.hpp:308
constexpr SERAC_HOST_DEVICE auto symmetric_mat3_function_with_derivative(tensor< dual< Gradient >, 3, 3 > A, tensor< double, 3, 3 > f_A, vec3 lambda, mat3 Q, const Function &g)
Helper function for defining the derivative.
SERAC_HOST_DEVICE auto promote_to_dual_when(const T &x)
a function that optionally (decided at compile time) converts a value to its dual type
constexpr SERAC_HOST_DEVICE int size(const tensor< T, n... > &)
returns the total number of stored values in a tensor
Definition: tensor.hpp:1934
constexpr SERAC_HOST_DEVICE auto linear_solve(const LuFactorization< S, n > &lu_factors, const tensor< T, n, m... > &b)
Definition: tensor.hpp:1621
constexpr SERAC_HOST_DEVICE auto dev(const tensor< T, n, n > &A)
Calculates the deviator of a matrix (rank-2 tensor)
Definition: tensor.hpp:1195
SERAC_HOST_DEVICE tuple< T... > make_tuple(const T &... args)
helper function for combining a list of values into a tuple
Definition: tuple.hpp:267
int sgn(T val)
Signum, returns sign of input.
constexpr SERAC_HOST_DEVICE auto get_value(const T &arg)
return the "value" part from a given type. For non-dual types, this is just the identity function
Definition: dual.hpp:445
SERAC_HOST_DEVICE auto abs(dual< gradient_type > x)
Implementation of absolute value function for dual numbers.
Definition: dual.hpp:219
constexpr SERAC_HOST_DEVICE auto det(const isotropic_tensor< T, m, m > &I)
compute the determinant of an isotropic tensor
const ScalarSolverOptions default_solver_options
Default options for solve_scalar_equation.
SERAC_HOST_DEVICE auto acos(dual< gradient_type > a)
implementation of acos for dual numbers
Definition: dual.hpp:373
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
Definition: tuple.hpp:379
constexpr SERAC_HOST_DEVICE auto norm(const isotropic_tensor< T, m, m > &I)
compute the Frobenius norm (sqrt(tr(dot(transpose(I), I)))) of an isotropic tensor
auto cross(const tensor< T, 3, 2 > &A)
compute the cross product of the columns of A: A(:,1) x A(:,2)
Definition: tensor.hpp:961
constexpr SERAC_HOST_DEVICE auto inv(const isotropic_tensor< T, m, m > &I)
return the inverse of an isotropic tensor
constexpr SERAC_HOST_DEVICE auto transpose(const isotropic_tensor< T, m, m > &I)
return the transpose of an isotropic tensor
constexpr SERAC_HOST_DEVICE LuFactorization< T, n > factorize_lu(const tensor< T, n, n > &A)
Compute LU factorization of a matrix with partial pivoting.
auto eigenvalues(const serac::tensor< T, size, size > &A)
compute the eigenvalues of a symmetric matrix A
SERAC_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
Definition: tensor.hpp:1117
auto symmetric_mat3_function(tensor< T, 3, 3 > A, const Function &f, const EigvalSecantFunction &g)
Constructs an isotropic tensor-valued function of a symmetric 3x3 tensor from a scalar function.
SERAC_HOST_DEVICE auto log(dual< gradient_type > a)
implementation of the natural logarithm function for dual numbers
Definition: dual.hpp:389
typename one_hot< i, n, T >::type one_hot_t
a tuple type with n entries, all of which are of type serac::zero, except for the i^{th} entry,...
constexpr SERAC_HOST_DEVICE auto get_gradient(dual< gradient_type > arg)
return the "gradient" part from a dual number type
Definition: dual.hpp:459
auto sqrt_symm(tensor< T, 3, 3 > A)
Square root of a symmetric matrix.
constexpr SERAC_HOST_DEVICE auto operator/(const dual< gradient_type > &a, double b)
division of a dual number by a non-dual number
Definition: dual.hpp:129
SERAC_HOST_DEVICE auto exp(dual< gradient_type > a)
implementation of exponential function for dual numbers
Definition: dual.hpp:381
Representation of an LU factorization.
Definition: tensor.hpp:1542
Settings for solve_scalar_equation.
unsigned int max_iter
maximum allowed number of iterations
double xtol
absolute tolerance on Newton correction
double rtol
absolute tolerance on absolute value of residual
Status and diagnostics of nonlinear equation solvers.
bool converged
converged Flag indicating whether solver converged to a solution or aborted.
double residual
Final value of residual.
unsigned int iterations
Number of iterations taken.
Dual number struct (value plus gradient)
Definition: dual.hpp:28
double value
the actual numerical value
Definition: dual.hpp:29
class for checking if a type is a dual number or not
Definition: dual.hpp:466
class for checking if a type is a tensor of dual numbers or not
static constexpr bool value
whether or not type T is a dual number
checks if a type is zero
Definition: tensor.hpp:150
Arbitrary-rank tensor class.
Definition: tensor.hpp:28
Type that mimics std::tuple.
Definition: tuple.hpp:61
Type that mimics std::tuple.
Definition: tuple.hpp:48
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
Definition: tuple.hpp:28
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122
Implementation of the tensor class used by Functional.
Implements a std::tuple-like object that works in CUDA kernels.