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 -------------------------------------------------------------------------------
11 License
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 
27 Description
28  Finite-Area matrix basic solvers.
29 
30 \*---------------------------------------------------------------------------*/
31 
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<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 
53 template<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 =
72 
73  SolverPerformance<Type> solverPerfVec
74  (
75  "faMatrix<Type>::solve",
76  psi.name()
77  );
78 
79  scalarField saveDiag(diag());
80 
81  Field<Type> source(source_);
82  addBoundarySource(source);
83 
84  // Note: make a copy of interfaces: no longer a reference
85  lduInterfaceFieldPtrsList interfaces =
86  psi_.boundaryField().scalarInterfaces();
87 
88  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
89  {
90  // copy field and source
91 
92  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
93  addBoundaryDiag(diag(), cmpt);
94 
95  scalarField sourceCmpt(source.component(cmpt));
96 
97  FieldField<Field, scalar> bouCoeffsCmpt
98  (
99  boundaryCoeffs_.component(cmpt)
100  );
101 
102  FieldField<Field, scalar> intCoeffsCmpt
103  (
104  internalCoeffs_.component(cmpt)
105  );
106 
107  // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
108  // bouCoeffsCmpt for the explicit part of the coupled boundary
109  // conditions
110  {
111  PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
112  ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
113 
114  const label startRequest = UPstream::nRequests();
115 
116  initMatrixInterfaces
117  (
118  true,
119  bouCoeffsCmpt,
120  interfaces,
121  psiCmpt_ss(),
122  sourceCmpt_ss.ref(),
123  cmpt
124  );
125 
126  updateMatrixInterfaces
127  (
128  true,
129  bouCoeffsCmpt,
130  interfaces,
131  psiCmpt_ss(),
132  sourceCmpt_ss.ref(),
133  cmpt,
134  startRequest
135  );
136  }
137 
138  solverPerformance solverPerf;
139 
140  // Solver call
141  solverPerf = lduMatrix::solver::New
142  (
143  psi_.name() + pTraits<Type>::componentNames[cmpt],
144  *this,
145  bouCoeffsCmpt,
146  intCoeffsCmpt,
147  interfaces,
148  solverControls
149  )->solve(psiCmpt, sourceCmpt, cmpt);
150 
151  if (logLevel)
152  {
153  solverPerf.print(Info);
154  }
155 
156  solverPerfVec.replace(cmpt, solverPerf);
157  solverPerfVec.solverName() = solverPerf.solverName();
158 
159  psi.primitiveFieldRef().replace(cmpt, psiCmpt);
160  diag() = saveDiag;
161  }
162 
163  psi.correctBoundaryConditions();
164 
165  psi.mesh().data().setSolverPerformance(psi.name(), solverPerfVec);
167  return solverPerfVec;
168 }
169 
170 
171 template<class Type>
173 {
174  return solve(faMat_.solverDict());
175 }
176 
177 
178 template<class Type>
180 {
181  return this->solve(solverDict(name));
182 }
183 
184 
185 template<class Type>
187 {
188  return this->solve(solverDict());
189 }
190 
191 
192 template<class Type>
194 {
195  tmp<Field<Type>> tres(new Field<Type>(source_));
196  Field<Type>& res = tres().ref();
197 
198  addBoundarySource(res);
199 
200  // Loop over field components
201  for (direction cmpt = 0; cmpt < Type::nComponents; ++cmpt)
202  {
203  scalarField psiCmpt(psi_.internalField().component(cmpt));
204 
205  scalarField boundaryDiagCmpt(psi_.size(), Zero);
206  addBoundaryDiag(boundaryDiagCmpt, cmpt);
207 
208  FieldField<Field, scalar> bouCoeffsCmpt
209  (
210  boundaryCoeffs_.component(cmpt)
211  );
212 
213  res.replace
214  (
215  cmpt,
217  (
218  psiCmpt,
219  res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
220  bouCoeffsCmpt,
222  cmpt
223  )
224  );
225  }
226 
227  return tres;
228 }
229 
230 
231 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
const dictionary & solverDict() const
Return the solver dictionary for psi.
Definition: faMatrix.C:750
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
SolverPerformance< Type > solve()
Solve returning the solution statistics.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
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.
Definition: faMatrixSolve.C:32
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set...
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: faMatrix.C:111
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
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.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
A non-const Field/List wrapper with possible data conversion.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:620
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: faMatrix.C:145
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:51
CEqn solve()
const word & solverName() const noexcept
Return solver name.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SolverPerformance< Type > solve()
Solve returning the solution statistics.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
#define DebugInFunction
Report an information message using Foam::Info.
A const Field/List wrapper with possible data conversion.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:608
Container< Type > & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: refPtrI.H:230
void print(Ostream &os) const
Print summary of solver performance to the given stream.
tmp< Field< Type > > residual() const
Return the matrix residual.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const volScalarField & psi
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...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127