27 template <
typename T,
int... n>
31 template <
typename T,
int m,
int... n>
32 struct tensor<T, m, n...> {
33 template <
typename i_type>
38 template <
typename i_type>
43 template <
typename i_type,
typename... jklm_type>
46 return data[i](jklm...);
48 template <
typename i_type,
typename... jklm_type>
51 return data[i](jklm...);
55 SMITH_HOST_DEVICE constexpr
const auto& operator[](
int i)
const {
return data[i]; }
60 template <
typename T,
int m>
62 template <
typename i_type>
67 template <
typename i_type>
73 SMITH_HOST_DEVICE constexpr
const auto& operator[](
int i)
const {
return data[i]; }
99 template <
typename T,
int n1>
110 template <
typename T,
int n1,
int n2>
127 template <
typename T,
int... n>
134 template <
typename... T>
141 template <
typename T>
149 template <
typename T>
160 template <
typename T>
167 template <
typename T>
182 template <
typename T>
189 template <
typename T>
201 template <
typename T>
208 template <
typename T>
215 template <
typename T>
222 template <
typename T>
225 static_assert(::detail::always_false<T>{},
"Error: Can't divide by zero!");
249 template <
typename T>
256 template <
typename T>
263 template <
typename T,
int m,
int... n>
276 template <
typename T,
int n1,
int n2 = 1>
278 (n1 == 1 && n2 == 1),
double,
286 template <
typename T,
int... n>
303 template <
typename lambda_type>
306 using T = decltype(f());
322 template <
int n1,
typename lambda_type>
325 using T = decltype(f(n1));
327 for (
int i = 0; i < n1; i++) {
346 template <
int n1,
int n2,
typename lambda_type>
349 using T = decltype(f(n1, n2));
351 for (
int i = 0; i < n1; i++) {
352 for (
int j = 0; j < n2; j++) {
373 template <
int n1,
int n2,
int n3,
typename lambda_type>
376 using T = decltype(f(n1, n2, n3));
378 for (
int i = 0; i < n1; i++) {
379 for (
int j = 0; j < n2; j++) {
380 for (
int k = 0; k < n3; k++) {
381 A(i, j, k) = f(i, j, k);
403 template <
int n1,
int n2,
int n3,
int n4,
typename lambda_type>
406 using T = decltype(f(n1, n2, n3, n4));
408 for (
int i = 0; i < n1; i++) {
409 for (
int j = 0; j < n2; j++) {
410 for (
int k = 0; k < n3; k++) {
411 for (
int l = 0; l < n4; l++) {
412 A(i, j, k, l) = f(i, j, k, l);
428 template <
typename S,
typename T,
int m,
int... n>
431 tensor<decltype(S{} + T{}), m, n...> C{};
432 for (
int i = 0; i < m; i++) {
444 template <
typename T,
int m,
int... n>
448 for (
int i = 0; i < m; i++) {
462 template <
typename S,
typename T,
int m,
int... n>
465 tensor<decltype(S{} + T{}), m, n...> C{};
466 for (
int i = 0; i < m; i++) {
480 template <
typename S,
typename T,
int m,
int... n>
483 for (
int i = 0; i < m; i++) {
496 template <
typename T>
509 template <
typename T,
int n>
512 for (
int i = 0; i < n; i++) {
513 A.data[i][0] += B[i];
524 template <
typename T,
int n>
527 for (
int i = 0; i < n; i++) {
528 A.data[0][i] += B[i];
539 template <
typename T>
542 return A.data[0] += B;
551 template <
typename T>
554 return A.data[0][0] += B;
563 template <
typename T,
int... n>
577 template <
typename S,
typename T,
int m,
int... n>
580 for (
int i = 0; i < m; i++) {
592 template <
typename T,
int... n>
602 template <
typename T,
int n>
605 tensor<decltype(
double{} * T{}), n> AB{};
606 for (
int i = 0; i < n; i++) {
616 template <
typename T,
int m>
619 tensor<decltype(T{} *
double{}), m> AB{};
620 for (
int i = 0; i < m; i++) {
630 template <
typename T,
int n>
640 template <
typename T,
int n>
650 template <
typename S,
typename T,
int m,
int n>
653 tensor<decltype(S{} * T{}), m, n> AB{};
654 for (
int i = 0; i < m; i++) {
655 for (
int j = 0; j < n; j++) {
656 AB[i][j] = A[i] * B[j];
671 template <
typename S,
typename T,
int m,
int n>
674 decltype(S{} * T{}) sum{};
675 for (
int i = 0; i < m; i++) {
676 for (
int j = 0; j < n; j++) {
677 sum += A[i][j] * B[i][j];
687 template <
typename S,
typename T,
int m>
690 decltype(S{} * T{}) sum{};
691 for (
int i = 0; i < m; i++) {
709 template <
typename S,
int m,
int n>
719 template <
typename S,
int m>
737 template <
typename T,
int m,
int n>
747 template <
typename T,
int m>
767 template <
typename S,
typename T,
int m,
int n,
int p>
770 tensor<decltype(S{} * T{}), m, p> AB{};
771 for (
int i = 0; i < m; i++) {
772 for (
int j = 0; j < p; j++) {
773 for (
int k = 0; k < n; k++) {
774 AB[i][j] = AB[i][j] + A[i][k] * B[k][j];
785 template <
typename T,
int m>
795 template <
typename T,
int m>
805 template <
typename S,
typename T,
int m>
808 decltype(S{} * T{}) AB{};
809 for (
int i = 0; i < m; i++) {
810 AB = AB + A[i] * B[i];
819 template <
typename S,
typename T,
int m,
int n>
822 tensor<decltype(S{} * T{}), n> AB{};
823 for (
int i = 0; i < n; i++) {
824 for (
int j = 0; j < m; j++) {
825 AB[i] = AB[i] + A[j] * B[j][i];
835 template <
typename S,
typename T,
int m,
int n,
int p>
838 tensor<decltype(S{} * T{}), n, p> AB{};
839 for (
int j = 0; j < m; j++) {
840 AB = AB + A[j] * B[j];
855 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
858 tensor<decltype(S{} * T{}), n, p, q> AB{};
859 for (
int j = 0; j < m; j++) {
860 AB = AB + A[j] * B[j];
869 template <
typename S,
typename T,
int m,
int n>
872 tensor<decltype(S{} * T{}), m> AB{};
873 for (
int i = 0; i < m; i++) {
874 for (
int j = 0; j < n; j++) {
875 AB[i] = AB[i] + A[i][j] * B[j];
885 template <
typename S,
typename T,
int m,
int n,
int p,
int q,
int r>
888 tensor<decltype(S{} * T{}), m, p, q, r> AB{};
889 for (
int i = 0; i < m; i++) {
890 for (
int j = 0; j < n; j++) {
891 AB[i] = AB[i] + A[i][j] * B[j];
901 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
904 tensor<decltype(S{} * T{}), m, p, q> AB{};
905 for (
int i = 0; i < m; i++) {
906 for (
int j = 0; j < n; j++) {
907 AB[i] = AB[i] + A[i][j] * B[j];
917 template <
typename S,
typename T,
int m,
int n,
int p>
920 tensor<decltype(S{} * T{}), m, n> AB{};
921 for (
int i = 0; i < m; i++) {
922 for (
int j = 0; j < n; j++) {
923 for (
int k = 0; k < p; k++) {
924 AB[i][j] += A[i][j][k] * B[k];
935 template <
typename S,
typename T,
typename U,
int m,
int n>
938 decltype(S{} * T{} * U{}) uAv{};
939 for (
int i = 0; i < m; i++) {
940 for (
int j = 0; j < n; j++) {
941 uAv += u[i] * A[i][j] * v[j];
948 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
951 tensor<decltype(S{} * T{}), m, n, p> AB{};
952 for (
int i = 0; i < m; i++) {
953 for (
int j = 0; j < n; j++) {
954 for (
int k = 0; k < p; k++) {
955 for (
int l = 0; l < q; l++) {
956 AB[i][j][k] += A[i][j][k][l] * B[l];
965 template <
typename T>
968 return tensor<T, 3>{A(1, 0) * A(2, 1) - A(2, 0) * A(1, 1), A(2, 0) * A(0, 1) - A(0, 0) * A(2, 1),
969 A(0, 0) * A(1, 1) - A(1, 0) * A(0, 1)};
973 template <
typename T>
980 template <
typename T>
987 template <
typename S,
typename T>
990 return tensor<decltype(S{} * T{}), 3>{u(1) * v(2) - u(2) * v(1), u(2) * v(0) - u(0) * v(2),
991 u(0) * v(1) - u(1) * v(0)};
1005 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
1008 tensor<decltype(S{} * T{}), m, n> AB{};
1009 for (
int i = 0; i < m; i++) {
1010 for (
int j = 0; j < n; j++) {
1011 for (
int k = 0; k < p; k++) {
1012 for (
int l = 0; l < q; l++) {
1013 AB[i][j] += A[i][j][k][l] * B[k][l];
1025 template <
typename S,
typename T,
int m,
int n,
int p>
1028 tensor<decltype(S{} * T{}), m> AB{};
1029 for (
int i = 0; i < m; i++) {
1030 for (
int j = 0; j < n; j++) {
1031 for (
int k = 0; k < p; k++) {
1032 AB[i] += A[i][j][k] * B[j][k];
1043 template <
typename S,
typename T,
int m,
int n>
1046 decltype(S{} * T{}) AB{};
1047 for (
int i = 0; i < m; i++) {
1048 for (
int j = 0; j < n; j++) {
1049 AB += A[i][j] * B[i][j];
1058 template <
typename S,
typename T,
int... m,
int... n>
1068 template <
typename T,
int m>
1072 for (
int i = 0; i < m; i++) {
1073 total += A[i] * A[i];
1079 template <
typename T,
int m,
int n>
1083 for (
int i = 0; i < m; i++) {
1084 for (
int j = 0; j < n; j++) {
1085 total += A[i][j] * A[i][j];
1092 template <
typename T,
int... n>
1096 for_constexpr<n...>([&](
auto... i) { total += A(i...) * A(i...); });
1104 template <
typename T,
int... n>
1121 template <
typename T,
int... n>
1132 template <
typename T>
1136 output[0][0] = A[0][0];
1137 output[0][1] = A[0][1];
1138 output[1][0] = A[1][0];
1139 output[1][1] = A[1][1];
1148 template <
typename T,
int n>
1152 for (
int i = 0; i < n; i++) {
1153 trA = trA + A[i][i];
1163 template <
typename T,
int n>
1167 for (
int i = 0; i < n; i++) {
1168 for (
int j = 0; j < n; j++) {
1169 symA[i][j] = 0.5 * (A[i][j] + A[j][i]);
1180 template <
typename T,
int n>
1184 for (
int i = 0; i < n; i++) {
1185 for (
int j = 0; j < n; j++) {
1186 antisymA[i][j] = 0.5 * (A[i][j] - A[j][i]);
1199 template <
typename T,
int n>
1204 for (
int i = 0; i < n; i++) {
1205 devA[i][i] -= trA / n;
1216 template <
typename T,
int n>
1220 for (
int i = 0; i < n; i++) {
1230 template <
typename T,
int n>
1234 for (
int i = 0; i < n; i++) {
1244 template <
typename T,
int n>
1248 for (
int i = 0; i < n; i++) {
1262 for (
int i = 0; i < dim; i++) {
1263 for (
int j = 0; j < dim; j++) {
1274 template <
typename T,
int m,
int n>
1278 for (
int i = 0; i < n; i++) {
1279 for (
int j = 0; j < m; j++) {
1290 template <
typename T>
1293 return +A[0][0] * A[1][1] + A[1][1] * A[2][2] + A[2][2] * A[0][0] - A[0][1] * A[1][0] - A[1][2] * A[2][1] -
1301 template <
typename T>
1304 return A[0][0] * A[1][1] - A[0][1] * A[1][0];
1308 template <
typename T>
1311 return A[0][0] * A[1][1] * A[2][2] + A[0][1] * A[1][2] * A[2][0] + A[0][2] * A[1][0] * A[2][1] -
1312 A[0][0] * A[1][2] * A[2][1] - A[0][1] * A[1][0] * A[2][2] - A[0][2] * A[1][1] * A[2][0];
1325 template <
typename T>
1334 return A(0, 0) - A(0, 1) * A(1, 0) + A(1, 1) + A(0, 0) * A(1, 1);
1338 template <
typename T>
1345 return A(0, 0) + A(1, 1) + A(2, 2)
1346 - A(0, 1) * A(1, 0) * (1 + A(2, 2))
1347 + A(0, 0) * A(1, 1) * (1 + A(2, 2))
1348 - A(0, 2) * A(2, 0) * (1 + A(1, 1))
1349 - A(1, 2) * A(2, 1) * (1 + A(0, 0))
1352 + A(0, 1) * A(1, 2) * A(2, 0)
1353 + A(0, 2) * A(1, 0) * A(2, 1);
1365 template <
typename T,
int dim>
1369 for (
int i = 0; i < 15; i++) {
1370 B = 0.5 * (B +
dot(A,
inv(B)));
1404 template <
int i1,
int i2,
typename S,
int m,
int... n,
typename T,
int p,
int q>
1407 constexpr
int Adims[] = {m, n...};
1408 constexpr
int Bdims[] = {p, q};
1409 static_assert(
sizeof...(n) < 3);
1410 static_assert(Adims[i1] == Bdims[i2],
"error: incompatible tensor dimensions");
1413 constexpr
int new_dim = (i2 == 0) ? q : p;
1414 constexpr
int d1 = (i1 == 0) ? new_dim : Adims[0];
1415 constexpr
int d2 = (i1 == 1) ? new_dim : Adims[1];
1416 constexpr
int d3 =
sizeof...(n) == 1 ? 0 : ((i1 == 2) ? new_dim : Adims[2]);
1419 using U = decltype(S{} * T{});
1426 if constexpr (d3 == 0) {
1427 for (
int i = 0; i < d1; i++) {
1428 for (
int j = 0; j < d2; j++) {
1430 for (
int k = 0; k < Adims[i1]; k++) {
1431 if constexpr (i1 == 0 && i2 == 0) sum += A(k, j) * B(k, i);
1432 if constexpr (i1 == 1 && i2 == 0) sum += A(i, k) * B(k, j);
1433 if constexpr (i1 == 0 && i2 == 1) sum += A(k, j) * B(i, k);
1434 if constexpr (i1 == 1 && i2 == 1) sum += A(i, k) * B(j, k);
1440 for (
int i = 0; i < d1; i++) {
1441 for (
int j = 0; j < d2; j++) {
1442 for (
int k = 0; k < d3; k++) {
1444 for (
int l = 0; l < Adims[i1]; l++) {
1445 if constexpr (i1 == 0 && i2 == 0) sum += A(l, j, k) * B(l, i);
1446 if constexpr (i1 == 1 && i2 == 0) sum += A(i, l, k) * B(l, j);
1447 if constexpr (i1 == 2 && i2 == 0) sum += A(i, j, l) * B(l, k);
1448 if constexpr (i1 == 0 && i2 == 1) sum += A(l, j, k) * B(i, l);
1449 if constexpr (i1 == 1 && i2 == 1) sum += A(i, l, k) * B(j, l);
1450 if constexpr (i1 == 2 && i2 == 1) sum += A(i, j, l) * B(k, l);
1462 template <
int i1,
int i2,
typename T>
1477 template <
typename T,
int... n>
1494 for (
int i = 0; i < n; ++i) {
1495 for (
int j = i + 1; j < n; ++j) {
1496 if (
std::abs(A(i, j) - A(j, i)) > tolerance) {
1517 if (A(0, 0) < 0.0) {
1534 auto subtensor = make_tensor<2, 2>([A](
int i,
int j) {
return A(i, j); });
1546 template <
typename T,
int n>
1567 template <
typename T,
int n,
int... m>
1572 for (
int i = 0; i < n; i++) {
1574 for (
int j = 0; j < i; j++) {
1575 c -= L[i][j] * y[j];
1586 template <
typename T,
int n,
int... m>
1607 template <
typename T,
int n,
int... m>
1611 for (
int i = n - 1; i >= 0; i--) {
1613 for (
int j = i + 1; j < n; j++) {
1614 c -= U[i][j] * x[j];
1625 template <
typename S,
typename T,
int n,
int... m>
1641 template <
typename T,
int n>
1654 double inv_detA(1.0 /
det(A));
1658 invA[0][0] = A[1][1] * inv_detA;
1659 invA[0][1] = -A[0][1] * inv_detA;
1660 invA[1][0] = -A[1][0] * inv_detA;
1661 invA[1][1] = A[0][0] * inv_detA;
1672 double inv_detA(1.0 /
det(A));
1676 invA[0][0] = (A[1][1] * A[2][2] - A[1][2] * A[2][1]) * inv_detA;
1677 invA[0][1] = (A[0][2] * A[2][1] - A[0][1] * A[2][2]) * inv_detA;
1678 invA[0][2] = (A[0][1] * A[1][2] - A[0][2] * A[1][1]) * inv_detA;
1679 invA[1][0] = (A[1][2] * A[2][0] - A[1][0] * A[2][2]) * inv_detA;
1680 invA[1][1] = (A[0][0] * A[2][2] - A[0][2] * A[2][0]) * inv_detA;
1681 invA[1][2] = (A[0][2] * A[1][0] - A[0][0] * A[1][2]) * inv_detA;
1682 invA[2][0] = (A[1][0] * A[2][1] - A[1][1] * A[2][0]) * inv_detA;
1683 invA[2][1] = (A[0][1] * A[2][0] - A[0][0] * A[2][1]) * inv_detA;
1684 invA[2][2] = (A[0][0] * A[1][1] - A[0][1] * A[1][0]) * inv_detA;
1693 template <
typename T,
int n>
1696 auto I = DenseIdentity<n>();
1708 template <
typename T,
int m,
int... n>
1712 for (
int i = 1; i < m; i++) {
1713 out <<
", " << A[i];
1741 template <
int m,
int... n>
1746 for (
int i = 1; i < m; i++) {
1761 for (
int i = 0; i < n; i++) {
1762 if (copy[i] * copy[i] < 1.0e-20) {
1770 template <
int m,
int n>
1774 for (
int i = 0; i < m; i++) {
1775 for (
int j = 0; j < n; j++) {
1776 if (copy[i][j] * copy[i][j] < 1.0e-20) {
1787 template <
typename T1,
typename T2>
1790 template <
int... m,
int... n>
1791 struct outer_prod<
tensor<double, m...>,
tensor<double, n...>> {
1796 struct outer_prod<double,
tensor<double, n...>> {
1801 struct outer_prod<
tensor<double, n...>, double> {
1806 struct outer_prod<double, double> {
1807 using type = tensor<double>;
1810 template <
typename T>
1811 struct outer_prod<zero, T> {
1815 template <
typename T>
1816 struct outer_prod<T, zero> {
1828 template <
typename T1,
typename T2>
1857 template <
typename T>
1867 template <
typename T>
1897 for_constexpr<n...>([&](
auto... i) { total += df_dx(i...) * dx(i...); });
1905 template <
int m,
int... n>
1909 for (
int i = 0; i < m; i++) {
1919 template <
int m,
int n,
int... p>
1923 for (
int i = 0; i < m; i++) {
1924 for (
int j = 0; j < n; j++) {
1938 template <
typename T,
int... n>
1941 return (n * ... * 1);
1961 template <
int i,
typename T,
int... n>
1964 constexpr
int dimensions[] = {n...};
1965 return dimensions[i];
1976 template <
typename T,
int m,
int... n>
1983 template <
typename T,
int... n>
1986 bool found_nan =
false;
1998 inline float angle_between(
const vec < 2 > & a,
const vec < 2 > & b) {
2002 inline float angle_between(
const vec < 3 > & a,
const vec < 3 > & b) {
2007 inline float angle_between(
const mat < 3, 3 > & U,
const mat < 3, 3 > & V) {
2011 inline mat < 2, 2 > rotation(
const float theta) {
2013 {
cos(theta), -
sin(theta)},
2014 {
sin(theta),
cos(theta) }
2018 inline mat < 3, 3 > axis_to_rotation(
const vec < 3 > & omega) {
2020 float norm_omega =
norm(omega);
2022 if (fabs(norm_omega) < 0.000001f) {
2028 vec3 u = omega / norm_omega;
2030 float c =
cos(norm_omega);
2031 float s =
sin(norm_omega);
2033 return mat < 3, 3 >{
2035 u[0]*u[0]*(1.0f - c) + c,
2036 u[0]*u[1]*(1.0f - c) - u[2]*s,
2037 u[0]*u[2]*(1.0f - c) + u[1]*s
2039 u[1]*u[0]*(1.0f - c) + u[2]*s,
2040 u[1]*u[1]*(1.0f - c) + c,
2041 u[1]*u[2]*(1.0f - c) - u[0]*s
2043 u[2]*u[0]*(1.0f - c) - u[1]*s,
2044 u[2]*u[1]*(1.0f - c) + u[0]*s,
2045 u[2]*u[2]*(1.0f - c) + c
2054 inline vec < 3 > rotation_to_axis(
const mat < 3, 3 > & R) {
2056 float theta =
acos(clip(0.5f * (
tr(R) - 1.0f), -1.0f, 1.0f));
2061 if (fabs(theta) < 0.00001f) {
2062 scale = 0.5f + theta * theta / 12.0f;
2065 scale = 0.5f * theta /
sin(theta);
2068 return vec3{ R(2,1) - R(1,2), R(0,2) - R(2,0), R(1,0) - R(0,1) } *scale;
2072 inline mat < 3, 3 > look_at(
const vec < 3 > & direction,
const vec < 3 > & up =
vec3{ 0.0f, 0.0f, 1.0f }) {
2084 inline mat < 2, 2 > look_at(
const vec < 2 > & direction) {
2094 inline mat < 3, 3 > R3_basis(
const vec3 & n) {
2095 float sign = (n[2] >= 0.0f) ? 1.0f : -1.0f;
2096 float a = -1.0f / (sign + n[2]);
2097 float b = n[0] * n[1] * a;
2099 return mat < 3, 3 >{
2101 1.0f + sign * n[0] * n[0] * a,
2106 sign + n[1] * n[1] * a,
2119 #include "smith/numerics/functional/tuple_tensor_dual_functions.hpp"
This file contains the interface used for initializing/terminating any hardware accelerator-related f...
#define SMITH_SUPPRESS_NVCC_HOSTDEVICE_WARNING
Macro to turn off specific nvcc warnings.
#define SMITH_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc or amdclang and does nothing on ...
Implementation of isotropic tensor classes.
Accelerator functionality.
SMITH_HOST_DEVICE bool is_symmetric(tensor< double, n, n > A, double tolerance=1.0e-8)
Return whether a square rank 2 tensor is symmetric.
SMITH_HOST_DEVICE auto acos(dual< gradient_type > a)
implementation of acos for dual numbers
SMITH_HOST_DEVICE bool is_symmetric_and_positive_definite(tensor< double, 3, 3 > A)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SMITH_HOST_DEVICE auto antisym(const tensor< T, n, n > &A)
Returns the antisymmetric part of a square matrix.
constexpr SMITH_HOST_DEVICE auto & operator+=(dual< gradient_type > &a, const dual< gradient_type > &b)
compound assignment (+) for dual numbers
constexpr SMITH_HOST_DEVICE auto & operator-=(tensor< T, n... > &A, zero)
compound assignment (-) between a tensor and zero (no-op)
SMITH_SUPPRESS_NVCC_HOSTDEVICE_WARNING constexpr SMITH_HOST_DEVICE auto make_tensor(lambda_type f)
Creates a tensor of requested dimension by subsequent calls to a functor Can be thought of as analogo...
constexpr T & get(variant< T0, T1 > &v)
Returns the variant member of specified type.
constexpr SMITH_HOST_DEVICE int size(zero)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SMITH_HOST_DEVICE auto det(const tensor< T, 3, 3 > &A)
This is an overloaded member function, provided for convenience. It differs from the above function o...
double relative_error(tensor< T, n... > A, tensor< T, n... > B)
computes the relative error (in the frobenius norm) between two tensors of the same shape
constexpr SMITH_HOST_DEVICE int leading_dimension(tensor< T, m, n... >)
a function for querying the first dimension of a tensor
constexpr SMITH_HOST_DEVICE auto get_gradient(const tensor< double, n... > &)
get the gradient of type tensor (note: since its stored type is not a dual number,...
constexpr SMITH_HOST_DEVICE auto operator*(const tensor< S, m... > &A, const tensor< T, n... > &B)
this is a shorthand for dot(A, B)
SMITH_HOST_DEVICE auto contract(const zero &, const T &)
This is an overloaded member function, provided for convenience. It differs from the above function o...
FieldStateWeightedSum operator+(const FieldState &x, const FieldState &y)
add two FieldState
SMITH_HOST_DEVICE auto cross(const tensor< S, 3 > &u, const tensor< T, 3 > &v)
compute the (right handed) cross product of two 3-vectors
SMITH_HOST_DEVICE auto chain_rule(const tensor< double, m, n, p... > &df_dx, const tensor< double, p... > &dx)
constexpr SMITH_HOST_DEVICE auto I2(const tensor< T, 3, 3 > &A)
Returns the second invariant of a 3x3 matrix.
SMITH_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
constexpr SMITH_HOST_DEVICE auto inner(zero, double)
constexpr SMITH_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
constexpr SMITH_HOST_DEVICE auto solve_lower_triangular(const tensor< T, n, n > &L, const tensor< T, n, m... > &b)
constexpr SMITH_HOST_DEVICE auto norm(zero)
overload of Frobenius norm for zero type
constexpr SMITH_HOST_DEVICE auto squared_norm(const tensor< T, n... > &A)
This is an overloaded member function, provided for convenience. It differs from the above function o...
tensor< double, 3 > vec3
statically sized vector of 3 doubles
auto & operator<<(std::ostream &out, zero)
Write a zero out to an output stream.
SMITH_HOST_DEVICE auto cos(dual< gradient_type > a)
implementation of cosine for dual numbers
std::conditional_t<(n1==1 &&n2==1), double, std::conditional_t< n1==1, tensor< T, n2 >, std::conditional_t< n2==1, tensor< T, n1 >, tensor< T, n1, n2 > >> > reduced_tensor
Removes 1s from tensor dimensions For example, a tensor<T, 1, 10> is equivalent to a tensor<T,...
constexpr SMITH_HOST_DEVICE auto chop(const tensor< double, m, n > &A)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SMITH_HOST_DEVICE tensor< T, n > diag(const tensor< T, n, n > &D)
Returns an array containing the diagonal entries of a square matrix.
bool isnan(const zero &)
This is an overloaded member function, provided for convenience. It differs from the above function o...
tensor(const T(&data)[n1]) -> tensor< T, n1 >
class template argument deduction guide for type tensor.
constexpr auto double_dot(const tensor< S, m, n > &A, const tensor< T, m, n > &B)
constexpr SMITH_HOST_DEVICE auto inv(const tensor< T, n, n > &A)
constexpr SMITH_HOST_DEVICE auto dev(const tensor< T, n, n > &A)
Calculates the deviator of a matrix (rank-2 tensor)
auto matrix_sqrt(const tensor< T, dim, dim > &A)
compute the matrix square root of a square, real-valued, symmetric matrix i.e. given A,...
SMITH_HOST_DEVICE auto abs(dual< gradient_type > x)
Implementation of absolute value function for dual numbers.
tensor(const T(&data)[n1][n2]) -> tensor< T, n1, n2 >
class template argument deduction guide for type tensor.
SMITH_HOST_DEVICE tensor< T, 3, 3 > to_3x3(const tensor< T, 2, 2 > &A)
promotes a 2x2 matrix to a 3x3 matrix, by populating the upper left block, leaving zeroes in the thir...
SMITH_HOST_DEVICE consteval int first_dim(const tensor< T, m, n... > &)
return the size of the leftmost tensor dimension
SMITH_HOST_DEVICE void print(const tensor< double, m, n... > &A)
print a tensor using printf, so that it is suitable for use inside cuda kernels.
constexpr SMITH_HOST_DEVICE auto sym(const tensor< T, n, n > &A)
Returns the symmetric part of a square matrix.
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 SMITH_HOST_DEVICE auto tensor_with_shape(std::integer_sequence< int, n... >)
Creates a tensor given the dimensions in a std::integer_sequence.
constexpr SMITH_HOST_DEVICE auto dot(const tensor< S, m, n, p, q > &A, const tensor< T, q > &B)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SMITH_HOST_DEVICE auto outer(const tensor< S, m > &A, const tensor< T, n > &B)
constexpr SMITH_HOST_DEVICE auto & operator+=(tensor< T, n... > &A, zero)
compound assignment (+) between a tensor and zero (no-op)
constexpr SMITH_HOST_DEVICE auto linear_solve(const LuFactorization< T, n > &, const zero)
tensor< double, 2 > vec2
statically sized vector of 2 doubles
constexpr SMITH_HOST_DEVICE int dimension(const tensor< T, n... > &)
a function for querying the ith dimension of a tensor
FieldStateWeightedSum operator*(double a, const FieldState &b)
multiply scalar by a FieldState to get a temporary FieldStateWeightedSum which can cast back to a Fie...
tensor< double, 2, 2 > mat2
statically sized 2x2 matrix of doubles
constexpr SMITH_HOST_DEVICE auto transpose(const tensor< T, m, n > &A)
Returns the transpose of the matrix.
constexpr SMITH_HOST_DEVICE tensor< double, dim, dim > DenseIdentity()
Obtains the identity matrix of the specified dimension.
constexpr SMITH_HOST_DEVICE auto & operator-=(dual< gradient_type > &a, const dual< gradient_type > &b)
compound assignment (-) for dual numbers
SMITH_HOST_DEVICE auto sin(dual< gradient_type > a)
implementation of sine for dual numbers
SMITH_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
constexpr SMITH_HOST_DEVICE auto detApIm1(const tensor< T, 3, 3 > &A)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SMITH_HOST_DEVICE auto tr(const tensor< T, n, n > &A)
Returns the trace of a square matrix.
constexpr SMITH_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
constexpr SMITH_HOST_DEVICE auto operator-(const tensor< S, m, n... > &A, const tensor< T, m, n... > &B)
return the difference of two tensors
constexpr SMITH_HOST_DEVICE auto solve_upper_triangular(const tensor< T, n, n > &U, const tensor< T, n, m... > &y)
Solves an upper triangular system Ux = y.
FieldStateWeightedSum operator-(const FieldState &x, const FieldState &y)
subtract two FieldState
constexpr SMITH_HOST_DEVICE auto operator/(const dual< gradient_type > &a, double b)
division of a dual number by a non-dual number
tensor< double, 3, 3 > mat3
statically sized 3x3 matrix of doubles
constexpr SMITH_HOST_DEVICE auto diagonal_matrix(const tensor< T, n, n > &A)
Returns a square matrix (rank-2 tensor) containing the diagonal entries of the input square matrix wi...
Representation of an LU factorization.
tensor< T, n, n > U
Upper triangular factor.
tensor< int, n > P
Row permutation indices due to partial pivoting.
tensor< T, n, n > L
Lower triangular factor. Has ones on diagonal.
Arbitrary-rank tensor class.
A sentinel struct for eliding no-op tensor operations.
SMITH_HOST_DEVICE auto operator()(T...) const
zero can be accessed like a multidimensional array
SMITH_HOST_DEVICE auto operator=(T)
anything assigned to zero does not change its value and returns zero