Loading...
Searching...
No Matches
scalarTransport.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2023 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 "scalarTransport.H"
30#include "CMULES.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37namespace Foam
38{
39namespace functionObjects
40{
42
44 (
48 );
49}
50}
51
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
55Foam::volScalarField& Foam::functionObjects::scalarTransport::transportedField()
56{
57 if (!foundObject<volScalarField>(fieldName_))
58 {
59 auto tfldPtr = tmp<volScalarField>::New
60 (
62 (
63 fieldName_,
65 mesh_,
69 ),
70 mesh_
71 );
72 store(fieldName_, tfldPtr);
73
74 if (phaseName_ != "none")
75 {
76 mesh_.setFluxRequired(fieldName_);
77 }
78 }
79
80 return lookupObjectRef<volScalarField>(fieldName_);
81}
82
83
84Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
85(
86 const volScalarField& s,
88) const
89{
90 const word Dname("D" + s.name());
91
92 if (constantD_)
93 {
95 (
96 Dname,
98 mesh_,
99 dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
100 );
101 }
102
103 if (nutName_ != "none")
104 {
105 const volScalarField& nutMean =
106 mesh_.lookupObject<volScalarField>(nutName_);
107
108 return tmp<volScalarField>::New(Dname, nutMean);
109 }
110
111 // Incompressible
112 {
113 const auto* turb =
114 findObject<incompressible::turbulenceModel>
115 (
117 );
118
119 if (turb)
120 {
122 (
123 Dname,
125 alphaD_ * turb->nu() + alphaDt_ * turb->nut()
126 );
127 }
128 }
129
130 // Compressible
131 {
132 const auto* turb =
133 findObject<compressible::turbulenceModel>
134 (
136 );
137
138 if (turb)
139 {
141 (
142 Dname,
144 alphaD_ * turb->mu() + alphaDt_ * turb->mut()
145 );
146 }
147 }
148
149
151 (
152 Dname,
154 mesh_,
156 );
157}
158
159
160// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161
163(
164 const word& name,
165 const Time& runTime,
166 const dictionary& dict
167)
168:
170 fvOptions_(mesh_),
171 fieldName_(dict.getOrDefault<word>("field", "s")),
172 schemesField_("unknown-schemesField"),
173 phiName_(dict.getOrDefault<word>("phi", "phi")),
174 rhoName_(dict.getOrDefault<word>("rho", "rho")),
175 nutName_(dict.getOrDefault<word>("nut", "none")),
176 phaseName_(dict.getOrDefault<word>("phase", "none")),
177 phasePhiCompressedName_
178 (
179 dict.getOrDefault<word>("phasePhiCompressed", "alphaPhiUn")
180 ),
181 D_(0),
182 alphaD_(1),
183 alphaDt_(1),
184 tol_(1),
185 nCorr_(0),
186 resetOnStartUp_(false),
187 constantD_(false),
188 bounded01_(dict.getOrDefault("bounded01", true))
189{
190 read(dict);
191
192 // Force creation of transported field so any BCs using it can
193 // look it up
194 volScalarField& s = transportedField();
195
196 if (resetOnStartUp_)
197 {
198 s == Zero;
199 }
200}
201
202
203// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204
206{
208 {
209 return false;
210 }
211
212 dict.readIfPresent("phi", phiName_);
213 dict.readIfPresent("rho", rhoName_);
214 dict.readIfPresent("nut", nutName_);
215 dict.readIfPresent("phase", phaseName_);
216 dict.readIfPresent("phasePhiCompressed", phasePhiCompressedName_);
217
218 schemesField_ = dict.getOrDefault("schemesField", fieldName_);
219
220 dict.readIfPresent("alphaD", alphaD_);
221 dict.readIfPresent("alphaDt", alphaDt_);
222 dict.readIfPresent("tolerance", tol_);
223 dict.readIfPresent("nCorr", nCorr_);
224 dict.readIfPresent("resetOnStartUp", resetOnStartUp_);
225 constantD_ = dict.readIfPresent("D", D_);
226 dict.readIfPresent("bounded01", bounded01_);
227
228 if (dict.found("fvOptions"))
229 {
230 fvOptions_.reset(dict.subDict("fvOptions"));
231 }
232
233 return true;
234}
235
236
238{
239 volScalarField& s = transportedField();
240
241 Log << type() << " execute: " << s.name() << endl;
242
243 const surfaceScalarField& phi =
244 mesh_.lookupObject<surfaceScalarField>(phiName_);
245
246 // Calculate the diffusivity
247 volScalarField D(this->D(s, phi));
248
249 word divScheme("div(phi," + schemesField_ + ")");
250 word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
251
252 // Set under-relaxation coeff
253 scalar relaxCoeff = 0;
254 mesh_.relaxEquation(schemesField_, relaxCoeff);
255
256 // Convergence monitor parameters
257 bool converged = false;
258 int iter = 0;
259
260 // Two phase scalar transport
261 if (phaseName_ != "none")
262 {
263 const volScalarField& alpha =
264 mesh_.lookupObject<volScalarField>(phaseName_);
265
266 const surfaceScalarField& limitedPhiAlpha =
267 mesh_.lookupObject<surfaceScalarField>(phasePhiCompressedName_);
268
269 D *= pos(alpha - 0.99);
270
271 // Reset D dimensions consistent with limitedPhiAlpha
272 D.dimensions().reset(limitedPhiAlpha.dimensions()/dimLength);
273
274 // Solve
276 for (int i = 0; i <= nCorr_; ++i)
277 {
278 fvScalarMatrix sEqn
279 (
280 fvm::ddt(s)
281 + fvm::div(limitedPhiAlpha, s, divScheme)
282 - fvm::laplacian(D, s, laplacianScheme)
283 ==
284 alpha*fvOptions_(s)
285 );
286
287 sEqn.relax(relaxCoeff);
288 fvOptions_.constrain(sEqn);
289
290 ++iter;
291 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
292
293 tTPhiUD = sEqn.flux();
294
295 if (converged) break;
296 }
297
298 if (bounded01_)
299 {
301 (
302 geometricOneField(),
303 s,
304 phi,
305 tTPhiUD.ref(),
306 oneField(),
307 zeroField()
308 );
309 }
310 }
311 else if (phi.dimensions() == dimMass/dimTime)
312 {
313 const volScalarField& rho = lookupObject<volScalarField>(rhoName_);
314
315 for (int i = 0; i <= nCorr_; ++i)
316 {
317 fvScalarMatrix sEqn
318 (
319 fvm::ddt(rho, s)
320 + fvm::div(phi, s, divScheme)
321 - fvm::laplacian(D, s, laplacianScheme)
322 ==
323 fvOptions_(rho, s)
324 );
325
326 sEqn.relax(relaxCoeff);
327
328 fvOptions_.constrain(sEqn);
329
330 ++iter;
331 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
332 if (converged) break;
333 }
334 }
335 else if (phi.dimensions() == dimVolume/dimTime)
336 {
337 for (int i = 0; i <= nCorr_; ++i)
338 {
339 fvScalarMatrix sEqn
340 (
341 fvm::ddt(s)
342 + fvm::div(phi, s, divScheme)
343 - fvm::laplacian(D, s, laplacianScheme)
344 ==
345 fvOptions_(s)
346 );
347
348 sEqn.relax(relaxCoeff);
349
350 fvOptions_.constrain(sEqn);
351
352 ++iter;
353 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
354 if (converged) break;
355 }
356 }
357 else
358 {
360 << "Incompatible dimensions for phi: " << phi.dimensions() << nl
361 << "Dimensions should be " << dimMass/dimTime << " or "
363 }
364
365 if (converged)
366 {
367 Log << type() << ": " << name() << ": "
368 << s.name() << " is converged." << nl
369 << tab << "initial-residual tolerance: " << tol_ << nl
370 << tab << "outer iteration: " << iter << nl;
371 }
373 Log << endl;
374
375 return true;
376}
377
378
380{
381 Log << type() << " write: " << name() << nl
382 << tab << fieldName_ << nl
383 << endl;
384
385 volScalarField& s = transportedField();
386
387 s.write();
388
389 return true;
390}
391
392
393// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
compressible::turbulenceModel & turb
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
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.
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.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
ObjectType & lookupObjectRef(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
Computes the transport equation for a passive scalar in single-phase or two-phase flow,...
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function object output.
virtual bool read(const dictionary &)
Read the function-object dictionary.
scalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition fvMatrix.C:1415
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition oneField.H:50
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable).
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
volScalarField & alpha
dictionary dict
const dimensionedScalar & D