Loading...
Searching...
No Matches
LLTMatrix.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) 2016 OpenFOAM Foundation
9 Copyright (C) 2019 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 "LLTMatrix.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35{
37}
38
39
40// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
41
42template<class Type>
44{
45 SquareMatrix<Type>& LLT = *this;
46
47 // Initialize the LLT decomposition matrix to M
48 LLT = mat;
49
50 const label m = LLT.m();
51
52 for (label i = 0; i < m; ++i)
53 {
54 for (label j = 0; j < m; ++j)
55 {
56 if (j > i)
57 {
58 LLT(i, j) = Zero;
59 continue;
60 }
61
62 Type sum = LLT(i, j);
63
64 for (label k = 0; k < j; ++k)
65 {
66 sum -= LLT(i, k)*LLT(j, k);
67 }
68
69 if (i > j)
70 {
71 LLT(i, j) = sum/LLT(j, j);
72 }
73 else if (sum > 0)
74 {
75 LLT(i, i) = sqrt(sum);
76 }
77 else
78 {
80 << "Cholesky decomposition failed, "
81 "matrix is not symmetric positive definite"
82 << abort(FatalError);
83 }
84 }
85 }
86}
87
88
89template<class Type>
90template<template<typename> class ListContainer>
91void Foam::LLTMatrix<Type>::solveImpl
92(
93 List<Type>& x,
94 const ListContainer<Type>& source
95) const
96{
97 // If x and source are different, copy initialize x = source
98 if (&x != &source)
99 {
100 x = source;
101 }
102
103 const SquareMatrix<Type>& LLT = *this;
104 const label m = LLT.m();
105
106 for (label i = 0; i < m; ++i)
107 {
108 Type sum = source[i];
109
110 for (label j = 0; j < i; ++j)
111 {
112 sum = sum - LLT(i, j)*x[j];
113 }
114
115 x[i] = sum/LLT(i, i);
116 }
117
118 for (label i = m - 1; i >= 0; --i)
119 {
120 Type sum = x[i];
121
122 for (label j = i + 1; j < m; ++j)
123 {
124 sum = sum - LLT(j, i)*x[j];
125 }
127 x[i] = sum/LLT(i, i);
128 }
129}
130
131
132template<class Type>
134(
135 List<Type>& x,
136 const UList<Type>& source
137) const
139 solveImpl(x, source);
140}
141
142
143template<class Type>
144template<class Addr>
146(
147 List<Type>& x,
149) const
150{
151 solveImpl(x, source);
152}
153
154template<class Type>
156(
157 const UList<Type>& source
158) const
159{
160 auto tresult(tmp<Field<Type>>::New(source.size()));
161
162 solve(tresult.ref(), source);
164 return tresult;
165}
166
167
168template<class Type>
169template<class Addr>
171(
172 const IndirectListBase<Type, Addr>& source
173) const
174{
175 auto tresult(tmp<Field<Type>>::New(source.size()));
176
177 solve(tresult.ref(), source);
178
179 return tresult;
180}
181
182
183// ************************************************************************* //
label k
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Base for lists with indirect addressing, templated on the list contents type and the addressing type....
label size() const noexcept
The number of elements in the list.
LLTMatrix()=default
Default construct.
void solve(List< Type > &x, const UList< Type > &source) const
Solve the linear system with the given source and return the solution in the argument x.
Definition LLTMatrix.C:127
void decompose(const SquareMatrix< Type > &mat)
Copy matrix and perform Cholesky decomposition.
Definition LLTMatrix.C:36
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
SquareMatrix()=default
Default construct.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A class for managing temporary objects.
Definition tmp.H:75
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
CEqn solve()