fvScalarMatrix.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-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 \*---------------------------------------------------------------------------*/
28 
29 #include "fvScalarMatrix.H"
31 #include "profiling.H"
32 #include "PrecisionAdaptor.H"
33 #include "jumpCyclicFvPatchField.H"
35 #include "cyclicAMIPolyPatch.H"
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 template<>
41 (
42  const label patchi,
43  const label facei,
44  const direction,
45  const scalar value
46 )
47 {
48  if (psi_.needReference())
49  {
50  if (Pstream::master())
51  {
52  internalCoeffs_[patchi][facei] +=
53  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
54 
55  boundaryCoeffs_[patchi][facei] +=
56  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
57  *value;
58  }
59  }
60 }
61 
62 
63 template<>
66 (
67  const dictionary& solverControls
68 )
69 {
70  word regionName;
71  if (psi_.mesh().name() != polyMesh::defaultRegion)
72  {
73  regionName = psi_.mesh().name() + "::";
74  }
75  addProfiling(solve, "fvMatrix::solve." + regionName + psi_.name());
76 
77  if (debug)
78  {
79  Info.masterStream(this->mesh().comm())
80  << "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
81  "solver for fvMatrix<scalar>"
82  << endl;
83  }
84 
85  scalarField saveDiag(diag());
86  addBoundaryDiag(diag(), 0);
87 
88  lduInterfaceFieldPtrsList interfaces =
89  psi_.boundaryField().scalarInterfaces();
90 
91 
92  autoPtr<fvMatrix<scalar>::fvSolver> solverPtr
93  (
94  new fvMatrix<scalar>::fvSolver
95  (
96  *this,
98  (
99  psi_.name(),
100  *this,
101  boundaryCoeffs_,
102  internalCoeffs_,
103  interfaces,
104  solverControls
105  )
106  )
107  );
108 
109  diag() = saveDiag;
111  return solverPtr;
112 }
113 
114 
115 template<>
117 (
118  const dictionary& solverControls
119 )
120 {
121  const int logLevel =
122  solverControls.getOrDefault<int>
123  (
124  "log",
126  );
127 
128  auto& psi =
130  (
131  fvMat_.psi()
132  );
133 
134  scalarField saveDiag(fvMat_.diag());
135  fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
136 
137  scalarField totalSource(fvMat_.source());
138  fvMat_.addBoundarySource(totalSource, false);
139 
140  // Assign new solver controls
141  solver_->read(solverControls);
142 
143  solverPerformance solverPerf = solver_->solve
144  (
145  psi.primitiveFieldRef(),
146  totalSource
147  );
148 
149  if (logLevel)
150  {
151  solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
152  }
153 
154  fvMat_.diag() = saveDiag;
155 
156  psi.correctBoundaryConditions();
157 
158  psi.mesh().data().setSolverPerformance(psi.name(), solverPerf);
160  return solverPerf;
161 }
162 
163 
164 template<>
166 (
167  const dictionary& solverControls
168 )
169 {
170  if (debug)
171  {
172  Info.masterStream(this->mesh().comm())
173  << "fvMatrix<scalar>::solveSegregated"
174  "(const dictionary& solverControls) : "
175  "solving fvMatrix<scalar>"
176  << endl;
177  }
178 
179  const int logLevel =
180  solverControls.getOrDefault<int>
181  (
182  "log",
184  );
185 
186  scalarField saveLower;
187  scalarField saveUpper;
188 
189  if (useImplicit_)
190  {
192 
193  if (psi_.mesh().fluxRequired(psi_.name()))
194  {
195  // Save lower/upper for flux calculation
196  if (asymmetric())
197  {
198  saveLower = lower();
199  }
200  saveUpper = upper();
201  }
202 
206  direction cmpt = 0;
207  manipulateMatrix(cmpt);
208  }
209 
210  scalarField saveDiag(diag());
211  addBoundaryDiag(diag(), 0);
212 
213  scalarField totalSource(source_);
214  addBoundarySource(totalSource, false);
215 
216  lduInterfaceFieldPtrsList interfaces;
217  PtrDynList<lduInterfaceField> newInterfaces;
218  if (!useImplicit_)
219  {
220  interfaces = this->psi(0).boundaryField().scalarInterfaces();
221  }
222  else
223  {
224  setInterfaces(interfaces, newInterfaces);
225  }
226 
227  tmp<scalarField> tpsi;
228  if (!useImplicit_)
229  {
230  tpsi.ref
231  (
232  const_cast<GeometricField<scalar, fvPatchField, volMesh>&>
233  (
234  psi_
236  );
237  }
238  else
239  {
240  tpsi = tmp<scalarField>::New(lduAddr().size(), Zero);
241  scalarField& psi = tpsi.ref();
242 
243  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
244  {
245  const label cellOffset = lduMeshPtr()->cellOffsets()[fieldi];
246  const auto& psiInternal = this->psi(fieldi).primitiveField();
247 
248  forAll(psiInternal, localCellI)
249  {
250  psi[cellOffset + localCellI] = psiInternal[localCellI];
251  }
252  }
253  }
254  scalarField& psi = tpsi.ref();
255 
256  // Solver call
258  (
259  this->psi(0).name(),
260  *this,
261  boundaryCoeffs_,
262  internalCoeffs_,
263  interfaces,
264  solverControls
265  )->solve(psi, totalSource);
266 
267  if (useImplicit_)
268  {
269  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
270  {
271  auto& psiInternal =
272  const_cast<GeometricField<scalar, fvPatchField, volMesh>&>
273  (
274  this->psi(fieldi)
275  ).primitiveFieldRef();
276 
277  const label cellOffset = lduMeshPtr()->cellOffsets()[fieldi];
278 
279  forAll(psiInternal, localCellI)
280  {
281  psiInternal[localCellI] = psi[localCellI + cellOffset];
282  }
283  }
284  }
285 
286  if (logLevel)
287  {
288  solverPerf.print(Info.masterStream(mesh().comm()));
289  }
290 
291  diag() = saveDiag;
292 
293  if (useImplicit_)
294  {
295  if (psi_.mesh().fluxRequired(psi_.name()))
296  {
297  // Restore lower/upper
298  if (asymmetric())
299  {
300  lower().setSize(saveLower.size());
301  lower() = saveLower;
302  }
303 
304  upper().setSize(saveUpper.size());
305  upper() = saveUpper;
306  }
307  // Set the original lduMesh
308  setLduMesh(psi_.mesh());
309  }
310 
311  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
312  {
313  auto& localPsi =
314  const_cast<GeometricField<scalar, fvPatchField, volMesh>&>
315  (
316  this->psi(fieldi)
317  );
318 
319  localPsi.correctBoundaryConditions();
320  localPsi.mesh().data().setSolverPerformance
321  (
322  localPsi.name(),
323  solverPerf
324  );
325  }
326 
327  return solverPerf;
328 }
329 
330 
331 template<>
333 {
334  scalarField boundaryDiag(psi_.size(), Zero);
335  addBoundaryDiag(boundaryDiag, 0);
336 
337  const scalarField& psif = psi_.primitiveField();
338  ConstPrecisionAdaptor<solveScalar, scalar> tpsi(psif);
339  const solveScalarField& psi = tpsi();
340 
341  tmp<solveScalarField> tres
342  (
344  (
345  psi,
346  source_ - boundaryDiag*psif,
347  boundaryCoeffs_,
349  0
350  )
351  );
352 
353  ConstPrecisionAdaptor<scalar, solveScalar> tres_s(tres);
354  addBoundarySource(tres_s.ref());
355 
356  return tres_s;
357 }
358 
359 
360 template<>
362 {
363  auto tHphi = volScalarField::New
364  (
365  "H(" + psi_.name() + ')',
366  psi_.mesh(),
367  dimensions_/dimVol,
369  );
370  auto& Hphi = tHphi.ref();
371 
372  Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
373  addBoundarySource(Hphi.primitiveFieldRef());
374 
375  Hphi.primitiveFieldRef() /= psi_.mesh().V();
376  Hphi.correctBoundaryConditions();
377 
378  return tHphi;
379 }
380 
381 
382 template<>
384 {
385  auto tH1 = volScalarField::New
386  (
387  "H(1)",
388  psi_.mesh(),
389  dimensions_/(dimVol*psi_.dimensions()),
391  );
392  auto& H1_ = tH1.ref();
393 
394  H1_.primitiveFieldRef() = lduMatrix::H1();
395  //addBoundarySource(Hphi.primitiveField());
396 
397  H1_.primitiveFieldRef() /= psi_.mesh().V();
398  H1_.correctBoundaryConditions();
399 
400  return tH1;
401 }
402 
403 
404 // ************************************************************************* //
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:117
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition: fvMatrix.C:914
Field< solveScalar > solveScalarField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const labelList & cellOffsets() const
Return cellOffsets.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void residual(solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set...
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition: fvMatrix.C:798
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
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.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1384
scalarField & upper()
Definition: lduMatrix.C:208
tmp< scalarField > H1() const
conserve primitiveFieldRef()+
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
CEqn solve()
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
label nMatrices() const
Definition: fvMatrix.H:392
#define addProfiling(Name, Descr)
Define profiling trigger with specified name and description string.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
SolverPerformance< Type > solve()
Solve returning the solution statistics.
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition: fvMatrix.C:478
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition: fvMatrix.C:821
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
bool asymmetric() const noexcept
Definition: lduMatrix.H:791
const Mesh & mesh() const noexcept
Return mesh.
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.
int debug
Static debugging option.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition: fvMatrix.C:888
OSstream & masterStream(const label communicator)
Return OSstream for output operations on the master process only, Snull on other processes.
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:718
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:734
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
scalarField & lower()
Definition: lduMatrix.C:179
A scalar instance of fvMatrix.
autoPtr< fvSolver > solver()
Construct and return the solver.
void setLduMesh(const lduMesh &m)
Set the LDU mesh containing the addressing is obtained.
Definition: lduMatrix.H:726
messageStream Info
Information stream (stdout output on master, null elsewhere)
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:170
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
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: fvMatrixSolve.C:31
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
tmp< Field< Type > > H(const Field< Type > &) const
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
scalarField & diag()
Definition: lduMatrix.C:197
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
tmp< Field< Type > > residual() const
Return the matrix residual.
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition: fvMatrix.C:668
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127