All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
scalarMatrices.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::scalarMatrices
26 
27 Description
28  Scalar matrices
29 
30  LUDecompose for scalarSymmetricSquareMatrix implements the Cholesky
31  decomposition method from JAMA, a public-domain library developed at NIST,
32  available at http://math.nist.gov/tnt/index.html
33 
34 SourceFiles
35  scalarMatrices.C
36  scalarMatricesTemplates.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef scalarMatrices_H
41 #define scalarMatrices_H
42 
43 #include "RectangularMatrix.H"
44 #include "SquareMatrix.H"
45 #include "SymmetricSquareMatrix.H"
46 #include "DiagonalMatrix.H"
47 #include "scalarField.H"
48 #include "labelList.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
59 
60 //- Solve the matrix using Gaussian elimination with pivoting,
61 // returning the solution in the source
62 template<class Type>
63 void solve(scalarSquareMatrix& matrix, List<Type>& source);
64 
65 //- Solve the matrix using Gaussian elimination with pivoting
66 // and return the solution
67 template<class Type>
68 void solve
69 (
70  List<Type>& psi,
71  const scalarSquareMatrix& matrix,
72  const List<Type>& source
73 );
74 
75 //- LU decompose the matrix with pivoting
76 void LUDecompose
77 (
78  scalarSquareMatrix& matrix,
79  labelList& pivotIndices
80 );
81 
82 //- LU decompose the matrix with pivoting.
83 // sign is -1 for odd number of row interchanges and 1 for even number.
84 void LUDecompose
85 (
86  scalarSquareMatrix& matrix,
87  labelList& pivotIndices,
88  label& sign
89 );
90 
91 //- LU decompose the matrix into a lower (L) and upper (U) part. U = L.T()
92 void LUDecompose(scalarSymmetricSquareMatrix& matrix);
93 
94 //- LU back-substitution with given source, returning the solution
95 // in the source
96 template<class Type>
98 (
99  const scalarSquareMatrix& luMmatrix,
100  const labelList& pivotIndices,
101  List<Type>& source
102 );
103 
104 //- LU back-substitution with given source, returning the solution
105 // in the source. Specialised for symmetric square matrices that have been
106 // decomposed into LU where U = L.T() as pivoting is not required
107 template<class Type>
108 void LUBacksubstitute
109 (
110  const scalarSymmetricSquareMatrix& luMmatrix,
111  List<Type>& source
112 );
113 
114 //- Solve the matrix using LU decomposition with pivoting
115 // returning the LU form of the matrix and the solution in the source
116 template<class Type>
117 void LUsolve(scalarSquareMatrix& matrix, List<Type>& source);
118 
119 //- Solve the matrix using LU decomposition returning the LU form of the matrix
120 // and the solution in the source, where U = L.T()
121 template<class Type>
122 void LUsolve(scalarSymmetricSquareMatrix& matrix, List<Type>& source);
123 
124 template<class Form, class Type>
125 void multiply
126 (
127  Matrix<Form, Type>& answer, // value changed in return
128  const Matrix<Form, Type>& A,
129  const Matrix<Form, Type>& B
130 );
131 
132 void multiply
133 (
134  scalarRectangularMatrix& answer, // value changed in return
135  const scalarRectangularMatrix& A,
136  const scalarRectangularMatrix& B,
137  const scalarRectangularMatrix& C
138 );
139 
140 void multiply
141 (
142  scalarRectangularMatrix& answer, // value changed in return
143  const scalarRectangularMatrix& A,
144  const DiagonalMatrix<scalar>& B,
145  const scalarRectangularMatrix& C
146 );
147 
148 void multiply
149 (
150  scalarSquareMatrix& answer, // value changed in return
151  const scalarSquareMatrix& A,
152  const DiagonalMatrix<scalar>& B,
153  const scalarSquareMatrix& C
154 );
155 
156 //- Return the inverse of matrix A using SVD
157 scalarRectangularMatrix SVDinv
158 (
159  const scalarRectangularMatrix& A,
160  scalar minCondition = 0
161 );
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #ifdef NoRepository
171  #include "scalarMatricesTemplates.C"
172 #endif
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 #endif
177 
178 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
Graphite solid properties.
Definition: C.H:48
DiagonalMatrix< scalar > scalarDiagonalMatrix
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
A templated (m x n) matrix of objects of <T>.
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
A templated 2D square symmetric matrix of objects of <T>, where the n x n matrix dimension is known a...
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
const volScalarField & psi
SquareMatrix< scalar > scalarSquareMatrix
Namespace for OpenFOAM.
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
RectangularMatrix< scalar > scalarRectangularMatrix