age.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) 2018-2021 OpenFOAM Foundation
9  Copyright (C) 2021-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 "age.H"
30 #include "fvmDiv.H"
31 #include "fvmLaplacian.H"
32 #include "fvOptions.H"
35 #include "turbulenceModel.H"
37 #include "wallFvPatch.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace functionObjects
46 {
47  defineTypeNameAndDebug(age, 0);
48  addToRunTimeSelectionTable(functionObject, age, dictionary);
49 }
50 }
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 Foam::wordList Foam::functionObjects::age::patchTypes() const
55 {
56  wordList result
57  (
58  mesh_.boundary().size(),
60  );
61 
62  forAll(mesh_.boundary(), patchi)
63  {
64  if (isA<wallFvPatch>(mesh_.boundary()[patchi]))
65  {
66  result[patchi] = fvPatchFieldBase::zeroGradientType();
67  }
68  }
69 
70  return result;
71 }
72 
73 
74 bool Foam::functionObjects::age::converged
75 (
76  const int nCorr,
77  const scalar initialResidual
78 ) const
79 {
80  if (initialResidual < tolerance_)
81  {
82  Info<< "Field " << typeName
83  << " converged in " << nCorr << " correctors"
84  << nl << endl;
85 
86  return true;
87  }
88 
89  return false;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
98  const Time& runTime,
99  const dictionary& dict
100 )
101 :
103 {
104  read(dict);
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
113  {
114  phiName_ = dict.getOrDefault<word>("phi", "phi");
115  rhoName_ = dict.getOrDefault<word>("rho", "rho");
116  schemesField_ = dict.getOrDefault<word>("schemesField", typeName);
117  tolerance_ = dict.getOrDefault<scalar>("tolerance", 1e-5);
118  nCorr_ = dict.getOrDefault<int>("nCorr", 5);
119  diffusion_ = dict.getOrDefault<bool>("diffusion", false);
120 
121  return true;
122  }
123 
124  return false;
125 }
126 
127 
129 {
130  auto tage = tmp<volScalarField>::New
131  (
132  IOobject
133  (
134  typeName,
135  mesh_.time().timeName(),
136  mesh_,
140  ),
141  mesh_,
143  patchTypes()
144  );
145  volScalarField& age = tage.ref();
146 
147  const word divScheme("div(phi," + schemesField_ + ")");
148 
149  // Set under-relaxation coeff
150  scalar relaxCoeff = 0;
151  mesh_.relaxEquation(schemesField_, relaxCoeff);
152 
154 
155 
156  // This only works because the null constructed inletValue for an
157  // inletOutletFvPatchField is zero. If we needed any other value we would
158  // have to loop over the inletOutlet patches and explicitly set the
159  // inletValues. We would need to change the interface of inletOutlet in
160  // order to do this.
161 
162  const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
163 
164  if (phi.dimensions() == dimMass/dimTime)
165  {
166  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
167 
168  tmp<volScalarField> tmuEff;
169  word laplacianScheme;
170 
171  if (diffusion_)
172  {
173  tmuEff =
174  mesh_.lookupObject<compressible::turbulenceModel>
175  (
177  ).muEff();
178 
179  laplacianScheme =
180  "laplacian(" + tmuEff().name() + ',' + schemesField_ + ")";
181  }
182 
183  for (int i = 0; i <= nCorr_; ++i)
184  {
185  fvScalarMatrix ageEqn
186  (
187  fvm::div(phi, age, divScheme) == rho //+ fvOptions(rho, age)
188  );
189 
190  if (diffusion_)
191  {
192  ageEqn -= fvm::laplacian(tmuEff(), age, laplacianScheme);
193  }
194 
195  ageEqn.relax(relaxCoeff);
196 
197  fvOptions.constrain(ageEqn);
198 
199  if (converged(i, ageEqn.solve().initialResidual()))
200  {
201  break;
202  };
203 
204  fvOptions.correct(age);
205  }
206  }
207  else
208  {
209  tmp<volScalarField> tnuEff;
210  word laplacianScheme;
211 
212  if (diffusion_)
213  {
214  tnuEff =
215  mesh_.lookupObject<incompressible::turbulenceModel>
216  (
218  ).nuEff();
219 
220  laplacianScheme =
221  "laplacian(" + tnuEff().name() + ',' + schemesField_ + ")";
222  }
223 
224  for (int i = 0; i <= nCorr_; ++i)
225  {
226  fvScalarMatrix ageEqn
227  (
228  fvm::div(phi, age, divScheme)
229  == dimensionedScalar(1) + fvOptions(age)
230  );
231 
232  if (diffusion_)
233  {
234  ageEqn -= fvm::laplacian(tnuEff(), age, laplacianScheme);
235  }
236 
237  ageEqn.relax(relaxCoeff);
238 
239  fvOptions.constrain(ageEqn);
240 
241  if (converged(i, ageEqn.solve().initialResidual()))
242  {
243  break;
244  }
245 
246  fvOptions.correct(age);
247  }
248  }
249 
250  Info<< "Min/max age:"
251  << min(age).value() << ' '
252  << max(age).value()
253  << endl;
254 
255  // Workaround
256  word fieldName = typeName;
257 
258  return store(fieldName, tage);
259 }
260 
261 
263 {
264  return writeObject(typeName);
265 }
266 
267 
268 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
virtual bool execute()
Execute.
Definition: age.C:121
virtual bool read(const dictionary &)
Read the data.
Definition: age.C:103
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
age(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: age.C:89
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
wordList patchTypes(nPatches)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Calculate the matrix for the laplacian of the field.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Finite-volume options.
Definition: fvOptions.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const char *const typeName
Typename for Field.
Definition: Field.H:86
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Reading is optional [identical to LAZY_READ].
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
List of word.
Definition: fileName.H:59
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
virtual bool read(const dictionary &dict)
Read optional controls.
virtual bool write()
Write.
Definition: age.C:255
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127