20 static constexpr
bool value =
false;
24 template <
typename T,
int... n>
26 static constexpr
bool value =
true;
37 template <
typename S,
typename T,
int m,
int... n,
41 tensor<decltype(S{} * T{}), m, n...> C{};
42 for (
int i = 0; i < m; i++) {
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>>
60 tensor<decltype(T{} * S{}), m, n...> C{};
61 for (
int i = 0; i < m; i++) {
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>>
79 tensor<decltype(S{} * T{}), n...> C{};
80 for (
int i = 0; i < m; i++) {
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>>
98 tensor<decltype(T{} * S{}), m, n...> C{};
99 for (
int i = 0; i < m; i++) {
106 template <
int i,
typename S,
typename T>
107 struct one_hot_helper;
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>...>;
114 template <
int i,
int n,
typename T>
115 struct one_hot :
public one_hot_helper<i, std::make_integer_sequence<int, n>, T> {
125 template <
int i,
int n,
typename T>
129 template <
int i,
int N>
142 template <
int i,
int N>
147 arg_dual.
value = arg;
148 serac::get<i>(arg_dual.gradient) = 1.0;
159 template <
int i,
int N,
typename T,
int... n>
165 arg_dual(j...).value = arg(j...);
166 serac::get<i>(arg_dual(j...).gradient)(j...) = 1.0;
188 template <
typename T0,
typename T1>
191 return tuple{make_dual_helper<0, 2>(get<0>(args)), make_dual_helper<1, 2>(get<1>(args))};
195 template <
typename T0,
typename T1,
typename T2>
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))};
209 template <
bool dualify,
typename T>
212 if constexpr (dualify) {
215 if constexpr (!dualify) {
228 template <
bool dualify,
typename T,
int n>
231 if constexpr (dualify) {
232 using return_type = decltype(
make_dual(T{}));
234 for (
int i = 0; i < n; i++) {
239 if constexpr (!dualify) {
245 template <
int n,
typename... T,
int... i>
274 template <
int n,
typename... T>
277 return make_dual_helper<n>(args, std::make_integer_sequence<
int,
static_cast<int>(
sizeof...(T))>{});
289 template <
typename T1,
typename T2,
int n>
293 for (
int i = 0; i < n; i++) {
304 template <
typename... T>
315 template <
typename... T>
322 template <
typename... T,
int... n>
327 for_constexpr<
sizeof...(T)>([&](
auto j) { serac::get<j>(g)(i...) = serac::get<j>(arg(i...).gradient); });
333 template <
typename... T>
349 A_dual(i...).value = A(i...);
350 A_dual(i...).gradient(i...) = 1.0;
364 template <
typename T,
int n>
367 constexpr
auto abs = [](
double x) {
return (x < 0) ? -x : x; };
368 constexpr
auto swap = [](
auto& x,
auto& y) {
379 for (
int i = 0; i < n; i++) {
388 for (
int i = 0; i < n; i++) {
393 for (
int j = i + 1; j < n; j++) {
395 if (
abs(U_ji) > max_val) {
401 swap(P[max_row], P[i]);
402 swap(U[max_row], U[i]);
405 for (
int i = 0; i < n; i++) {
408 for (
int j = i + 1; j < n; j++) {
409 auto c = U[j][i] / U[i][i];
425 template <
typename S,
typename T,
int n,
int... m>
445 if constexpr (
is_zero<decltype(dx)>{}) {
455 template <
typename T,
int n>
458 return make_tensor<n>([&](
int i) {
return dual<T>{x[i], dx[i]}; });
464 template <
typename T,
int m,
int n>
467 return make_tensor<m, n>([&](
int i,
int j) {
return dual<T>{x[i][j], dx[i][j]}; });
479 template <
typename gradient_type,
int n>
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];
499 template <
typename T,
int... n>
502 tensor<double, n...> value{};
503 for_constexpr<n...>([&](
auto... i) { value(i...) = arg(i...).value; });
515 for_constexpr<n...>([&](
auto... i) { g(i...) = arg(i...).gradient; });
520 template <
int... n,
int... m>
523 tensor<double, n..., m...> g{};
524 for_constexpr<n...>([&](
auto... i) { g(i...) = arg(i...).gradient; });
576 template <
typename function,
typename... ParamTypes>
581 double fl = f(lower_bound,
get_value(params)...);
582 double fh = f(upper_bound,
get_value(params)...);
584 SLIC_ERROR_ROOT_IF(fl * fh > 0,
"solve_scalar_equation: root not bracketed by input bounds.");
586 unsigned int iterations = 0;
587 bool converged =
false;
593 }
else if (fh == 0) {
603 double xl = lower_bound;
604 double xh = upper_bound;
611 if (x0 < lower_bound || x0 > upper_bound) {
612 x0 = 0.5 * (lower_bound + upper_bound);
616 double delta_x_old =
std::abs(upper_bound - lower_bound);
617 double delta_x = delta_x_old;
623 if (iterations == options.
max_iter) {
624 SLIC_WARNING(
"solve_scalar_equation failed to converge in allotted iterations.");
629 if ((x - xh) * df_dx - fval > 0 || (x - xl) * df_dx - fval < 0 ||
631 delta_x_old = delta_x;
632 delta_x = 0.5 * (xh - xl);
634 converged = (x == xl);
636 delta_x_old = delta_x;
637 delta_x = fval / df_dx;
640 converged = (x == temp);
670 constexpr
bool contains_duals =
672 if constexpr (contains_duals) {
673 auto [fval, df] = f(x, params...);
674 auto dx = -df / df_dx;
678 if constexpr (!contains_duals) {
679 auto fval = f(x, params...);
681 return tuple{x, status};
697 template <
typename function,
int n>
701 "error: f(x) must have the same number of equations as unknowns");
703 double epsilon = 1.0e-8;
704 int max_iterations = 10;
708 for (
int k = 0; k < max_iterations; k++) {
711 if (
norm(r) < epsilon)
break;
727 template <
typename T,
int size>
732 for (
int i = 0; i <
size; i++) {
733 for (
int j = 0; j <
size; j++) {
735 matA(i, j) = A[i][j].value;
737 matA(i, j) = A[i][j];
743 mfem::DenseMatrixEigensystem eig_sys(matA);
748 for (
int k = 0; k <
size; k++) {
750 output[k] = eig_sys.Eigenvalue(k);
755 auto dA = make_tensor<size, size>([&](
int i,
int j) {
return A(i, j).gradient; });
756 output[k].gradient =
dot(phi, dA, phi);
769 template <
typename T>
774 return (T(0) < val) - (val < T(0));
783 template <
typename T>
786 auto swap = [](
int& first,
int& second) {
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]);
817 double J2 = 0.5 *
inner(A_dev, A_dev);
818 double J3 =
det(A_dev);
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;
826 if (6.0 * alpha < M_PI) {
836 double norm_max = -1.0;
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);
843 double norm_r =
norm(r[i]);
844 if (norm_max < norm_r) {
850 vec3 s0, s1, t1, t2, v0, v1, v2, w;
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;
859 for (
int i = 0; i < 3; i++) {
865 auto A_dev_s0 =
dot(A_dev, s0);
866 auto A_dev_s1 =
dot(A_dev, s1);
868 double A11 =
dot(s0, A_dev_s0);
869 double A12 =
dot(s0, A_dev_s1);
871 double A22 =
dot(s1, A_dev_s1);
873 double delta = 0.5 *
std::sqrt((A11 - A22) * (A11 - A22) + 4 * A12 * A21);
875 eta(1) = 0.5 * (A11 + A22) - delta;
876 eta(2) = 0.5 * (A11 + A22) + delta;
881 if (fabs(delta) <= 1.0e-15) {
882 for (
int i = 0; i < 3; i++) {
889 t1 = A_dev_s0 - eta(1) * s0;
890 t2 = A_dev_s1 - eta(1) * s1;
895 for (
int i = 0; i < 3; i++) Q[i][1] = v1(i);
901 for (
int i = 0; i < 3; i++) Q[i][2] = v2(i);
906 for (
int i = 0; i < 3; i++) eta[i] +=
tr(A) / 3.0;
910 vec3 eigvals{{eta[order[0]], eta[order[1]], eta[order[2]]}};
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]]}}};
917 return {eigvals, eigvecs};
968 template <
typename T,
typename Function,
typename EigvalSecantFunction>
973 for (
int i = 0; i < 3; i++) {
986 template <
typename Gradient,
typename Function>
991 return make_tensor<3, 3>([&](
int i,
int j) {
992 auto value = f_A[i][j];
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;
1013 template <
typename T>
1016 auto g = [](
double lam1,
double lam2) {
1020 double y = lam1 / lam2;
1021 return (
std::log(y) / (y - 1.0)) / lam2;
1025 A, [](
double x) {
return std::log(x); }, g);
1034 template <
typename T>
1037 auto g = [](
double lam1,
double lam2) {
1041 double arg = lam1 - lam2;
1042 return std::exp(lam2) * std::expm1(arg) / arg;
1046 A, [](
double x) {
return std::exp(x); }, g);
1055 template <
typename T>
1058 auto g = [](
double lam1,
double lam2) {
return 1.0 / (
std::sqrt(lam1) +
std::sqrt(lam2)); };
1060 A, [](
double x) {
return std::sqrt(x); }, g);
#define SERAC_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
This file contains the declaration of a dual number class.
Accelerator functionality.
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)
constexpr SERAC_HOST_DEVICE auto make_dual(double x)
promote a value to a dual number of the appropriate type
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
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
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
SERAC_HOST_DEVICE auto cos(dual< gradient_type > a)
implementation of cosine for dual numbers
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
constexpr SERAC_HOST_DEVICE tensor< T, n, n > diag(const tensor< T, n > &d)
Returns a square diagonal matrix by specifying the diagonal entries.
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
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
constexpr SERAC_HOST_DEVICE auto linear_solve(const LuFactorization< S, n > &lu_factors, const tensor< T, n, m... > &b)
constexpr SERAC_HOST_DEVICE auto dev(const tensor< T, n, n > &A)
Calculates the deviator of a matrix (rank-2 tensor)
SERAC_HOST_DEVICE tuple< T... > make_tuple(const T &... args)
helper function for combining a list of values into a tuple
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
SERAC_HOST_DEVICE auto abs(dual< gradient_type > x)
Implementation of absolute value function for dual numbers.
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
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
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)
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,...
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
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
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
SERAC_HOST_DEVICE auto exp(dual< gradient_type > a)
implementation of exponential function for dual numbers
Representation of an LU factorization.
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)
double value
the actual numerical value
class for checking if a type is a dual number or not
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
Arbitrary-rank tensor class.
Type that mimics std::tuple.
Type that mimics std::tuple.
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
A sentinel struct for eliding no-op tensor operations.
Implementation of the tensor class used by Functional.
Implements a std::tuple-like object that works in CUDA kernels.