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 -------------------------------------------------------------------------------
11 License
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 
42 namespace Foam
43 {
44 namespace functionObjects
45 {
46  defineTypeNameAndDebug(scalarTransport, 0);
47 
49  (
50  functionObject,
51  scalarTransport,
52  dictionary
53  );
54 }
55 }
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 Foam::volScalarField& Foam::functionObjects::scalarTransport::transportedField()
61 {
62  if (!foundObject<volScalarField>(fieldName_))
63  {
64  auto tfldPtr = tmp<volScalarField>::New
65  (
66  IOobject
67  (
68  fieldName_,
69  mesh_.time().timeName(),
70  mesh_,
74  ),
75  mesh_
76  );
77  store(fieldName_, tfldPtr);
78 
79  if (phaseName_ != "none")
80  {
81  mesh_.setFluxRequired(fieldName_);
82  }
83  }
84 
85  return lookupObjectRef<volScalarField>(fieldName_);
86 }
87 
88 
89 Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
90 (
91  const volScalarField& s,
92  const surfaceScalarField& phi
93 ) const
94 {
95  const word Dname("D" + s.name());
96 
97  if (constantD_)
98  {
99  return volScalarField::New
100  (
101  Dname,
103  mesh_,
104  dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
105  );
106  }
107 
108  if (nutName_ != "none")
109  {
110  const volScalarField& nutMean =
111  mesh_.lookupObject<volScalarField>(nutName_);
112 
113  return tmp<volScalarField>::New(Dname, nutMean);
114  }
115 
116  // Incompressible
117  {
118  const auto* turb =
119  findObject<incompressible::turbulenceModel>
120  (
122  );
123 
124  if (turb)
125  {
126  return volScalarField::New
127  (
128  Dname,
130  alphaD_ * turb->nu() + alphaDt_ * turb->nut()
131  );
132  }
133  }
134 
135  // Compressible
136  {
137  const auto* turb =
138  findObject<compressible::turbulenceModel>
139  (
141  );
142 
143  if (turb)
144  {
145  return volScalarField::New
146  (
147  Dname,
149  alphaD_ * turb->mu() + alphaDt_ * turb->mut()
150  );
151  }
152  }
153 
154 
155  return volScalarField::New
156  (
157  Dname,
159  mesh_,
160  dimensionedScalar(phi.dimensions()/dimLength, Zero)
161  );
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
166 
167 Foam::functionObjects::scalarTransport::scalarTransport
168 (
169  const word& name,
170  const Time& runTime,
171  const dictionary& dict
172 )
173 :
175  fieldName_(dict.getOrDefault<word>("field", "s")),
176  phiName_(dict.getOrDefault<word>("phi", "phi")),
177  rhoName_(dict.getOrDefault<word>("rho", "rho")),
178  nutName_(dict.getOrDefault<word>("nut", "none")),
179  phaseName_(dict.getOrDefault<word>("phase", "none")),
180  phasePhiCompressedName_
181  (
182  dict.getOrDefault<word>("phasePhiCompressed", "alphaPhiUn")
183  ),
184  D_(0),
185  constantD_(false),
186  nCorr_(0),
187  resetOnStartUp_(false),
188  schemesField_("unknown-schemesField"),
189  fvOptions_(mesh_),
190  bounded01_(dict.getOrDefault("bounded01", true))
191 {
192  read(dict);
193 
194  // Force creation of transported field so any BCs using it can
195  // look it up
196  volScalarField& s = transportedField();
197 
198  if (resetOnStartUp_)
199  {
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
214 {
216 
217  dict.readIfPresent("phi", phiName_);
218  dict.readIfPresent("rho", rhoName_);
219  dict.readIfPresent("nut", nutName_);
220  dict.readIfPresent("phase", phaseName_);
221  dict.readIfPresent("bounded01", bounded01_);
222 
223  schemesField_ = dict.getOrDefault("schemesField", fieldName_);
224  constantD_ = dict.readIfPresent("D", D_);
225  alphaD_ = dict.getOrDefault<scalar>("alphaD", 1);
226  alphaDt_ = dict.getOrDefault<scalar>("alphaDt", 1);
227 
228  dict.readIfPresent("nCorr", nCorr_);
229  dict.readIfPresent("resetOnStartUp", resetOnStartUp_);
230 
231  if (dict.found("fvOptions"))
232  {
233  fvOptions_.reset(dict.subDict("fvOptions"));
234  }
235 
236  return true;
237 }
238 
239 
241 {
242  volScalarField& s = transportedField();
243 
244  Log << type() << " execute: " << s.name() << endl;
245 
246  const surfaceScalarField& phi =
247  mesh_.lookupObject<surfaceScalarField>(phiName_);
248 
249  // Calculate the diffusivity
250  volScalarField D(this->D(s, phi));
251 
252  word divScheme("div(phi," + schemesField_ + ")");
253  word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
254 
255  // Set under-relaxation coeff
256  scalar relaxCoeff = 0;
257  mesh_.relaxEquation(schemesField_, relaxCoeff);
258 
259  // Two phase scalar transport
260  if (phaseName_ != "none")
261  {
262  const volScalarField& alpha =
263  mesh_.lookupObject<volScalarField>(phaseName_);
264 
265  const surfaceScalarField& limitedPhiAlpha =
266  mesh_.lookupObject<surfaceScalarField>(phasePhiCompressedName_);
267 
268  D *= pos(alpha - 0.99);
269 
270  // Reset D dimensions consistent with limitedPhiAlpha
271  D.dimensions().reset(limitedPhiAlpha.dimensions()/dimLength);
272 
273  // Solve
274  tmp<surfaceScalarField> tTPhiUD;
275  for (label i = 0; i <= nCorr_; i++)
276  {
277  fvScalarMatrix sEqn
278  (
279  fvm::ddt(s)
280  + fvm::div(limitedPhiAlpha, s, divScheme)
281  - fvm::laplacian(D, s, laplacianScheme)
282  ==
283  alpha*fvOptions_(s)
284  );
285 
286  sEqn.relax(relaxCoeff);
287  fvOptions_.constrain(sEqn);
288  sEqn.solve(schemesField_);
289 
290  tTPhiUD = sEqn.flux();
291  }
292 
293  if (bounded01_)
294  {
296  (
297  geometricOneField(),
298  s,
299  phi,
300  tTPhiUD.ref(),
301  oneField(),
302  zeroField()
303  );
304  }
305  }
306  else if (phi.dimensions() == dimMass/dimTime)
307  {
308  const volScalarField& rho = lookupObject<volScalarField>(rhoName_);
309 
310  for (label i = 0; i <= nCorr_; i++)
311  {
312 
313  fvScalarMatrix sEqn
314  (
315  fvm::ddt(rho, s)
316  + fvm::div(phi, s, divScheme)
317  - fvm::laplacian(D, s, laplacianScheme)
318  ==
319  fvOptions_(rho, s)
320  );
321 
322  sEqn.relax(relaxCoeff);
323 
324  fvOptions_.constrain(sEqn);
325 
326  sEqn.solve(schemesField_);
327  }
328  }
329  else if (phi.dimensions() == dimVolume/dimTime)
330  {
331  for (label i = 0; i <= nCorr_; i++)
332  {
333  fvScalarMatrix sEqn
334  (
335  fvm::ddt(s)
336  + fvm::div(phi, s, divScheme)
337  - fvm::laplacian(D, s, laplacianScheme)
338  ==
339  fvOptions_(s)
340  );
341 
342  sEqn.relax(relaxCoeff);
343 
344  fvOptions_.constrain(sEqn);
345 
346  sEqn.solve(schemesField_);
347  }
348  }
349  else
350  {
352  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
353  << "Dimensions should be " << dimMass/dimTime << " or "
355  }
357  Log << endl;
358 
359  return true;
360 }
361 
362 
364 {
365  Log << type() << " write: " << name() << nl
366  << tab << fieldName_ << nl
367  << endl;
368 
369  volScalarField& s = transportedField();
370 
371  s.write();
372 
373  return true;
374 }
375 
376 
377 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Calculate the matrix for the laplacian of the field.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
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.
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable)
virtual bool execute()
Calculate the scalarTransport.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
Calculate the matrix for the first temporal derivative.
virtual bool read(const dictionary &)
Read the scalarTransport data.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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:64
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
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
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1095
Calculate the matrix for the divergence of the given field and flux.
virtual bool write()
Do nothing.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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:41
#define Log
Definition: PDRblock.C:28
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual bool read(const dictionary &dict)
Read optional controls.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1416
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const dimensionedScalar & D
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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))
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127