Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::scalarMatrices
28
29Description
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
36SourceFiles
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"
48#include "DiagonalMatrix.H"
49#include "scalarField.H"
50#include "labelList.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54namespace Foam
56
61
62//- Solve the matrix using Gaussian elimination with pivoting,
63//- returning the solution in the source
64template<class Type>
65void solve(scalarSquareMatrix& matrix, List<Type>& source);
66
67//- Solve the matrix using Gaussian elimination with pivoting
68//- and return the solution
69template<class Type>
70void solve
71(
73 const scalarSquareMatrix& matrix,
74 const List<Type>& source
75);
76
77//- LU decompose the matrix with pivoting.
78void LUDecompose
79(
80 scalarSquareMatrix& matrix,
82 labelList& pivotIndices
83);
84
85//- LU decompose the matrix with pivoting.
86void 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
100template<class Type>
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
111template<class Type>
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
120template<class Type>
121void 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()
125template<class Type>
126void LUsolve(scalarSymmetricSquareMatrix& matrix, List<Type>& source);
127
128template<class Form, class Type>
129void 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
136void multiply
137(
138 scalarRectangularMatrix& answer, // value changed in return
142);
143
144void multiply
145(
146 scalarRectangularMatrix& answer, // value changed in return
150);
151
152void multiply
153(
154 scalarSquareMatrix& answer, // value changed in return
155 const scalarSquareMatrix& A,
157 const scalarSquareMatrix& C
158);
159
160//- Return the inverse of matrix A using SVD
162(
164 scalar minCondition = 0
165);
166
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170} // End namespace Foam
171
172// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173
174#ifdef NoRepository
176#endif
177
178// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179
180#endif
181
182// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
Definition C.H:49
A templated (N x N) diagonal matrix of objects of <Type>, effectively containing N elements,...
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 (m x n) matrix of objects of <T>. The layout is (mRows x nCols) - row-major order:
Definition Matrix.H:77
A templated (M x N) rectangular matrix of objects of <Type>, containing M*N elements,...
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
const volScalarField & psi
Namespace for OpenFOAM.
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
RectangularMatrix< scalar > scalarRectangularMatrix
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
SquareMatrix< scalar > scalarSquareMatrix
DiagonalMatrix< scalar > scalarDiagonalMatrix
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 multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
CEqn solve()