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 SERAC_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 SERAC_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>
162 template <
typename T>
169 template <
typename T>
184 template <
typename T>
191 template <
typename T>
203 template <
typename T>
210 template <
typename T>
217 template <
typename T>
224 template <
typename T>
227 static_assert(::detail::always_false<T>{},
"Error: Can't divide by zero!");
251 template <
typename T>
258 template <
typename T>
271 template <
typename T,
int n1,
int n2 = 1>
273 (n1 == 1 && n2 == 1),
double,
281 template <
typename T,
int... n>
298 template <
typename lambda_type>
301 using T = decltype(f());
317 template <
int n1,
typename lambda_type>
320 using T = decltype(f(n1));
322 for (
int i = 0; i < n1; i++) {
341 template <
int n1,
int n2,
typename lambda_type>
344 using T = decltype(f(n1, n2));
346 for (
int i = 0; i < n1; i++) {
347 for (
int j = 0; j < n2; j++) {
368 template <
int n1,
int n2,
int n3,
typename lambda_type>
371 using T = decltype(f(n1, n2, n3));
373 for (
int i = 0; i < n1; i++) {
374 for (
int j = 0; j < n2; j++) {
375 for (
int k = 0; k < n3; k++) {
376 A(i, j, k) = f(i, j, k);
398 template <
int n1,
int n2,
int n3,
int n4,
typename lambda_type>
401 using T = decltype(f(n1, n2, n3, n4));
403 for (
int i = 0; i < n1; i++) {
404 for (
int j = 0; j < n2; j++) {
405 for (
int k = 0; k < n3; k++) {
406 for (
int l = 0; l < n4; l++) {
407 A(i, j, k, l) = f(i, j, k, l);
423 template <
typename S,
typename T,
int m,
int... n>
426 tensor<decltype(S{} + T{}), m, n...> C{};
427 for (
int i = 0; i < m; i++) {
439 template <
typename T,
int m,
int... n>
443 for (
int i = 0; i < m; i++) {
457 template <
typename S,
typename T,
int m,
int... n>
460 tensor<decltype(S{} + T{}), m, n...> C{};
461 for (
int i = 0; i < m; i++) {
475 template <
typename S,
typename T,
int m,
int... n>
478 for (
int i = 0; i < m; i++) {
491 template <
typename T>
504 template <
typename T,
int n>
507 for (
int i = 0; i < n; i++) {
508 A.data[i][0] += B[i];
519 template <
typename T,
int n>
522 for (
int i = 0; i < n; i++) {
523 A.data[0][i] += B[i];
534 template <
typename T>
537 return A.data[0] += B;
546 template <
typename T>
549 return A.data[0][0] += B;
558 template <
typename T,
int... n>
572 template <
typename S,
typename T,
int m,
int... n>
575 for (
int i = 0; i < m; i++) {
587 template <
typename T,
int... n>
597 template <
typename T,
int n>
600 tensor<decltype(
double{} * T{}), n> AB{};
601 for (
int i = 0; i < n; i++) {
611 template <
typename T,
int m>
614 tensor<decltype(T{} *
double{}), m> AB{};
615 for (
int i = 0; i < m; i++) {
625 template <
typename T,
int n>
635 template <
typename T,
int n>
645 template <
typename S,
typename T,
int m,
int n>
648 tensor<decltype(S{} * T{}), m, n> AB{};
649 for (
int i = 0; i < m; i++) {
650 for (
int j = 0; j < n; j++) {
651 AB[i][j] = A[i] * B[j];
666 template <
typename S,
typename T,
int m,
int n>
669 decltype(S{} * T{}) sum{};
670 for (
int i = 0; i < m; i++) {
671 for (
int j = 0; j < n; j++) {
672 sum += A[i][j] * B[i][j];
682 template <
typename S,
typename T,
int m>
685 decltype(S{} * T{}) sum{};
686 for (
int i = 0; i < m; i++) {
704 template <
typename S,
int m,
int n>
714 template <
typename S,
int m>
732 template <
typename T,
int m,
int n>
742 template <
typename T,
int m>
762 template <
typename S,
typename T,
int m,
int n,
int p>
765 tensor<decltype(S{} * T{}), m, p> AB{};
766 for (
int i = 0; i < m; i++) {
767 for (
int j = 0; j < p; j++) {
768 for (
int k = 0; k < n; k++) {
769 AB[i][j] = AB[i][j] + A[i][k] * B[k][j];
780 template <
typename T,
int m>
790 template <
typename T,
int m>
800 template <
typename S,
typename T,
int m>
803 decltype(S{} * T{}) AB{};
804 for (
int i = 0; i < m; i++) {
805 AB = AB + A[i] * B[i];
814 template <
typename S,
typename T,
int m,
int n>
817 tensor<decltype(S{} * T{}), n> AB{};
818 for (
int i = 0; i < n; i++) {
819 for (
int j = 0; j < m; j++) {
820 AB[i] = AB[i] + A[j] * B[j][i];
830 template <
typename S,
typename T,
int m,
int n,
int p>
833 tensor<decltype(S{} * T{}), n, p> AB{};
834 for (
int j = 0; j < m; j++) {
835 AB = AB + A[j] * B[j];
850 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
853 tensor<decltype(S{} * T{}), n, p, q> AB{};
854 for (
int j = 0; j < m; j++) {
855 AB = AB + A[j] * B[j];
864 template <
typename S,
typename T,
int m,
int n>
867 tensor<decltype(S{} * T{}), m> AB{};
868 for (
int i = 0; i < m; i++) {
869 for (
int j = 0; j < n; j++) {
870 AB[i] = AB[i] + A[i][j] * B[j];
880 template <
typename S,
typename T,
int m,
int n,
int p,
int q,
int r>
883 tensor<decltype(S{} * T{}), m, p, q, r> AB{};
884 for (
int i = 0; i < m; i++) {
885 for (
int j = 0; j < n; j++) {
886 AB[i] = AB[i] + A[i][j] * B[j];
896 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
899 tensor<decltype(S{} * T{}), m, p, q> AB{};
900 for (
int i = 0; i < m; i++) {
901 for (
int j = 0; j < n; j++) {
902 AB[i] = AB[i] + A[i][j] * B[j];
912 template <
typename S,
typename T,
int m,
int n,
int p>
915 tensor<decltype(S{} * T{}), m, n> AB{};
916 for (
int i = 0; i < m; i++) {
917 for (
int j = 0; j < n; j++) {
918 for (
int k = 0; k < p; k++) {
919 AB[i][j] += A[i][j][k] * B[k];
930 template <
typename S,
typename T,
typename U,
int m,
int n>
933 decltype(S{} * T{} * U{}) uAv{};
934 for (
int i = 0; i < m; i++) {
935 for (
int j = 0; j < n; j++) {
936 uAv += u[i] * A[i][j] * v[j];
943 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
946 tensor<decltype(S{} * T{}), m, n, p> AB{};
947 for (
int i = 0; i < m; i++) {
948 for (
int j = 0; j < n; j++) {
949 for (
int k = 0; k < p; k++) {
950 for (
int l = 0; l < q; l++) {
951 AB[i][j][k] += A[i][j][k][l] * B[l];
960 template <
typename T>
963 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),
964 A(0, 0) * A(1, 1) - A(1, 0) * A(0, 1)};
968 template <
typename T>
975 template <
typename T>
982 template <
typename S,
typename T>
985 return tensor<decltype(S{} * T{}), 3>{u(1) * v(2) - u(2) * v(1), u(2) * v(0) - u(0) * v(2),
986 u(0) * v(1) - u(1) * v(0)};
1000 template <
typename S,
typename T,
int m,
int n,
int p,
int q>
1003 tensor<decltype(S{} * T{}), m, n> AB{};
1004 for (
int i = 0; i < m; i++) {
1005 for (
int j = 0; j < n; j++) {
1006 for (
int k = 0; k < p; k++) {
1007 for (
int l = 0; l < q; l++) {
1008 AB[i][j] += A[i][j][k][l] * B[k][l];
1020 template <
typename S,
typename T,
int m,
int n,
int p>
1023 tensor<decltype(S{} * T{}), m> AB{};
1024 for (
int i = 0; i < m; i++) {
1025 for (
int j = 0; j < n; j++) {
1026 for (
int k = 0; k < p; k++) {
1027 AB[i] += A[i][j][k] * B[j][k];
1038 template <
typename S,
typename T,
int m,
int n>
1041 decltype(S{} * T{}) AB{};
1042 for (
int i = 0; i < m; i++) {
1043 for (
int j = 0; j < n; j++) {
1044 AB += A[i][j] * B[i][j];
1053 template <
typename S,
typename T,
int... m,
int... n>
1063 template <
typename T,
int m>
1067 for (
int i = 0; i < m; i++) {
1068 total += A[i] * A[i];
1074 template <
typename T,
int m,
int n>
1078 for (
int i = 0; i < m; i++) {
1079 for (
int j = 0; j < n; j++) {
1080 total += A[i][j] * A[i][j];
1087 template <
typename T,
int... n>
1091 for_constexpr<n...>([&](
auto... i) { total += A(i...) * A(i...); });
1099 template <
typename T,
int... n>
1116 template <
typename T,
int... n>
1127 template <
typename T>
1131 output[0][0] = A[0][0];
1132 output[0][1] = A[0][1];
1133 output[1][0] = A[1][0];
1134 output[1][1] = A[1][1];
1143 template <
typename T,
int n>
1147 for (
int i = 0; i < n; i++) {
1148 trA = trA + A[i][i];
1158 template <
typename T,
int n>
1162 for (
int i = 0; i < n; i++) {
1163 for (
int j = 0; j < n; j++) {
1164 symA[i][j] = 0.5 * (A[i][j] + A[j][i]);
1175 template <
typename T,
int n>
1179 for (
int i = 0; i < n; i++) {
1180 for (
int j = 0; j < n; j++) {
1181 antisymA[i][j] = 0.5 * (A[i][j] - A[j][i]);
1194 template <
typename T,
int n>
1199 for (
int i = 0; i < n; i++) {
1200 devA[i][i] -= trA / n;
1211 template <
typename T,
int n>
1215 for (
int i = 0; i < n; i++) {
1225 template <
typename T,
int n>
1229 for (
int i = 0; i < n; i++) {
1239 template <
typename T,
int n>
1243 for (
int i = 0; i < n; i++) {
1257 for (
int i = 0; i < dim; i++) {
1258 for (
int j = 0; j < dim; j++) {
1269 template <
typename T,
int m,
int n>
1273 for (
int i = 0; i < n; i++) {
1274 for (
int j = 0; j < m; j++) {
1285 template <
typename T>
1288 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] -
1296 template <
typename T>
1299 return A[0][0] * A[1][1] - A[0][1] * A[1][0];
1303 template <
typename T>
1306 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] -
1307 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];
1320 template <
typename T>
1329 return A(0, 0) - A(0, 1) * A(1, 0) + A(1, 1) + A(0, 0) * A(1, 1);
1333 template <
typename T>
1340 return A(0, 0) + A(1, 1) + A(2, 2)
1341 - A(0, 1) * A(1, 0) * (1 + A(2, 2))
1342 + A(0, 0) * A(1, 1) * (1 + A(2, 2))
1343 - A(0, 2) * A(2, 0) * (1 + A(1, 1))
1344 - A(1, 2) * A(2, 1) * (1 + A(0, 0))
1347 + A(0, 1) * A(1, 2) * A(2, 0)
1348 + A(0, 2) * A(1, 0) * A(2, 1);
1360 template <
typename T,
int dim>
1364 for (
int i = 0; i < 15; i++) {
1365 B = 0.5 * (B +
dot(A,
inv(B)));
1399 template <
int i1,
int i2,
typename S,
int m,
int... n,
typename T,
int p,
int q>
1402 constexpr
int Adims[] = {m, n...};
1403 constexpr
int Bdims[] = {p, q};
1404 static_assert(
sizeof...(n) < 3);
1405 static_assert(Adims[i1] == Bdims[i2],
"error: incompatible tensor dimensions");
1408 constexpr
int new_dim = (i2 == 0) ? q : p;
1409 constexpr
int d1 = (i1 == 0) ? new_dim : Adims[0];
1410 constexpr
int d2 = (i1 == 1) ? new_dim : Adims[1];
1411 constexpr
int d3 =
sizeof...(n) == 1 ? 0 : ((i1 == 2) ? new_dim : Adims[2]);
1414 using U = decltype(S{} * T{});
1421 if constexpr (d3 == 0) {
1422 for (
int i = 0; i < d1; i++) {
1423 for (
int j = 0; j < d2; j++) {
1425 for (
int k = 0; k < Adims[i1]; k++) {
1426 if constexpr (i1 == 0 && i2 == 0) sum += A(k, j) * B(k, i);
1427 if constexpr (i1 == 1 && i2 == 0) sum += A(i, k) * B(k, j);
1428 if constexpr (i1 == 0 && i2 == 1) sum += A(k, j) * B(i, k);
1429 if constexpr (i1 == 1 && i2 == 1) sum += A(i, k) * B(j, k);
1435 for (
int i = 0; i < d1; i++) {
1436 for (
int j = 0; j < d2; j++) {
1437 for (
int k = 0; k < d3; k++) {
1439 for (
int l = 0; l < Adims[i1]; l++) {
1440 if constexpr (i1 == 0 && i2 == 0) sum += A(l, j, k) * B(l, i);
1441 if constexpr (i1 == 1 && i2 == 0) sum += A(i, l, k) * B(l, j);
1442 if constexpr (i1 == 2 && i2 == 0) sum += A(i, j, l) * B(l, k);
1443 if constexpr (i1 == 0 && i2 == 1) sum += A(l, j, k) * B(i, l);
1444 if constexpr (i1 == 1 && i2 == 1) sum += A(i, l, k) * B(j, l);
1445 if constexpr (i1 == 2 && i2 == 1) sum += A(i, j, l) * B(k, l);
1457 template <
int i1,
int i2,
typename T>
1472 template <
typename T,
int... n>
1489 for (
int i = 0; i < n; ++i) {
1490 for (
int j = i + 1; j < n; ++j) {
1491 if (
std::abs(A(i, j) - A(j, i)) > tolerance) {
1512 if (A(0, 0) < 0.0) {
1529 auto subtensor = make_tensor<2, 2>([A](
int i,
int j) {
return A(i, j); });
1541 template <
typename T,
int n>
1562 template <
typename T,
int n,
int... m>
1567 for (
int i = 0; i < n; i++) {
1569 for (
int j = 0; j < i; j++) {
1570 c -= L[i][j] * y[j];
1581 template <
typename T,
int n,
int... m>
1602 template <
typename T,
int n,
int... m>
1606 for (
int i = n - 1; i >= 0; i--) {
1608 for (
int j = i + 1; j < n; j++) {
1609 c -= U[i][j] * x[j];
1620 template <
typename S,
typename T,
int n,
int... m>
1636 template <
typename T,
int n>
1649 double inv_detA(1.0 /
det(A));
1653 invA[0][0] = A[1][1] * inv_detA;
1654 invA[0][1] = -A[0][1] * inv_detA;
1655 invA[1][0] = -A[1][0] * inv_detA;
1656 invA[1][1] = A[0][0] * inv_detA;
1667 double inv_detA(1.0 /
det(A));
1671 invA[0][0] = (A[1][1] * A[2][2] - A[1][2] * A[2][1]) * inv_detA;
1672 invA[0][1] = (A[0][2] * A[2][1] - A[0][1] * A[2][2]) * inv_detA;
1673 invA[0][2] = (A[0][1] * A[1][2] - A[0][2] * A[1][1]) * inv_detA;
1674 invA[1][0] = (A[1][2] * A[2][0] - A[1][0] * A[2][2]) * inv_detA;
1675 invA[1][1] = (A[0][0] * A[2][2] - A[0][2] * A[2][0]) * inv_detA;
1676 invA[1][2] = (A[0][2] * A[1][0] - A[0][0] * A[1][2]) * inv_detA;
1677 invA[2][0] = (A[1][0] * A[2][1] - A[1][1] * A[2][0]) * inv_detA;
1678 invA[2][1] = (A[0][1] * A[2][0] - A[0][0] * A[2][1]) * inv_detA;
1679 invA[2][2] = (A[0][0] * A[1][1] - A[0][1] * A[1][0]) * inv_detA;
1688 template <
typename T,
int n>
1691 auto I = DenseIdentity<n>();
1703 template <
typename T,
int m,
int... n>
1707 for (
int i = 1; i < m; i++) {
1708 out <<
", " << A[i];
1736 template <
int m,
int... n>
1741 for (
int i = 1; i < m; i++) {
1756 for (
int i = 0; i < n; i++) {
1757 if (copy[i] * copy[i] < 1.0e-20) {
1765 template <
int m,
int n>
1769 for (
int i = 0; i < m; i++) {
1770 for (
int j = 0; j < n; j++) {
1771 if (copy[i][j] * copy[i][j] < 1.0e-20) {
1782 template <
typename T1,
typename T2>
1785 template <
int... m,
int... n>
1786 struct outer_prod<
tensor<double, m...>,
tensor<double, n...>> {
1791 struct outer_prod<double,
tensor<double, n...>> {
1796 struct outer_prod<
tensor<double, n...>, double> {
1801 struct outer_prod<double, double> {
1802 using type = tensor<double>;
1805 template <
typename T>
1806 struct outer_prod<zero, T> {
1810 template <
typename T>
1811 struct outer_prod<T, zero> {
1823 template <
typename T1,
typename T2>
1852 template <
typename T>
1862 template <
typename T>
1892 for_constexpr<n...>([&](
auto... i) { total += df_dx(i...) * dx(i...); });
1900 template <
int m,
int... n>
1904 for (
int i = 0; i < m; i++) {
1914 template <
int m,
int n,
int... p>
1918 for (
int i = 0; i < m; i++) {
1919 for (
int j = 0; j < n; j++) {
1933 template <
typename T,
int... n>
1936 return (n * ... * 1);
1956 template <
int i,
typename T,
int... n>
1959 constexpr
int dimensions[] = {n...};
1960 return dimensions[i];
1971 template <
typename T,
int m,
int... n>
1978 template <
typename T,
int... n>
1981 bool found_nan =
false;
1993 inline float angle_between(
const vec < 2 > & a,
const vec < 2 > & b) {
1997 inline float angle_between(
const vec < 3 > & a,
const vec < 3 > & b) {
2002 inline float angle_between(
const mat < 3, 3 > & U,
const mat < 3, 3 > & V) {
2006 inline mat < 2, 2 > rotation(
const float theta) {
2008 {
cos(theta), -
sin(theta)},
2009 {
sin(theta),
cos(theta) }
2013 inline mat < 3, 3 > axis_to_rotation(
const vec < 3 > & omega) {
2015 float norm_omega =
norm(omega);
2017 if (fabs(norm_omega) < 0.000001f) {
2023 vec3 u = omega / norm_omega;
2025 float c =
cos(norm_omega);
2026 float s =
sin(norm_omega);
2028 return mat < 3, 3 >{
2030 u[0]*u[0]*(1.0f - c) + c,
2031 u[0]*u[1]*(1.0f - c) - u[2]*s,
2032 u[0]*u[2]*(1.0f - c) + u[1]*s
2034 u[1]*u[0]*(1.0f - c) + u[2]*s,
2035 u[1]*u[1]*(1.0f - c) + c,
2036 u[1]*u[2]*(1.0f - c) - u[0]*s
2038 u[2]*u[0]*(1.0f - c) - u[1]*s,
2039 u[2]*u[1]*(1.0f - c) + u[0]*s,
2040 u[2]*u[2]*(1.0f - c) + c
2049 inline vec < 3 > rotation_to_axis(
const mat < 3, 3 > & R) {
2051 float theta =
acos(clip(0.5f * (
tr(R) - 1.0f), -1.0f, 1.0f));
2056 if (fabs(theta) < 0.00001f) {
2057 scale = 0.5f + theta * theta / 12.0f;
2060 scale = 0.5f * theta /
sin(theta);
2063 return vec3{ R(2,1) - R(1,2), R(0,2) - R(2,0), R(1,0) - R(0,1) } *scale;
2067 inline mat < 3, 3 > look_at(
const vec < 3 > & direction,
const vec < 3 > & up =
vec3{ 0.0f, 0.0f, 1.0f }) {
2079 inline mat < 2, 2 > look_at(
const vec < 2 > & direction) {
2089 inline mat < 3, 3 > R3_basis(
const vec3 & n) {
2090 float sign = (n[2] >= 0.0f) ? 1.0f : -1.0f;
2091 float a = -1.0f / (sign + n[2]);
2092 float b = n[0] * n[1] * a;
2094 return mat < 3, 3 >{
2096 1.0f + sign * n[0] * n[0] * a,
2101 sign + n[1] * n[1] * a,
2114 #include "serac/numerics/functional/tuple_tensor_dual_functions.hpp"
This file contains the interface used for initializing/terminating any hardware accelerator-related f...
#define SERAC_SUPPRESS_NVCC_HOSTDEVICE_WARNING
Macro to turn off specific nvcc warnings.
#define SERAC_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
Implementation of isotropic tensor classes.
Accelerator functionality.
constexpr SERAC_HOST_DEVICE auto inv(const tensor< T, n, n > &A)
constexpr SERAC_HOST_DEVICE tensor< double, dim, dim > DenseIdentity()
Obtains the identity matrix of the specified dimension.
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 tensor< T, n > diag(const tensor< T, n, n > &D)
Returns an array containing the diagonal entries of a square matrix.
SERAC_HOST_DEVICE auto sin(dual< gradient_type > a)
implementation of sine for dual numbers
constexpr SERAC_HOST_DEVICE auto operator-(const tensor< S, m, n... > &A, const tensor< T, m, n... > &B)
return the difference of two tensors
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 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 bool is_symmetric(tensor< double, n, n > A, double tolerance=1.0e-8)
Return whether a square rank 2 tensor is symmetric.
constexpr SERAC_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 SERAC_HOST_DEVICE auto operator*(const dual< gradient_type > &a, double b)
multiplication of a dual number and a non-dual number
tensor(const T(&data)[n1]) -> tensor< T, n1 >
class template argument deduction guide for type tensor.
SERAC_HOST_DEVICE auto cos(dual< gradient_type > a)
implementation of cosine for dual numbers
Domain operator-(const Domain &a, const Domain &b)
create a new domain that is the set difference of a and b
constexpr SERAC_HOST_DEVICE auto & operator+=(tensor< T, n... > &A, zero)
compound assignment (+) between a tensor and zero (no-op)
constexpr SERAC_HOST_DEVICE int size(zero)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SERAC_HOST_DEVICE int leading_dimension(tensor< T, m, n... >)
a function for querying the first dimension of a tensor
constexpr T & get(variant< T0, T1 > &v)
Returns the variant member of specified type.
constexpr SERAC_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 SERAC_HOST_DEVICE auto transpose(const tensor< T, m, n > &A)
Returns the transpose of the matrix.
constexpr SERAC_HOST_DEVICE auto I2(const tensor< T, 3, 3 > &A)
Returns the second invariant of a 3x3 matrix.
SERAC_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...
SERAC_SUPPRESS_NVCC_HOSTDEVICE_WARNING constexpr SERAC_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 SERAC_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,...
SERAC_HOST_DEVICE auto chain_rule(const tensor< double, m, n, p... > &df_dx, const tensor< double, p... > &dx)
SERAC_HOST_DEVICE auto contract(const zero &, const T &)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SERAC_HOST_DEVICE auto linear_solve(const LuFactorization< T, n > &, const zero)
SERAC_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
constexpr SERAC_HOST_DEVICE auto operator+(dual< gradient_type > a, double b)
addition of a dual number and a non-dual number
constexpr SERAC_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...
constexpr SERAC_HOST_DEVICE auto dev(const tensor< T, n, n > &A)
Calculates the deviator of a matrix (rank-2 tensor)
tensor< double, 2, 2 > mat2
statically sized 2x2 matrix of doubles
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,...
constexpr SERAC_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.
constexpr SERAC_HOST_DEVICE auto & operator-=(dual< gradient_type > &a, const dual< gradient_type > &b)
compound assignment (-) for dual numbers
constexpr SERAC_HOST_DEVICE auto tensor_with_shape(std::integer_sequence< int, n... >)
Creates a tensor given the dimensions in a std::integer_sequence.
constexpr SERAC_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
bool isnan(const zero &)
This is an overloaded member function, provided for convenience. It differs from the above function o...
constexpr SERAC_HOST_DEVICE auto sym(const tensor< T, n, n > &A)
Returns the symmetric part of a square matrix.
constexpr SERAC_HOST_DEVICE auto antisym(const tensor< T, n, n > &A)
Returns the antisymmetric part of a square matrix.
SERAC_HOST_DEVICE auto abs(dual< gradient_type > x)
Implementation of absolute value function for dual numbers.
constexpr SERAC_HOST_DEVICE auto inner(zero, double)
constexpr SERAC_HOST_DEVICE auto outer(const tensor< S, m > &A, const tensor< T, n > &B)
SERAC_HOST_DEVICE auto acos(dual< gradient_type > a)
implementation of acos for dual numbers
constexpr SERAC_HOST_DEVICE int dimension(const tensor< T, n... > &)
a function for querying the ith dimension of a tensor
auto & operator<<(std::ostream &out, zero)
Write a zero out to an output stream.
constexpr auto double_dot(const tensor< S, m, n > &A, const tensor< T, m, n > &B)
auto cross(const tensor< S, 3 > &u, const tensor< T, 3 > &v)
compute the (right handed) cross product of two 3-vectors
constexpr SERAC_HOST_DEVICE auto operator*(const tensor< S, m... > &A, const tensor< T, n... > &B)
this is a shorthand for dot(A, B)
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 & operator-=(tensor< T, n... > &A, zero)
compound assignment (-) between a tensor and zero (no-op)
constexpr SERAC_HOST_DEVICE auto & operator+=(dual< gradient_type > &a, const dual< gradient_type > &b)
compound assignment (+) for dual numbers
constexpr SERAC_HOST_DEVICE auto tr(const tensor< T, n, n > &A)
Returns the trace of a square matrix.
SERAC_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
tensor< double, 2 > vec2
statically sized vector of 2 doubles
constexpr SERAC_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 SERAC_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...
tensor< double, 3, 3 > mat3
statically sized 3x3 matrix of doubles
tensor(const T(&data)[n1][n2]) -> tensor< T, n1, n2 >
class template argument deduction guide for type tensor.
SERAC_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.
SERAC_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...
constexpr SERAC_HOST_DEVICE auto operator/(const dual< gradient_type > &a, double b)
division of a dual number by a non-dual number
constexpr SERAC_HOST_DEVICE auto solve_lower_triangular(const tensor< T, n, n > &L, const tensor< T, n, m... > &b)
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 SERAC_HOST_DEVICE auto norm(zero)
overload of Frobenius norm for zero type
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.
SERAC_HOST_DEVICE auto operator=(T)
anything assigned to zero does not change its value and returns zero
SERAC_HOST_DEVICE auto operator()(T...) const
zero can be accessed like a multidimensional array