Loading...
Searching...
No Matches
faMatrixSolve.C
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) 2016-2017 Wikki Ltd
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Description
28 Finite-Area matrix basic solvers.
29
30\*---------------------------------------------------------------------------*/
32#include "PrecisionAdaptor.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const label patchi,
40 const label facei,
41 const direction cmpt,
42 const scalar value
43)
44{
45 internalCoeffs_[patchi][facei].component(cmpt) +=
46 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
47
48 boundaryCoeffs_[patchi][facei].component(cmpt) +=
49 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]*value;
50}
51
52
53template<class Type>
55(
56 const dictionary& solverControls
57)
58{
60 << "solving faMatrix<Type>"
61 << endl;
62
63 const int logLevel =
64 solverControls.getOrDefault<int>
65 (
66 "log",
68 );
69
70 auto& psi = psi_.constCast();
71
72 SolverPerformance<Type> solverPerfVec
73 (
74 "faMatrix<Type>::solve",
75 psi.name()
76 );
77
78 scalarField saveDiag(diag());
79
80 Field<Type> source(source_);
81 addBoundarySource(source);
82
83 // Note: make a copy of interfaces: no longer a reference
84 lduInterfaceFieldPtrsList interfaces =
85 psi_.boundaryField().scalarInterfaces();
86
87 for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
88 {
89 // copy field and source
90
91 scalarField psiCmpt(psi_.primitiveField().component(cmpt));
92 addBoundaryDiag(diag(), cmpt);
93
94 scalarField sourceCmpt(source.component(cmpt));
95
96 FieldField<Field, scalar> bouCoeffsCmpt
97 (
98 boundaryCoeffs_.component(cmpt)
99 );
100
101 FieldField<Field, scalar> intCoeffsCmpt
102 (
103 internalCoeffs_.component(cmpt)
104 );
105
106 // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
107 // bouCoeffsCmpt for the explicit part of the coupled boundary
108 // conditions
109 {
110 PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
112
113 const label startRequest = UPstream::nRequests();
114
115 initMatrixInterfaces
116 (
117 true,
118 bouCoeffsCmpt,
119 interfaces,
120 psiCmpt_ss(),
121 sourceCmpt_ss.ref(),
122 cmpt
123 );
124
125 updateMatrixInterfaces
126 (
127 true,
128 bouCoeffsCmpt,
129 interfaces,
130 psiCmpt_ss(),
131 sourceCmpt_ss.ref(),
132 cmpt,
133 startRequest
134 );
135 }
136
137 solverPerformance solverPerf;
138
139 // Solver call
140 solverPerf = lduMatrix::solver::New
141 (
142 psi_.name() + pTraits<Type>::componentNames[cmpt],
143 *this,
144 bouCoeffsCmpt,
145 intCoeffsCmpt,
146 interfaces,
147 solverControls
148 )->solve(psiCmpt, sourceCmpt, cmpt);
149
150 if (logLevel)
151 {
152 solverPerf.print(Info);
153 }
154
155 solverPerfVec.replace(cmpt, solverPerf);
156 solverPerfVec.solverName() = solverPerf.solverName();
157
158 psi.primitiveFieldRef().replace(cmpt, psiCmpt);
159 diag() = saveDiag;
160 }
161
162 psi.correctBoundaryConditions();
163
164 psi.mesh().data().setSolverPerformance(psi.name(), solverPerfVec);
166 return solverPerfVec;
167}
168
169
170template<class Type>
173 return solve(faMat_.solverDict());
174}
175
176
177template<class Type>
180 return this->solve(solverDict(name));
181}
182
183
184template<class Type>
187 return this->solve(solverDict());
188}
189
190
191template<class Type>
193{
194 auto tres = tmp<Field<Type>>::New(source_);
195 auto& res = tres().ref();
196
198
199 // Loop over field components
200 for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
201 {
202 scalarField psiCmpt(psi_.internalField().component(cmpt));
203
204 scalarField boundaryDiagCmpt(psi_.size(), Zero);
205 addBoundaryDiag(boundaryDiagCmpt, cmpt);
206
207 FieldField<Field, scalar> bouCoeffsCmpt
208 (
209 boundaryCoeffs_.component(cmpt)
210 );
211
212 res.replace
213 (
214 cmpt,
216 (
217 psiCmpt,
218 res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
219 bouCoeffsCmpt,
221 cmpt
222 )
223 );
224 }
225
226 return tres;
227}
228
229
230// ************************************************************************* //
A const Field/List wrapper with possible data conversion.
tmp< DimensionedField< cmptType, GeoMesh > > component(const direction d) const
Return a component field of the field.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
A non-const Field/List wrapper with possible data conversion.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void replace(const label cmpt, const SolverPerformance< typename pTraits< Type >::cmptType > &sp)
Replace component based on the minimal SolverPerformance.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
const word & solverName() const noexcept
Return solver name.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
SolverPerformance< Type > solve()
Solve returning the solution statistics.
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition faMatrix.H:348
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition faMatrix.C:142
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution on a given patch face.
tmp< Field< Type > > residual() const
Return the matrix residual.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition faMatrix.C:108
Field< Type > & source() noexcept
Definition faMatrix.H:358
const dictionary & solverDict() const
Return the solver dictionary for psi.
Definition faMatrix.C:732
static autoPtr< solver > New(const word &solverName, const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new solver of given type.
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
const scalarField & diag() const
Definition lduMatrix.C:195
void initMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
Initialise the update of interfaced interfaces for matrix operations.
void updateMatrixInterfaces(const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt, const label startRequest) const
Update interfaced interfaces for matrix operations.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
#define DebugInFunction
Report an information message using Foam::Info.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.