age.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "age.H"
27 #include "fvmDiv.H"
28 #include "fvmLaplacian.H"
29 #include "fvModels.H"
30 #include "fvConstraints.H"
31 #include "momentumTransportModel.H"
33 #include "wallFvPatch.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(age, 0);
44  addToRunTimeSelectionTable(functionObject, age, dictionary);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 Foam::wordList Foam::functionObjects::age::patchTypes() const
52 {
53  wordList result
54  (
55  mesh_.boundary().size(),
57  );
58 
59  forAll(mesh_.boundary(), patchi)
60  {
61  if (isA<wallFvPatch>(mesh_.boundary()[patchi]))
62  {
64  }
65  }
66 
67  return result;
68 }
69 
70 
71 bool Foam::functionObjects::age::converged
72 (
73  const int nCorr,
74  const scalar initialResidual
75 ) const
76 {
77  if (initialResidual < tolerance_)
78  {
79  Info<< "Field " << typeName
80  << " converged in " << nCorr << " correctors\n" << endl;
81 
82  return true;
83  }
84  else
85  {
86  return false;
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const word& name,
96  const Time& runTime,
97  const dictionary& dict
98 )
99 :
100  fvMeshFunctionObject(name, runTime, dict)
101 {
102  read(dict);
103 }
104 
105 
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
107 
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
117  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
118  nCorr_ = dict.lookupOrDefault<int>("nCorr", 5);
119  schemesField_ = dict.lookupOrDefault<word>("schemesField", typeName);
120  diffusion_ = dict.lookupOrDefault<Switch>("diffusion", false);
121  tolerance_ = dict.lookupOrDefault<scalar>("tolerance", 1e-5);
122 
123  return true;
124 }
125 
126 
128 {
129  return wordList{phiName_, rhoName_};
130 }
131 
132 
134 {
136  (
137  new volScalarField
138  (
139  IOobject
140  (
141  typeName,
142  mesh_.time().timeName(),
143  mesh_,
146  false
147  ),
148  mesh_,
150  patchTypes()
151  )
152  );
153  volScalarField& age = tage.ref();
154 
155  const word divScheme("div(phi," + schemesField_ + ")");
156 
157  // Set under-relaxation coeff
158  scalar relaxCoeff = 0.0;
159  if (mesh_.solution().relaxEquation(schemesField_))
160  {
161  relaxCoeff = mesh_.solution().equationRelaxationFactor(schemesField_);
162  }
163 
166  (
168  );
169 
170  // This only works because the null constructed inletValue for an
171  // inletOutletFvPatchField is zero. If we needed any other value we would
172  // have to loop over the inletOutlet patches and explicitly set the
173  // inletValues. We would need to change the interface of inletOutlet in
174  // order to do this.
175 
176  const surfaceScalarField& phi =
177  mesh_.lookupObject<surfaceScalarField>(phiName_);
178 
179  if (phi.dimensions() == dimMass/dimTime)
180  {
181  const volScalarField& rho =
182  mesh_.lookupObject<volScalarField>(rhoName_);
183 
184  const word laplacianScheme("laplacian(muEff," + schemesField_ + ")");
185 
186  tmp<volScalarField> tnuEff;
187  if (diffusion_)
188  {
189  tnuEff =
190  mesh_.lookupObject<momentumTransportModel>
191  (
192  momentumTransportModel::typeName
193  ).nuEff();
194  }
195 
196  for (int i=0; i<=nCorr_; i++)
197  {
198  fvScalarMatrix ageEqn
199  (
200  fvm::div(phi, age, divScheme) == rho + fvModels.source(rho, age)
201  );
202 
203  if (diffusion_)
204  {
205  ageEqn -= fvm::laplacian(rho*tnuEff(), age, laplacianScheme);
206  }
207 
208  ageEqn.relax(relaxCoeff);
209 
210  fvConstraints.constrain(ageEqn);
211 
212  if (converged(i, ageEqn.solve(schemesField_).initialResidual()))
213  {
214  break;
215  };
216 
217  fvConstraints.constrain(age);
218  }
219  }
220  else
221  {
222  tmp<volScalarField> tnuEff;
223  word laplacianScheme;
224 
225  if (diffusion_)
226  {
227  tnuEff =
228  mesh_.lookupObject<momentumTransportModel>
229  (
230  momentumTransportModel::typeName
231  ).nuEff();
232 
233  laplacianScheme =
234  "laplacian(" + tnuEff().name() + ',' + schemesField_ + ")";
235  }
236 
237  for (int i=0; i<=nCorr_; i++)
238  {
239  fvScalarMatrix ageEqn
240  (
241  fvm::div(phi, age, divScheme)
242  == dimensionedScalar(1) + fvModels.source(age)
243  );
244 
245  if (diffusion_)
246  {
247  ageEqn -= fvm::laplacian(tnuEff(), age, laplacianScheme);
248  }
249 
250  ageEqn.relax(relaxCoeff);
251 
252  fvConstraints.constrain(ageEqn);
253 
254  if (converged(i, ageEqn.solve(schemesField_).initialResidual()))
255  {
256  break;
257  }
258 
259  fvConstraints.constrain(age);
260  }
261  }
262 
263  Info<< "Min/max age:" << min(age).value()
264  << ' ' << max(age).value() << endl;
265 
266  return store(tage);
267 }
268 
269 
271 {
272  return writeObject(typeName);
273 }
274 
275 
276 // ************************************************************************* //
virtual bool execute()
Execute.
Definition: age.C:133
virtual bool read(const dictionary &)
Read the data.
Definition: age.C:114
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static const char *const typeName
Definition: Field.H:105
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
age(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: age.C:94
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
wordList patchTypes(nPatches)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Calculate the matrix for the laplacian of the field.
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.
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
Finite volume models.
Definition: fvModels.H:60
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
const dimensionSet dimMass
defineTypeNameAndDebug(Qdot, 0)
Calculates and writes out the time taken for a particle to travel from an inlet to the location...
Definition: age.H:129
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Abstract base class for turbulence models (RAS, LES and laminar).
Foam::fvModels & fvModels
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
virtual ~age()
Destructor.
Definition: age.C:108
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Finite volume constraints.
Definition: fvConstraints.H:57
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:46
virtual wordList fields() const
Return the list of fields required.
Definition: age.C:127
messageStream Info
virtual bool write()
Write.
Definition: age.C:270
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.