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-------------------------------------------------------------------------------
11License
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
43namespace Foam
44{
45namespace functionObjects
46{
49}
50}
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
54Foam::wordList Foam::functionObjects::age::patchTypes() const
55{
56 wordList result
57 (
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
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
93template<class GeoField>
96(
97 const word& baseName,
98 const wordList patches
99) const
100{
102 (
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 (
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool read()
Re-read model coefficients if they have changed.
bool converged() const noexcept
Has the solver converged?
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const complexVectorField & newField()
Definition: UOprocess.C:95
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Generic dimensioned Type class.
Abstract base-class for Time/database function objects.
Calculates and writes out the time taken for a particle to travel from an inlet to the location....
Definition: age.H:177
virtual bool execute()
Execute.
Definition: age.C:153
virtual bool write()
Write.
Definition: age.C:290
virtual bool read(const dictionary &)
Read the data.
Definition: age.C:135
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
engineTime & runTime
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< word > wordList
A List of words.
Definition: fileName.H:63
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
wordList patchTypes(nPatches)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.