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