Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/db1/math_8h_source.html
DOLFINx 0.6.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
math.h
1// Copyright (C) 2021 Igor Baratta
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <array>
10#include <basix/mdspan.hpp>
11#include <cmath>
12#include <string>
13#include <type_traits>
14
15namespace dolfinx::math
16{
17
22template <typename U, typename V>
23std::array<typename U::value_type, 3> cross(const U& u, const V& v)
24{
25 assert(u.size() == 3);
26 assert(v.size() == 3);
27 return {u[1] * v[2] - u[2] * v[1], u[2] * v[0] - u[0] * v[2],
28 u[0] * v[1] - u[1] * v[0]};
29}
30
33template <typename T>
34T difference_of_products(T a, T b, T c, T d) noexcept
35{
36 T w = b * c;
37 T err = std::fma(-b, c, w);
38 T diff = std::fma(a, d, -w);
39 return diff + err;
40}
41
48template <typename T>
49auto det(const T* A, std::array<std::size_t, 2> shape)
50{
51 assert(shape[0] == shape[1]);
52
53 // const int nrows = shape[0];
54 switch (shape[0])
55 {
56 case 1:
57 return *A;
58 case 2:
59 /* A(0, 0), A(0, 1), A(1, 0), A(1, 1) */
60 return difference_of_products(A[0], A[1], A[2], A[3]);
61 case 3:
62 {
63 // Leibniz formula combined with Kahan’s method for accurate
64 // computation of 3 x 3 determinants
65 T w0 = difference_of_products(A[3 + 1], A[3 + 2], A[3 * 2 + 1],
66 A[2 * 3 + 2]);
67 T w1 = difference_of_products(A[3], A[3 + 2], A[3 * 2], A[3 * 2 + 2]);
68 T w2 = difference_of_products(A[3], A[3 + 1], A[3 * 2], A[3 * 2 + 1]);
69 T w3 = difference_of_products(A[0], A[1], w1, w0);
70 T w4 = std::fma(A[2], w2, w3);
71 return w4;
72 }
73 default:
74 throw std::runtime_error("math::det is not implemented for "
75 + std::to_string(A[0]) + "x" + std::to_string(A[1])
76 + " matrices.");
77 }
78}
79
84template <typename Matrix>
85auto det(Matrix A)
86{
87 static_assert(Matrix::rank() == 2, "Must be rank 2");
88 assert(A.extent(0) == A.extent(1));
89
90 using value_type = typename Matrix::value_type;
91 const int nrows = A.extent(0);
92 switch (nrows)
93 {
94 case 1:
95 return A(0, 0);
96 case 2:
97 return difference_of_products(A(0, 0), A(0, 1), A(1, 0), A(1, 1));
98 case 3:
99 {
100 // Leibniz formula combined with Kahan’s method for accurate
101 // computation of 3 x 3 determinants
102 value_type w0 = difference_of_products(A(1, 1), A(1, 2), A(2, 1), A(2, 2));
103 value_type w1 = difference_of_products(A(1, 0), A(1, 2), A(2, 0), A(2, 2));
104 value_type w2 = difference_of_products(A(1, 0), A(1, 1), A(2, 0), A(2, 1));
105 value_type w3 = difference_of_products(A(0, 0), A(0, 1), w1, w0);
106 value_type w4 = std::fma(A(0, 2), w2, w3);
107 return w4;
108 }
109 default:
110 throw std::runtime_error("math::det is not implemented for "
111 + std::to_string(A.extent(0)) + "x"
112 + std::to_string(A.extent(1)) + " matrices.");
113 }
114}
115
122template <typename U, typename V>
123void inv(U A, V B)
124{
125 static_assert(U::rank() == 2, "Must be rank 2");
126 static_assert(V::rank() == 2, "Must be rank 2");
127
128 using value_type = typename U::value_type;
129 const std::size_t nrows = A.extent(0);
130 switch (nrows)
131 {
132 case 1:
133 B(0, 0) = 1 / A(0, 0);
134 break;
135 case 2:
136 {
137 value_type idet = 1. / det(A);
138 B(0, 0) = idet * A(1, 1);
139 B(0, 1) = -idet * A(0, 1);
140 B(1, 0) = -idet * A(1, 0);
141 B(1, 1) = idet * A(0, 0);
142 break;
143 }
144 case 3:
145 {
146 value_type w0 = difference_of_products(A(1, 1), A(1, 2), A(2, 1), A(2, 2));
147 value_type w1 = difference_of_products(A(1, 0), A(1, 2), A(2, 0), A(2, 2));
148 value_type w2 = difference_of_products(A(1, 0), A(1, 1), A(2, 0), A(2, 1));
149 value_type w3 = difference_of_products(A(0, 0), A(0, 1), w1, w0);
150 value_type det = std::fma(A(0, 2), w2, w3);
151 assert(det != 0.);
152 value_type idet = 1 / det;
153
154 B(0, 0) = w0 * idet;
155 B(1, 0) = -w1 * idet;
156 B(2, 0) = w2 * idet;
157 B(0, 1) = difference_of_products(A(0, 2), A(0, 1), A(2, 2), A(2, 1)) * idet;
158 B(0, 2) = difference_of_products(A(0, 1), A(0, 2), A(1, 1), A(1, 2)) * idet;
159 B(1, 1) = difference_of_products(A(0, 0), A(0, 2), A(2, 0), A(2, 2)) * idet;
160 B(1, 2) = difference_of_products(A(1, 0), A(0, 0), A(1, 2), A(0, 2)) * idet;
161 B(2, 1) = difference_of_products(A(2, 0), A(0, 0), A(2, 1), A(0, 1)) * idet;
162 B(2, 2) = difference_of_products(A(0, 0), A(1, 0), A(0, 1), A(1, 1)) * idet;
163 break;
164 }
165 default:
166 throw std::runtime_error("math::inv is not implemented for "
167 + std::to_string(A.extent(0)) + "x"
168 + std::to_string(A.extent(1)) + " matrices.");
169 }
170}
171
178template <typename U, typename V, typename P>
179void dot(U A, V B, P C, bool transpose = false)
180{
181 static_assert(U::rank() == 2, "Must be rank 2");
182 static_assert(V::rank() == 2, "Must be rank 2");
183 static_assert(P::rank() == 2, "Must be rank 2");
184
185 if (transpose)
186 {
187 assert(A.extent(0) == B.extent(1));
188 for (std::size_t i = 0; i < A.extent(1); i++)
189 for (std::size_t j = 0; j < B.extent(0); j++)
190 for (std::size_t k = 0; k < A.extent(0); k++)
191 C(i, j) += A(k, i) * B(j, k);
192 }
193 else
194 {
195 assert(A.extent(1) == B.extent(0));
196 for (std::size_t i = 0; i < A.extent(0); i++)
197 for (std::size_t j = 0; j < B.extent(1); j++)
198 for (std::size_t k = 0; k < A.extent(1); k++)
199 C(i, j) += A(i, k) * B(k, j);
200 }
201}
202
209template <typename U, typename V>
210void pinv(U A, V P)
211{
212 static_assert(U::rank() == 2, "Must be rank 2");
213 static_assert(V::rank() == 2, "Must be rank 2");
214
215 assert(A.extent(0) > A.extent(1));
216 assert(P.extent(1) == A.extent(0));
217 assert(P.extent(0) == A.extent(1));
218 using T = typename U::value_type;
219 if (A.extent(1) == 2)
220 {
221 namespace stdex = std::experimental;
222 std::array<T, 6> ATb;
223 std::array<T, 4> ATAb, Invb;
224 stdex::mdspan<T, stdex::extents<std::size_t, 2, 3>> AT(ATb.data(), 2, 3);
225 stdex::mdspan<T, stdex::extents<std::size_t, 2, 2>> ATA(ATAb.data(), 2, 2);
226 stdex::mdspan<T, stdex::extents<std::size_t, 2, 2>> Inv(Invb.data(), 2, 2);
227
228 for (std::size_t i = 0; i < AT.extent(0); ++i)
229 for (std::size_t j = 0; j < AT.extent(1); ++j)
230 AT(i, j) = A(j, i);
231
232 std::fill(ATAb.begin(), ATAb.end(), 0.0);
233 for (std::size_t i = 0; i < P.extent(0); ++i)
234 for (std::size_t j = 0; j < P.extent(1); ++j)
235 P(i, j) = 0;
236
237 // pinv(A) = (A^T * A)^-1 * A^T
238 dot(AT, A, ATA);
239 inv(ATA, Inv);
240 dot(Inv, AT, P);
241 }
242 else if (A.extent(1) == 1)
243 {
244 T res = 0;
245 for (std::size_t i = 0; i < A.extent(0); ++i)
246 for (std::size_t j = 0; j < A.extent(1); ++j)
247 res += A(i, j) * A(i, j);
248
249 for (std::size_t i = 0; i < A.extent(0); ++i)
250 for (std::size_t j = 0; j < A.extent(1); ++j)
251 P(j, i) = (1 / res) * A(i, j);
252 }
253 else
254 {
255 throw std::runtime_error("math::pinv is not implemented for "
256 + std::to_string(A.extent(0)) + "x"
257 + std::to_string(A.extent(1)) + " matrices.");
258 }
259}
260
261} // namespace dolfinx::math