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 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(),
59  inletOutletFvPatchField<scalar>::typeName
60  );
61 
62  forAll(mesh_.boundary(), patchi)
63  {
64  if (isA<wallFvPatch>(mesh_.boundary()[patchi]))
65  {
66  result[patchi] = zeroGradientFvPatchField<scalar>::typeName;
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 template<class GeoField>
95 Foam::functionObjects::age::newField
96 (
97  const word& baseName,
98  const wordList patches
99 ) const
100 {
102  (
103  IOobject
104  (
105  scopedName(baseName),
106  time_.timeName(),
107  mesh_,
110  ),
111  mesh_,
113  patches
114  );
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
121 (
122  const word& name,
123  const Time& runTime,
124  const dictionary& dict
125 )
126 :
128 {
129  read(dict);
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 {
138  {
139  phiName_ = dict.getOrDefault<word>("phi", "phi");
140  rhoName_ = dict.getOrDefault<word>("rho", "rho");
141  schemesField_ = dict.getOrDefault<word>("schemesField", typeName);
142  tolerance_ = dict.getOrDefault<scalar>("tolerance", 1e-5);
143  nCorr_ = dict.getOrDefault<int>("nCorr", 5);
144  diffusion_ = dict.getOrDefault<bool>("diffusion", false);
145 
146  return true;
147  }
148 
149  return false;
150 }
151 
152 
154 {
155  auto tage = tmp<volScalarField>::New
156  (
157  IOobject
158  (
159  typeName,
160  mesh_.time().timeName(),
161  mesh_,
164  false
165  ),
166  mesh_,
168  patchTypes()
169  );
170  volScalarField& age = tage.ref();
171 
172  const word divScheme("div(phi," + schemesField_ + ")");
173 
174  // Set under-relaxation coeff
175  scalar relaxCoeff = 0;
176  if (mesh_.relaxEquation(schemesField_))
177  {
178  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
179  }
180 
182 
183 
184  // This only works because the null constructed inletValue for an
185  // inletOutletFvPatchField is zero. If we needed any other value we would
186  // have to loop over the inletOutlet patches and explicitly set the
187  // inletValues. We would need to change the interface of inletOutlet in
188  // order to do this.
189 
190  const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
191 
192  if (phi.dimensions() == dimMass/dimTime)
193  {
194  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
195 
196  tmp<volScalarField> tmuEff;
197  word laplacianScheme;
198 
199  if (diffusion_)
200  {
201  tmuEff =
202  mesh_.lookupObject<compressible::turbulenceModel>
203  (
205  ).muEff();
206 
207  laplacianScheme =
208  "laplacian(" + tmuEff().name() + ',' + schemesField_ + ")";
209  }
210 
211  for (int i = 0; i <= nCorr_; ++i)
212  {
213  fvScalarMatrix ageEqn
214  (
215  fvm::div(phi, age, divScheme) == rho //+ fvOptions(rho, age)
216  );
217 
218  if (diffusion_)
219  {
220  ageEqn -= fvm::laplacian(tmuEff(), age, laplacianScheme);
221  }
222 
223  ageEqn.relax(relaxCoeff);
224 
225  fvOptions.constrain(ageEqn);
226 
227  if (converged(i, ageEqn.solve().initialResidual()))
228  {
229  break;
230  };
231 
232  fvOptions.correct(age);
233  }
234  }
235  else
236  {
237  tmp<volScalarField> tnuEff;
238  word laplacianScheme;
239 
240  if (diffusion_)
241  {
242  tnuEff =
243  mesh_.lookupObject<incompressible::turbulenceModel>
244  (
246  ).nuEff();
247 
248  laplacianScheme =
249  "laplacian(" + tnuEff().name() + ',' + schemesField_ + ")";
250  }
251 
252  for (int i = 0; i <= nCorr_; ++i)
253  {
254  fvScalarMatrix ageEqn
255  (
256  fvm::div(phi, age, divScheme)
258  );
259 
260  if (diffusion_)
261  {
262  ageEqn -= fvm::laplacian(tnuEff(), age, laplacianScheme);
263  }
264 
265  ageEqn.relax(relaxCoeff);
266 
267  fvOptions.constrain(ageEqn);
268 
269  if (converged(i, ageEqn.solve().initialResidual()))
270  {
271  break;
272  }
273 
274  fvOptions.correct(age);
275  }
276  }
277 
278  Info<< "Min/max age:"
279  << min(age).value() << ' '
280  << max(age).value()
281  << endl;
282 
283  // Workaround
284  word fieldName = typeName;
285 
286  return store(fieldName, tage);
287 }
288 
289 
291 {
292  return writeObject(typeName);
293 }
294 
295 
296 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
nCorr
const int nCorr
Definition: readFluidMultiRegionPIMPLEControls.H:3
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
age.H
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
turbulentTransportModel.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
wallFvPatch.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
zeroGradientFvPatchField.H
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
patchTypes
wordList patchTypes(nPatches)
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::functionObjects::age::read
virtual bool read(const dictionary &)
Read the data.
Definition: age.C:135
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::functionObjects::age::write
virtual bool write()
Write.
Definition: age.C:290
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::functionObjects::age::age
age(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: age.C:121
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::functionObjects::age
Calculates and writes out the time taken for a particle to travel from an inlet to the location....
Definition: age.H:174
Foam::functionObjects::age::execute
virtual bool execute()
Execute.
Definition: age.C:153
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
inletOutletFvPatchField.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::List< word >
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
turbulenceModel.H
turbulentFluidThermoModel.H