fvMatrixSolve.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2019 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 \*---------------------------------------------------------------------------*/
28 
29 #include "LduMatrix.H"
30 #include "diagTensorField.H"
31 #include "profiling.H"
32 #include "PrecisionAdaptor.H"
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  if (psi_.needReference())
46  {
47  if (Pstream::master())
48  {
49  internalCoeffs_[patchi][facei].component(cmpt) +=
50  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
51 
52  boundaryCoeffs_[patchi][facei].component(cmpt) +=
53  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
54  *value;
55  }
56  }
57 }
58 
59 
60 template<class Type>
62 (
63  const dictionary& solverControls
64 )
65 {
67  if (psi_.mesh().name() != polyMesh::defaultRegion)
68  {
69  regionName = psi_.mesh().name() + "::";
70  }
71  addProfiling(solve, "fvMatrix::solve." + regionName + psi_.name());
72 
73  if (debug)
74  {
75  Info.masterStream(this->mesh().comm())
76  << "fvMatrix<Type>::solveSegregatedOrCoupled"
77  "(const dictionary& solverControls) : "
78  "solving fvMatrix<Type>"
79  << endl;
80  }
81 
82  label maxIter = -1;
83  if (solverControls.readIfPresent("maxIter", maxIter))
84  {
85  if (maxIter == 0)
86  {
87  return SolverPerformance<Type>();
88  }
89  }
90 
91  word type(solverControls.lookupOrDefault<word>("type", "segregated"));
92 
93  if (type == "segregated")
94  {
95  return solveSegregated(solverControls);
96  }
97  else if (type == "coupled")
98  {
99  return solveCoupled(solverControls);
100  }
101  else
102  {
103  FatalIOErrorInFunction(solverControls)
104  << "Unknown type " << type
105  << "; currently supported solver types are segregated and coupled"
106  << exit(FatalIOError);
107 
108  return SolverPerformance<Type>();
109  }
110 }
111 
112 
113 template<class Type>
115 (
116  const dictionary& solverControls
117 )
118 {
119  if (debug)
120  {
121  Info.masterStream(this->mesh().comm())
122  << "fvMatrix<Type>::solveSegregated"
123  "(const dictionary& solverControls) : "
124  "solving fvMatrix<Type>"
125  << endl;
126  }
127 
130 
131  SolverPerformance<Type> solverPerfVec
132  (
133  "fvMatrix<Type>::solveSegregated",
134  psi.name()
135  );
136 
137  scalarField saveDiag(diag());
138 
139  Field<Type> source(source_);
140 
141  // At this point include the boundary source from the coupled boundaries.
142  // This is corrected for the implicit part by updateMatrixInterfaces within
143  // the component loop.
144  addBoundarySource(source);
145 
146  typename Type::labelType validComponents
147  (
148  psi.mesh().template validComponents<Type>()
149  );
150 
151  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
152  {
153  if (validComponents[cmpt] == -1) continue;
154 
155  // copy field and source
156 
157  scalarField psiCmpt(psi.primitiveField().component(cmpt));
158  addBoundaryDiag(diag(), cmpt);
159 
160  scalarField sourceCmpt(source.component(cmpt));
161 
162  FieldField<Field, scalar> bouCoeffsCmpt
163  (
164  boundaryCoeffs_.component(cmpt)
165  );
166 
167  FieldField<Field, scalar> intCoeffsCmpt
168  (
169  internalCoeffs_.component(cmpt)
170  );
171 
172  lduInterfaceFieldPtrsList interfaces =
173  psi.boundaryField().scalarInterfaces();
174 
175  // Use the initMatrixInterfaces and updateMatrixInterfaces to correct
176  // bouCoeffsCmpt for the explicit part of the coupled boundary
177  // conditions
178  {
179  PrecisionAdaptor<solveScalar, scalar> sourceCmpt_ss(sourceCmpt);
180  ConstPrecisionAdaptor<solveScalar, scalar> psiCmpt_ss(psiCmpt);
181 
182  initMatrixInterfaces
183  (
184  true,
185  bouCoeffsCmpt,
186  interfaces,
187  psiCmpt_ss(),
188  sourceCmpt_ss.ref(),
189  cmpt
190  );
191 
192  updateMatrixInterfaces
193  (
194  true,
195  bouCoeffsCmpt,
196  interfaces,
197  psiCmpt_ss(),
198  sourceCmpt_ss.ref(),
199  cmpt
200  );
201  }
202 
203  solverPerformance solverPerf;
204 
205  // Solver call
206  solverPerf = lduMatrix::solver::New
207  (
208  psi.name() + pTraits<Type>::componentNames[cmpt],
209  *this,
210  bouCoeffsCmpt,
211  intCoeffsCmpt,
212  interfaces,
213  solverControls
214  )->solve(psiCmpt, sourceCmpt, cmpt);
215 
217  {
218  solverPerf.print(Info.masterStream(this->mesh().comm()));
219  }
220 
221  solverPerfVec.replace(cmpt, solverPerf);
222  solverPerfVec.solverName() = solverPerf.solverName();
223 
224  psi.primitiveFieldRef().replace(cmpt, psiCmpt);
225  diag() = saveDiag;
226  }
227 
228  psi.correctBoundaryConditions();
229 
230  psi.mesh().setSolverPerformance(psi.name(), solverPerfVec);
231 
232  return solverPerfVec;
233 }
234 
235 
236 template<class Type>
238 (
239  const dictionary& solverControls
240 )
241 {
242  if (debug)
243  {
244  Info.masterStream(this->mesh().comm())
245  << "fvMatrix<Type>::solveCoupled"
246  "(const dictionary& solverControls) : "
247  "solving fvMatrix<Type>"
248  << endl;
249  }
250 
253 
254  LduMatrix<Type, scalar, scalar> coupledMatrix(psi.mesh());
255  coupledMatrix.diag() = diag();
256  coupledMatrix.upper() = upper();
257  coupledMatrix.lower() = lower();
258  coupledMatrix.source() = source();
259 
260  addBoundaryDiag(coupledMatrix.diag(), 0);
261  addBoundarySource(coupledMatrix.source(), false);
262 
263  coupledMatrix.interfaces() = psi.boundaryFieldRef().interfaces();
264  coupledMatrix.interfacesUpper() = boundaryCoeffs().component(0);
265  coupledMatrix.interfacesLower() = internalCoeffs().component(0);
266 
268  coupledMatrixSolver
269  (
271  (
272  psi.name(),
273  coupledMatrix,
274  solverControls
275  )
276  );
277 
278  SolverPerformance<Type> solverPerf
279  (
280  coupledMatrixSolver->solve(psi)
281  );
282 
284  {
285  solverPerf.print(Info.masterStream(this->mesh().comm()));
286  }
287 
288  psi.correctBoundaryConditions();
289 
290  psi.mesh().setSolverPerformance(psi.name(), solverPerf);
291 
292  return solverPerf;
293 }
294 
295 
296 template<class Type>
298 (
299  const dictionary& solverControls
300 )
301 {
302  return psi_.mesh().solve(*this, solverControls);
303 }
304 
305 
306 template<class Type>
309 {
310  return solver(solverDict());
311 }
312 
313 
314 template<class Type>
316 {
317  return solve(fvMat_.solverDict());
318 }
319 
320 
321 template<class Type>
323 {
324  return this->solve(solverDict());
325 }
326 
327 
328 template<class Type>
330 {
331  tmp<Field<Type>> tres(new Field<Type>(source_));
332  Field<Type>& res = tres.ref();
333 
334  addBoundarySource(res);
335 
336  // Loop over field components
337  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
338  {
339  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
340 
341  scalarField boundaryDiagCmpt(psi_.size(), Zero);
342  addBoundaryDiag(boundaryDiagCmpt, cmpt);
343 
344  FieldField<Field, scalar> bouCoeffsCmpt
345  (
346  boundaryCoeffs_.component(cmpt)
347  );
348 
349  res.replace
350  (
351  cmpt,
353  (
354  psiCmpt,
355  res.component(cmpt) - boundaryDiagCmpt*psiCmpt,
356  bouCoeffsCmpt,
357  psi_.boundaryField().scalarInterfaces(),
358  cmpt
359  )
360  );
361  }
362 
363  return tres;
364 }
365 
366 
367 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvMatrix::residual
tmp< Field< Type > > residual() const
Return the matrix residual.
Definition: fvMatrixSolve.C:329
profiling.H
Foam::fvMatrix::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:322
Foam::fvMatrix::addBoundarySource
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:148
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::fvMatrix::solveSegregatedOrCoupled
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
Foam::messageStream::masterStream
OSstream & masterStream(const label communicator)
Convert to OSstream.
Definition: messageStream.C:69
Foam::fvMatrix::fvSolver::solve
SolverPerformance< Type > solve()
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:315
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::LduMatrix
LduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: LduMatrix.H:72
PrecisionAdaptor.H
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::fvMatrix::addBoundaryDiag
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:114
Foam::stringOps::lower
string lower(const std::string &str)
Return string transformed with std::tolower on each character.
Definition: stringOps.C:1199
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
LduMatrix.H
Foam::SolverPerformance::solverName
const word & solverName() const
Return solver name.
Definition: SolverPerformance.H:149
solve
CEqn solve()
diagTensorField.H
Foam::SolverPerformance::print
void print(Ostream &os) const
Print summary of solver performance to the given stream.
Definition: SolverPerformance.C:99
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::fvMatrix::solveCoupled
SolverPerformance< Type > solveCoupled(const dictionary &)
Solve coupled returning the solution statistics.
Definition: fvMatrixSolve.C:238
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
Foam::solver
Base class for solution control classes.
Definition: solver.H:51
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::stringOps::upper
string upper(const std::string &str)
Return string transformed with std::toupper on each character.
Definition: stringOps.C:1219
Foam::UPtrList< const lduInterfaceField >
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:574
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:114
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::LduMatrix::diag
Field< DType > & diag()
Definition: LduMatrix.C:192
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:562
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fvMatrix::solver
autoPtr< fvSolver > solver()
Construct and return the solver.
Definition: fvMatrixSolve.C:308
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::PrecisionAdaptor
A Field wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:158
Foam::PrecisionAdaptor::ref
FieldType & ref()
Allow modification without const-ref check.
Definition: PrecisionAdaptor.H:219
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::fvMatrix::solverDict
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:1014
Foam::ConstPrecisionAdaptor
A const Field wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:50
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::fvMatrix::solveSegregated
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
Definition: fvMatrixSolve.C:115
Foam::GeometricField< Type, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::lduMatrix::residual
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
Definition: lduMatrixATmul.C:211
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
Foam::fvMatrix::setComponentReference
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution.
Definition: fvMatrixSolve.C:38
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62