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-2022 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 "surfaceFields.H"
31#include "fvmDdt.H"
32#include "fvmDiv.H"
33#include "fvmLaplacian.H"
34#include "fvmSup.H"
35#include "CMULES.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace functionObjects
45{
47
49 (
53 );
54}
55}
56
57
58// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59
60Foam::volScalarField& Foam::functionObjects::scalarTransport::transportedField()
61{
62 if (!foundObject<volScalarField>(fieldName_))
63 {
64 auto tfldPtr = tmp<volScalarField>::New
65 (
66 IOobject
67 (
68 fieldName_,
70 mesh_,
73 ),
74 mesh_
75 );
76 store(fieldName_, tfldPtr);
77
78 if (phaseName_ != "none")
79 {
80 mesh_.setFluxRequired(fieldName_);
81 }
82 }
83
84 return lookupObjectRef<volScalarField>(fieldName_);
85}
86
87
89(
90 const volScalarField& s,
92) const
93{
94 const word Dname("D" + s.name());
95
96 if (constantD_)
97 {
99 (
100 IOobject
101 (
102 Dname,
103 mesh_.time().timeName(),
104 mesh_.time(),
107 ),
108 mesh_,
109 dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
110 );
111 }
112
113 if (nutName_ != "none")
114 {
115 const volScalarField& nutMean =
116 mesh_.lookupObject<volScalarField>(nutName_);
117
118 return tmp<volScalarField>::New(Dname, nutMean);
119 }
120
121 // Incompressible
122 {
123 const auto* turb =
124 findObject<incompressible::turbulenceModel>
125 (
127 );
128
129 if (turb)
130 {
132 (
133 Dname,
134 alphaD_ * turb->nu() + alphaDt_ * turb->nut()
135 );
136 }
137 }
138
139 // Compressible
140 {
141 const auto* turb =
142 findObject<compressible::turbulenceModel>
143 (
145 );
146
147 if (turb)
148 {
150 (
151 Dname,
152 alphaD_ * turb->mu() + alphaDt_ * turb->mut()
153 );
154 }
155 }
156
157
159 (
160 IOobject
161 (
162 Dname,
163 mesh_.time().timeName(),
164 mesh_.time(),
167 ),
168 mesh_,
169 dimensionedScalar(phi.dimensions()/dimLength, Zero)
170 );
171}
172
173
174// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
175
177(
178 const word& name,
179 const Time& runTime,
180 const dictionary& dict
181)
182:
184 fieldName_(dict.getOrDefault<word>("field", "s")),
185 phiName_(dict.getOrDefault<word>("phi", "phi")),
186 rhoName_(dict.getOrDefault<word>("rho", "rho")),
187 nutName_(dict.getOrDefault<word>("nut", "none")),
188 phaseName_(dict.getOrDefault<word>("phase", "none")),
189 phasePhiCompressedName_
190 (
191 dict.getOrDefault<word>("phasePhiCompressed", "alphaPhiUn")
192 ),
193 D_(0),
194 constantD_(false),
195 nCorr_(0),
196 resetOnStartUp_(false),
197 schemesField_("unknown-schemesField"),
198 fvOptions_(mesh_),
199 bounded01_(dict.getOrDefault("bounded01", true))
200{
201 read(dict);
202
203 // Force creation of transported field so any BCs using it can
204 // look it up
205 volScalarField& s = transportedField();
206
207 if (resetOnStartUp_)
208 {
210 }
211}
212
213
214// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
215
217{}
218
219
220// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221
223{
225
226 dict.readIfPresent("phi", phiName_);
227 dict.readIfPresent("rho", rhoName_);
228 dict.readIfPresent("nut", nutName_);
229 dict.readIfPresent("phase", phaseName_);
230 dict.readIfPresent("bounded01", bounded01_);
231
232 schemesField_ = dict.getOrDefault("schemesField", fieldName_);
233 constantD_ = dict.readIfPresent("D", D_);
234 alphaD_ = dict.getOrDefault<scalar>("alphaD", 1);
235 alphaDt_ = dict.getOrDefault<scalar>("alphaDt", 1);
236
237 dict.readIfPresent("nCorr", nCorr_);
238 dict.readIfPresent("resetOnStartUp", resetOnStartUp_);
239
240 if (dict.found("fvOptions"))
241 {
242 fvOptions_.reset(dict.subDict("fvOptions"));
243 }
244
245 return true;
246}
247
248
250{
251 volScalarField& s = transportedField();
252
253 Log << type() << " execute: " << s.name() << endl;
254
255 const surfaceScalarField& phi =
256 mesh_.lookupObject<surfaceScalarField>(phiName_);
257
258 // Calculate the diffusivity
259 volScalarField D(this->D(s, phi));
260
261 word divScheme("div(phi," + schemesField_ + ")");
262 word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
263
264 // Set under-relaxation coeff
265 scalar relaxCoeff = 0.0;
266 if (mesh_.relaxEquation(schemesField_))
267 {
268 relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
269 }
270
271 // Two phase scalar transport
272 if (phaseName_ != "none")
273 {
274 const volScalarField& alpha =
275 mesh_.lookupObject<volScalarField>(phaseName_);
276
277 const surfaceScalarField& limitedPhiAlpha =
278 mesh_.lookupObject<surfaceScalarField>(phasePhiCompressedName_);
279
280 D *= pos(alpha - 0.99);
281
282 // Reset D dimensions consistent with limitedPhiAlpha
283 D.dimensions().reset(limitedPhiAlpha.dimensions()/dimLength);
284
285 // Solve
287 for (label i = 0; i <= nCorr_; i++)
288 {
289 fvScalarMatrix sEqn
290 (
291 fvm::ddt(s)
292 + fvm::div(limitedPhiAlpha, s, divScheme)
293 - fvm::laplacian(D, s, laplacianScheme)
294 ==
295 alpha*fvOptions_(s)
296 );
297
298 sEqn.relax(relaxCoeff);
299 fvOptions_.constrain(sEqn);
300 sEqn.solve(mesh_.solverDict(schemesField_));
301
302 tTPhiUD = sEqn.flux();
303 }
304
305 if (bounded01_)
306 {
308 (
310 s,
311 phi,
312 tTPhiUD.ref(),
313 oneField(),
314 zeroField()
315 );
316 }
317 }
318 else if (phi.dimensions() == dimMass/dimTime)
319 {
320 const volScalarField& rho = lookupObject<volScalarField>(rhoName_);
321
322 for (label i = 0; i <= nCorr_; i++)
323 {
324
325 fvScalarMatrix sEqn
326 (
327 fvm::ddt(rho, s)
328 + fvm::div(phi, s, divScheme)
329 - fvm::laplacian(D, s, laplacianScheme)
330 ==
331 fvOptions_(rho, s)
332 );
333
334 sEqn.relax(relaxCoeff);
335
336 fvOptions_.constrain(sEqn);
337
338 sEqn.solve(mesh_.solverDict(schemesField_));
339 }
340 }
341 else if (phi.dimensions() == dimVolume/dimTime)
342 {
343 for (label i = 0; i <= nCorr_; i++)
344 {
345 fvScalarMatrix sEqn
346 (
347 fvm::ddt(s)
348 + fvm::div(phi, s, divScheme)
349 - fvm::laplacian(D, s, laplacianScheme)
350 ==
351 fvOptions_(s)
352 );
353
354 sEqn.relax(relaxCoeff);
355
356 fvOptions_.constrain(sEqn);
357
358 sEqn.solve(mesh_.solverDict(schemesField_));
359 }
360 }
361 else
362 {
364 << "Incompatible dimensions for phi: " << phi.dimensions() << nl
365 << "Dimensions should be " << dimMass/dimTime << " or "
367 }
368
369 Log << endl;
370
371 return true;
372}
373
374
376{
377 Log << type() << " write: " << name() << nl
378 << tab << fieldName_ << nl
379 << endl;
380
381 volScalarField& s = transportedField();
382
383 s.write();
384
385 return true;
386}
387
388
389// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
compressible::turbulenceModel & turb
virtual bool read()
Re-read model coefficients if they have changed.
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
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
Evolves a passive scalar transport equation.
virtual bool execute()
Calculate the scalarTransport.
virtual bool read(const dictionary &)
Read the scalarTransport data.
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
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1443
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
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:56
void setFluxRequired(const word &name) const
Get flux-required for given name, or default.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:56
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
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)
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
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
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: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
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
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52
volScalarField & alpha
dictionary dict
const dimensionedScalar & D
Foam::surfaceFields.