LUscalarMatrix.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::LUscalarMatrix
28 
29 Description
30  Class to perform the LU decomposition on a symmetric matrix.
31 
32 SourceFiles
33  LUscalarMatrix.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef Foam_LUscalarMatrix_H
38 #define Foam_LUscalarMatrix_H
39 
40 #include "scalarMatrices.H"
41 #include "labelList.H"
42 #include "FieldField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward Declarations
51 class lduMatrix;
52 class procLduMatrix;
53 
54 /*---------------------------------------------------------------------------*\
55  Class LUscalarMatrix Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class LUscalarMatrix
59 :
60  public scalarSquareMatrix
61 {
62  // Private Data
63 
64  //- Communicator to use
65  const label comm_;
66 
67  //- Processor matrix offsets
68  labelList procOffsets_;
69 
70  //- The pivot indices used in the LU decomposition
71  labelList pivotIndices_;
72 
73 
74  // Private Member Functions
75 
76  //- Convert the given lduMatrix into this LUscalarMatrix
77  void convert
78  (
79  const lduMatrix& ldum,
80  const FieldField<Field, scalar>& interfaceCoeffs,
81  const lduInterfaceFieldPtrsList& interfaces
82  );
83 
84  //- Convert the given list of procLduMatrix into this LUscalarMatrix
85  //- on the master processor
86  void convert(const PtrList<procLduMatrix>& lduMatrices);
87 
88 
89  //- Print the ratio of the mag-sum of the off-diagonal coefficients
90  //- to the mag-diagonal
91  void printDiagonalDominance() const;
92 
93 
94 public:
95 
96  // Declare name of the class and its debug switch
97  ClassName("LUscalarMatrix");
98 
99 
100  // Constructors
101 
102  //- Default construct
104 
105  //- Construct from and perform LU decomposition of the given matrix
106  explicit LUscalarMatrix(const scalarSquareMatrix& mat);
107 
108  //- Construct from lduMatrix and perform LU decomposition.
109  //- In parallel it assembles the matrix on the master.
111  (
112  const lduMatrix& ldum,
113  const FieldField<Field, scalar>& interfaceCoeffs,
114  const lduInterfaceFieldPtrsList& interfaces
115  );
116 
117 
118  // Member Functions
119 
120  //- Perform the LU decomposition of the matrix
121  void decompose(const scalarSquareMatrix& mat);
122 
123  //- Solve the linear system with the given source
124  // and returning the solution in the Field argument x.
125  // This function may be called with the same field for x and source.
126  template<class Type>
127  void solve(List<Type>& x, const UList<Type>& source) const;
128 
129  //- Solve the linear system with the given source
130  // returning the solution
131  template<class Type>
132  tmp<Field<Type>> solve(const UList<Type>& source) const;
133 
134  //- Set M to the inverse of this square matrix
135  void inv(scalarSquareMatrix& M) const;
136 };
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 } // End namespace Foam
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 #ifdef NoRepository
146  #include "LUscalarMatrixTemplates.C"
147 #endif
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 #endif
152 
153 // ************************************************************************* //
void solve(List< Type > &x, const UList< Type > &source) const
Solve the linear system with the given source.
LUscalarMatrix() noexcept
Default construct.
Class to perform the LU decomposition on a symmetric matrix.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
const direction noexcept
Definition: Scalar.H:258
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:80
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
void decompose(const scalarSquareMatrix &mat)
Perform the LU decomposition of the matrix.
ClassName("LUscalarMatrix")
A class for managing temporary objects.
Definition: HashPtrTable.H:50
#define M(I)
Namespace for OpenFOAM.
void inv(scalarSquareMatrix &M) const
Set M to the inverse of this square matrix.