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,
82  labelList& pivotIndices
83 );
84 
85 //- LU decompose the matrix with pivoting.
86 void LUDecompose
87 (
88  scalarSquareMatrix& matrix,
90  labelList& pivotIndices,
92  label& sign
93 );
94 
95 //- LU decompose the matrix into a lower (L) and upper (U) part. U = L.T()
97 
98 //- LU back-substitution with given source, returning the solution
99 //- in the source
100 template<class Type>
101 void LUBacksubstitute
102 (
103  const scalarSquareMatrix& luMmatrix,
104  const labelList& pivotIndices,
105  List<Type>& source
106 );
107 
108 //- LU back-substitution with given source, returning the solution
109 //- in the source. Specialised for symmetric square matrices that have been
110 //- decomposed into LU where U = L.T() as pivoting is not required
111 template<class Type>
112 void LUBacksubstitute
113 (
114  const scalarSymmetricSquareMatrix& luMmatrix,
115  List<Type>& source
116 );
117 
118 //- Solve the matrix using LU decomposition with pivoting
119 //- returning the LU form of the matrix and the solution in the source
120 template<class Type>
121 void LUsolve(scalarSquareMatrix& matrix, List<Type>& source);
122 
123 //- Solve the matrix using LU decomposition returning the LU form of the matrix
124 //- and the solution in the source, where U = L.T()
125 template<class Type>
126 void LUsolve(scalarSymmetricSquareMatrix& matrix, List<Type>& source);
127 
128 template<class Form, class Type>
129 void multiply
130 (
131  Matrix<Form, Type>& answer, // value changed in return
132  const Matrix<Form, Type>& A,
133  const Matrix<Form, Type>& B
134 );
135 
136 void multiply
137 (
138  scalarRectangularMatrix& answer, // value changed in return
139  const scalarRectangularMatrix& A,
140  const scalarRectangularMatrix& B,
142 );
143 
144 void multiply
145 (
146  scalarRectangularMatrix& answer, // value changed in return
147  const scalarRectangularMatrix& A,
148  const DiagonalMatrix<scalar>& B,
150 );
151 
152 void multiply
153 (
154  scalarSquareMatrix& answer, // value changed in return
155  const scalarSquareMatrix& A,
156  const DiagonalMatrix<scalar>& B,
157  const scalarSquareMatrix& C
158 );
159 
160 //- Return the inverse of matrix A using SVD
162 (
163  const scalarRectangularMatrix& A,
164  scalar minCondition = 0
165 );
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 } // End namespace Foam
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 #ifdef NoRepository
175  #include "scalarMatricesTemplates.C"
176 #endif
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 #endif
181 
182 // ************************************************************************* //
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:46
DiagonalMatrix< scalar > scalarDiagonalMatrix
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
A templated (m x n) matrix of objects of <T>. The layout is (mRows x nCols) - row-major order: ...
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix...
const volScalarField & psi
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
SquareMatrix< scalar > scalarSquareMatrix
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Namespace for OpenFOAM.
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
RectangularMatrix< scalar > scalarRectangularMatrix