Loading...
Searching...
No Matches
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-2025 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"
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{
57
58 wordList result
59 (
60 patches.size(),
62 );
63
64 forAll(patches, patchi)
65 {
66 if (isA<wallFvPatch>(patches[patchi]))
67 {
68 result[patchi] = fvPatchFieldBase::zeroGradientType();
69 }
70 }
71
72 return result;
73}
74
75
76bool Foam::functionObjects::age::converged
77(
78 const int nCorr,
79 const scalar initialResidual
80) const
81{
82 if (initialResidual < tolerance_)
83 {
84 Info<< "Field " << typeName
85 << " converged in " << nCorr << " correctors"
86 << nl << endl;
87
88 return true;
89 }
91 return false;
92}
93
94
95// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96
98(
99 const word& name,
100 const Time& runTime,
101 const dictionary& dict
102)
103:
106 read(dict);
107}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
113{
115 {
116 return false;
117 }
118
119 phiName_ = dict.getOrDefault<word>("phi", "phi");
120 rhoName_ = dict.getOrDefault<word>("rho", "rho");
121 schemesField_ = dict.getOrDefault<word>("schemesField", typeName);
122 tolerance_ = dict.getOrDefault<scalar>("tolerance", 1e-5);
123 nCorr_ = dict.getOrDefault<int>("nCorr", 5);
124 diffusion_ = dict.getOrDefault<bool>("diffusion", false);
125
126
127 // Detect if compressible or incompressible
128 const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
129 isCompressible_ = (phi.dimensions() == dimMass/dimTime);
130
131 // Store the divergence scheme
132 divScheme_ = word("div(phi," + schemesField_ + ")");
133
134 // Store the Laplacian scheme
135 if (diffusion_)
136 {
137 if (isCompressible_)
138 {
139 laplacianScheme_ = "laplacian(muEff," + schemesField_ + ")";
140 }
141 else
142 {
143 laplacianScheme_ = "laplacian(nuEff," + schemesField_ + ")";
145 }
146
147 return true;
148}
149
150
152{
153 auto* agePtr = getObjectPtr<volScalarField>(typeName);
154 if (!agePtr)
155 {
156 agePtr = new volScalarField
157 (
158 IOobject
159 (
160 typeName,
161 obr().time().timeName(),
162 obr(),
166 ),
167 mesh_,
169 patchTypes()
170 );
171 regIOobject::store(agePtr);
172 }
173 auto& age = *agePtr;
174
175
176 // Set under-relaxation coeff
177 scalar relaxCoeff = 0;
178 mesh_.relaxEquation(schemesField_, relaxCoeff);
179
180 Foam::fv::options& fvOptions(Foam::fv::options::New(mesh_));
181
182
183 // This only works because the null constructed inletValue for an
184 // inletOutletFvPatchField is zero. If we needed any other value we would
185 // have to loop over the inletOutlet patches and explicitly set the
186 // inletValues. We would need to change the interface of inletOutlet in
187 // order to do this.
188
189 const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
190
191 if (isCompressible_)
192 {
193 const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
194
195 tmp<volScalarField> tmuEff;
196
197 if (diffusion_)
198 {
199 tmuEff =
200 mesh_.lookupObject<compressible::turbulenceModel>
201 (
203 ).muEff();
204 }
205
206 for (int i = 0; i <= nCorr_; ++i)
207 {
208 fvScalarMatrix ageEqn
209 (
210 fvm::div(phi, age, divScheme_) == rho //+ fvOptions(rho, age)
211 );
212
213 if (diffusion_)
214 {
215 ageEqn -= fvm::laplacian(tmuEff(), age, laplacianScheme_);
216 }
217
218 ageEqn.relax(relaxCoeff);
219
220 fvOptions.constrain(ageEqn);
221
222 if (converged(i, ageEqn.solve().initialResidual()))
223 {
224 break;
225 }
226
227 fvOptions.correct(age);
228 }
229 }
230 else
231 {
232 tmp<volScalarField> tnuEff;
233
234 if (diffusion_)
235 {
236 tnuEff =
237 mesh_.lookupObject<incompressible::turbulenceModel>
238 (
240 ).nuEff();
241 }
242
243 for (int i = 0; i <= nCorr_; ++i)
244 {
245 fvScalarMatrix ageEqn
246 (
247 fvm::div(phi, age, divScheme_)
249 );
250
251 if (diffusion_)
252 {
253 ageEqn -= fvm::laplacian(tnuEff(), age, laplacianScheme_);
254 }
255
256 ageEqn.relax(relaxCoeff);
257
258 fvOptions.constrain(ageEqn);
259
260 if (converged(i, ageEqn.solve().initialResidual()))
261 {
262 break;
263 }
264
265 fvOptions.correct(age);
266 }
267 }
268
269 Info<< "Min/max age:"
270 << min(age).value() << ' '
271 << max(age).value()
272 << endl;
273
274 return true;
275}
276
277
279{
280 return writeObject(typeName);
281}
282
283
284// ************************************************************************* //
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
static const char *const typeName
Typename for Field.
Definition Field.H:93
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Type & initialResidual() const noexcept
Return initial residual.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Calculates and writes out the time taken for a particle to travel from an inlet to the location....
Definition age.H:173
age(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition age.C:91
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition age.C:105
virtual bool execute()
Execute the function-object operations.
Definition age.C:144
virtual bool write()
Write the function-object results.
Definition age.C:271
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
ObjectType * getObjectPtr(const word &fieldName) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
bool writeObject(const word &fieldName)
Write field if present in the (sub) objectRegistry.
const Time & time() const
Return time database.
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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.
word timeName
Definition getTimeIndex.H:3
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchTypes(nPatches)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299