Loading...
Searching...
No Matches
SVD.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "SVD.H"
30#include "scalarList.H"
31#include "scalarMatrices.H"
32#include "ListOps.H"
33
34#include "tensor.H"
35#include "diagTensor.H"
37// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41
42// Same as implementation in scalarMatrices, but for tensor/diagTensor/tensor
44(
45 const tensor& A,
46 const diagTensor& B,
47 const tensor& C
48)
49{
50 constexpr direction size = 3;
51
52 tensor ans(Foam::zero{});
53
54 for (direction i = 0; i < size; ++i)
55 {
56 for (direction g = 0; g < size; ++g)
57 {
58 for (direction l = 0; l < size; ++l)
59 {
60 ans(i, g) += C(l, g)*A(i, l)*B[l];
61 }
62 }
63 }
64
65 return ans;
66}
67
68
69template<class MatrixType, class DiagMatrixType, class WorkArrayType>
70static std::pair<label, bool> SVDcomp
71(
72 // input
73 const MatrixType& A,
74
75 // input
76 const scalar minCondition,
77
78 // output
79 MatrixType& U_,
80
81 // output
82 MatrixType& V_,
83
84 // output
85 DiagMatrixType& S_,
86
87 // scratch space
88 WorkArrayType& rv1
89)
90{
91 label nZeros_(0);
92 bool converged_(true);
93
94 // SVDcomp to find U_, V_ and S_ - the singular values
95
96 U_ = A;
97
98 const label Un = U_.n();
99 const label Um = U_.m();
100
101 scalar g = 0;
102 scalar scale = 0;
103 scalar s = 0;
104 scalar anorm = 0;
105 label l = 0;
106
107
108 const auto sign = [](scalar a, scalar b) -> scalar
109 {
110 //return b >= 0 ? (a >= 0 ? a : -a) : (a >= 0 ? -a : a);
111 return b >= 0 ? a : -a;
112 };
113
114
115 for (label i=0; i<Un; i++)
116 {
117 l = i + 2;
118 rv1[i] = scale*g;
119 g = s = scale = 0;
120
121 if (i < Um)
122 {
123 for (label k=i; k<Um; k++)
124 {
125 scale += mag(U_(k, i));
126 }
127
128 if (scale != 0)
129 {
130 for (label k=i; k<Um; k++)
131 {
132 U_(k, i) /= scale;
133 s += U_(k, i)*U_(k, i);
134 }
135
136 scalar f = U_(i, i);
137 g = -sign(sqrt(s), f);
138 scalar h = f*g - s;
139 U_(i, i) = f - g;
140
141 for (label j=l-1; j<Un; j++)
142 {
143 s = 0;
144 for (label k=i; k<Um; k++)
145 {
146 s += U_(k, i)*U_(k, j);
147 }
148
149 f = s/h;
150 for (label k=i; k<Um; k++)
151 {
152 U_(k, j) += f*U_(k, i);
153 }
154 }
155
156 for (label k=i; k<Um; k++)
157 {
158 U_(k, i) *= scale;
159 }
160 }
161 }
162
163 S_[i] = scale*g;
164
165 g = s = scale = 0;
166
167 if (i+1 <= Um && i != Un)
168 {
169 for (label k=l-1; k<Un; k++)
170 {
171 scale += mag(U_(i, k));
172 }
173
174 if (scale != 0)
175 {
176 for (label k=l-1; k<Un; k++)
177 {
178 U_(i, k) /= scale;
179 s += U_(i, k)*U_(i, k);
180 }
181
182 scalar f = U_(i, l-1);
183 g = -sign(sqrt(s), f);
184 scalar h = f*g - s;
185 U_(i, l-1) = f - g;
186
187 for (label k=l-1; k<Un; k++)
188 {
189 rv1[k] = U_(i, k)/h;
190 }
191
192 for (label j=l-1; j<Um; j++)
193 {
194 s = 0;
195 for (label k=l-1; k<Un; k++)
196 {
197 s += U_(j, k)*U_(i, k);
198 }
199
200 for (label k=l-1; k<Un; k++)
201 {
202 U_(j, k) += s*rv1[k];
203 }
204 }
205 for (label k=l-1; k<Un; k++)
206 {
207 U_(i, k) *= scale;
208 }
209 }
210 }
211
212 anorm = max(anorm, mag(S_[i]) + mag(rv1[i]));
213 }
214
215 anorm *= SMALL;
216
217 for (label i=Un-1; i >= 0; i--)
218 {
219 if (i < Un-1)
220 {
221 if (g*U_(i, l) != 0)
222 {
223 for (label j=l; j<Un; j++)
224 {
225 V_(j, i) = (U_(i, j)/U_(i, l))/g;
226 }
227
228 for (label j=l; j<Un; j++)
229 {
230 s = 0;
231 for (label k=l; k<Un; k++)
232 {
233 s += U_(i, k)*V_(k, j);
234 }
235
236 for (label k=l; k<Un; k++)
237 {
238 V_(k, j) += s*V_(k, i);
239 }
240 }
241 }
242
243 for (label j=l; j<Un;j++)
244 {
245 V_(i, j) = V_(j, i) = 0;
246 }
247 }
248
249 V_(i, i) = 1;
250 g = rv1[i];
251 l = i;
252 }
253
254 for (label i=min(Un, Um) - 1; i>=0; i--)
255 {
256 l = i+1;
257 g = S_[i];
258
259 for (label j=l; j<Un; j++)
260 {
261 U_(i, j) = 0;
262 }
263
264 if (g*U_(i, i) != 0)
265 {
266 g = 1.0/g;
267
268 for (label j=l; j<Un; j++)
269 {
270 s = 0;
271 for (label k=l; k<Um; k++)
272 {
273 s += U_(k, i)*U_(k, j);
274 }
275
276 scalar f = (s/U_(i, i))*g;
277
278 for (label k=i; k<Um; k++)
279 {
280 U_(k, j) += f*U_(k, i);
281 }
282 }
283
284 for (label j=i; j<Um; j++)
285 {
286 U_(j, i) *= g;
287 }
288 }
289 else
290 {
291 for (label j=i; j<Um; j++)
292 {
293 U_(j, i) = 0;
294 }
295 }
296
297 ++U_(i, i);
298 }
299
300 for (label k=Un-1; k >= 0; k--)
301 {
302 for (label its = 0; its < 30; its++)
303 {
304 bool flag = true;
305
306 label mn;
307 for (l = k; l >= 0; l--)
308 {
309 mn = l - 1;
310
311 if (l == 0 || mag(rv1[l]) <= anorm)
312 {
313 flag = false;
314 break;
315 }
316
317 if (mag(S_[mn]) <= anorm)
318 {
319 break;
320 }
321 }
322
323 if (flag)
324 {
325 scalar c = 0;
326 s = 1;
327 for (label i=l; i<k+1; i++)
328 {
329 scalar f = s*rv1[i];
330 rv1[i] = c*rv1[i];
331
332 if (mag(f) <= anorm)
333 {
334 break;
335 }
336
337 g = S_[i];
338 scalar h = sqrtSumSqr(f, g);
339 S_[i] = h;
340 h = 1.0/h;
341 c = g*h;
342 s = -f*h;
343
344 for (label j=0; j<Um; j++)
345 {
346 scalar y = U_(j, mn);
347 scalar z = U_(j, i);
348 U_(j, mn) = y*c + z*s;
349 U_(j, i) = z*c - y*s;
350 }
351 }
352 }
353
354 scalar z = S_[k];
355
356 if (l == k)
357 {
358 if (z < 0)
359 {
360 S_[k] = -z;
361 for (label j=0; j<Un; j++)
362 {
363 V_(j, k) = -V_(j, k);
364 }
365 }
366 break;
367 }
368
369 if (its == 29)
370 {
371 converged_ = false;
372 }
373
374 scalar x = S_[l];
375 mn = k-1;
376 scalar y = S_[mn];
377 g = rv1[mn];
378 scalar h = rv1[k];
379 scalar f = ((y - z)*(y + z) + (g - h)*(g + h))/(2*h*y);
380 g = sqrtSumSqr(f, scalar(1));
381 f = ((x - z)*(x + z) + h*((y/(f + sign(g, f))) - h))/x;
382 scalar c = 1;
383 s = 1;
384
385 for (label j=l; j <= mn; j++)
386 {
387 label i = j + 1;
388 g = rv1[i];
389 y = S_[i];
390 h = s*g;
391 g = c*g;
392 scalar z = sqrtSumSqr(f, h);
393 rv1[j] = z;
394 c = f/z;
395 s = h/z;
396 f = x*c + g*s;
397 g = g*c - x*s;
398 h = y*s;
399 y *= c;
400
401 for (label jj = 0; jj < Un; jj++)
402 {
403 x = V_(jj, j);
404 z = V_(jj, i);
405 V_(jj, j) = x*c + z*s;
406 V_(jj, i) = z*c - x*s;
407 }
408
409 z = sqrtSumSqr(f, h);
410 S_[j] = z;
411 if (z)
412 {
413 z = 1.0/z;
414 c = f*z;
415 s = h*z;
416 }
417 f = c*g + s*y;
418 x = c*y - s*g;
419
420 for (label jj=0; jj < Um; jj++)
421 {
422 y = U_(jj, j);
423 z = U_(jj, i);
424 U_(jj, j) = y*c + z*s;
425 U_(jj, i) = z*c - y*s;
426 }
427 }
428 rv1[l] = 0;
429 rv1[k] = f;
430 S_[k] = x;
431 }
432 }
433
434 // Zero singular values that are less than minCondition*maxS
435 const scalar minS = minCondition*S_[findMax(S_)];
436 for (auto& val : S_)
437 {
438 if (val <= minS)
439 {
440 val = 0;
441 ++nZeros_;
442 }
443 }
444
445 return { nZeros_, converged_ };
447
448} // End namespace Foam
449
450
451// * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * //
452
453Foam::SVD::SVD(const scalarRectangularMatrix& A, const scalar minCondition)
454:
455 U_(),
456 V_(A.n(), A.n()),
457 S_(A.n()),
458 converged_(true),
459 nZeros_(0)
460{
461 scalarList rv1(A.n());
462
463 // SVDcomp to find U_, V_ and S_ - the singular values
464
465 auto status = SVDcomp(A, minCondition, U_, V_, S_, rv1);
466
467 nZeros_ = status.first;
468 converged_ = status.second;
469}
470
471
472// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
473
474Foam::scalar Foam::SVD::minNonZeroS() const
475{
476 scalar minS = S_[0];
477 for (label i=1; i<S_.size(); i++)
478 {
479 scalar s = S_[i];
480 if (s > VSMALL && s < minS) minS = s;
481 }
482 return minS;
483}
484
485
487{
489 multiply(VSinvUt, V_, Foam::inv(S_), U_.T());
490 return VSinvUt;
491}
492
493
494// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
495
497(
499 const scalar minCondition
500)
502 SVD svd(A, minCondition);
503 return svd.VSinvUt();
504}
505
506
508Foam::SVD::pinv(const Tensor<scalar>& A, const scalar minCondition)
509{
510 // SVDcomp to find U_, V_ and S_ - the singular values
511
512 tensor U_;
513 tensor V_;
514 diagTensor S_;
516
517 SVDcomp(A, minCondition, U_, V_, S_, rv1);
518
519 return do_multiply(V_, Foam::inv(S_), U_.T());
520}
521
522
523// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
scalar y
label k
Various functions to operate on Lists.
label n
const uniformDimensionedVectorField & g
Graphite solid properties.
Definition C.H:49
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
scalarRectangularMatrix VSinvUt() const
Return the matrix product V S^(-1) U^T (the pseudo inverse).
Definition SVD.C:479
scalar minNonZeroS() const
Return the minimum non-zero singular value.
Definition SVD.C:467
SVD(const SVD &)=delete
No copy construct.
static scalarRectangularMatrix pinv(const scalarRectangularMatrix &A, const scalar minCondition=0)
Return the pseudo inverse of the given matrix.
Definition SVD.C:490
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition Tensor.H:60
Tensor< Cmpt > T() const
Return non-Hermitian transpose.
Definition TensorI.H:419
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DiagTensor< scalar > diagTensor
DiagTensor of scalars, i.e. DiagTensor<scalar>.
Definition diagTensor.H:49
dimensionedScalar sign(const dimensionedScalar &ds)
RectangularMatrix< scalar > scalarRectangularMatrix
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static tensor do_multiply(const tensor &A, const diagTensor &B, const tensor &C)
Definition SVD.C:37
uint8_t direction
Definition direction.H:49
static std::pair< label, bool > SVDcomp(const MatrixType &A, const scalar minCondition, MatrixType &U_, MatrixType &V_, DiagMatrixType &S_, WorkArrayType &rv1)
Definition SVD.C:64
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
Scalar sqrtSumSqr(const Scalar a, const Scalar b)
Definition scalarImpl.H:504
labelList f(nPoints)
volScalarField & h
volScalarField & b