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{}
36
37
38template<class Type>
40{
41 decompose(mat);
42}
43
44
45// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
46
47template<class Type>
49{
50 SquareMatrix<Type>& LLT = *this;
51
52 // Initialize the LLT decomposition matrix to M
53 LLT = mat;
54
55 const label m = LLT.m();
56
57 for (label i = 0; i < m; ++i)
58 {
59 for (label j = 0; j < m; ++j)
60 {
61 if (j > i)
62 {
63 LLT(i, j) = Zero;
64 continue;
65 }
66
67 Type sum = LLT(i, j);
68
69 for (label k = 0; k < j; ++k)
70 {
71 sum -= LLT(i, k)*LLT(j, k);
72 }
73
74 if (i > j)
75 {
76 LLT(i, j) = sum/LLT(j, j);
77 }
78 else if (sum > 0)
79 {
80 LLT(i, i) = sqrt(sum);
81 }
82 else
83 {
85 << "Cholesky decomposition failed, "
86 "matrix is not symmetric positive definite"
87 << abort(FatalError);
88 }
89 }
90 }
91}
92
93
94template<class Type>
95template<template<typename> class ListContainer>
97(
99 const ListContainer<Type>& source
100) const
101{
102 // If x and source are different, copy initialize x = source
103 if (&x != &source)
104 {
105 x = source;
106 }
107
108 const SquareMatrix<Type>& LLT = *this;
109 const label m = LLT.m();
110
111 for (label i = 0; i < m; ++i)
112 {
113 Type sum = source[i];
114
115 for (label j = 0; j < i; ++j)
116 {
117 sum = sum - LLT(i, j)*x[j];
118 }
119
120 x[i] = sum/LLT(i, i);
121 }
122
123 for (label i = m - 1; i >= 0; --i)
124 {
125 Type sum = x[i];
126
127 for (label j = i + 1; j < m; ++j)
128 {
129 sum = sum - LLT(j, i)*x[j];
130 }
131
132 x[i] = sum/LLT(i, i);
133 }
134}
135
136
137template<class Type>
139(
140 List<Type>& x,
141 const UList<Type>& source
142) const
143{
144 solveImpl(x, source);
145}
146
147
148template<class Type>
149template<class Addr>
151(
152 List<Type>& x,
153 const IndirectListBase<Type, Addr>& source
154) const
155{
156 solveImpl(x, source);
157}
158
159template<class Type>
161(
162 const UList<Type>& source
163) const
164{
165 auto tresult(tmp<Field<Type>>::New(source.size()));
166
167 solve(tresult.ref(), source);
168
169 return tresult;
170}
171
172
173template<class Type>
174template<class Addr>
176(
177 const IndirectListBase<Type, Addr>& source
178) const
179{
180 auto tresult(tmp<Field<Type>>::New(source.size()));
181
182 solve(tresult.ref(), source);
183
184 return tresult;
185}
186
187
188// ************************************************************************* //
label k
Generic templated field type.
Definition: Field.H:82
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.
Templated class to perform the Cholesky decomposition on a symmetric positive-definite matrix.
Definition: LLTMatrix.H:60
LLTMatrix()
Construct null.
Definition: LLTMatrix.C:34
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
label m() const noexcept
The number of rows.
Definition: MatrixI.H:96
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
Definition: SquareMatrix.H:66
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:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A class for managing temporary objects.
Definition: tmp.H:65
bool decompose() const noexcept
Query the decompose flag (normally off)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
CEqn solve()