phaseScalarTransport.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2019-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fixedValueFvPatchField.H"
28 #include "fvcDdt.H"
29 #include "fvcDiv.H"
30 #include "fvmDdt.H"
31 #include "fvmDiv.H"
32 #include "fvmLaplacian.H"
33 #include "fvModels.H"
34 #include "fvConstraints.H"
36 #include "phaseScalarTransport.H"
37 #include "surfaceFields.H"
38 #include "momentumTransportModel.H"
39 #include "wallFvPatch.H"
41 
42 #define PhiDimensionErrorInFunction(phi) \
43  FatalErrorInFunction \
44  << "Incompatible dimensions for " << phi.name() << ": " \
45  << phi.dimensions() << nl \
46  << "Dimensions should be " << dimMass/dimTime << " or " \
47  << dimVolume/dimTime << exit(FatalError)
48 
49 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace functionObjects
54 {
56 
58  (
62  );
63 }
64 }
65 
66 
67 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
68 
69 Foam::volScalarField& Foam::functionObjects::phaseScalarTransport::Phi()
70 {
71  if (!PhiPtr_.valid())
72  {
73  const surfaceScalarField& phi =
75  const volScalarField& p =
77 
78  wordList PhiPatchFieldTypes(mesh_.boundaryMesh().size());
79  forAll(p.boundaryField(), patchi)
80  {
81  PhiPatchFieldTypes[patchi] =
82  p.boundaryField()[patchi].fixesValue()
85  }
86 
87  PhiPtr_.set
88  (
89  new volScalarField
90  (
91  IOobject
92  (
93  typedName(IOobject::groupName("Phi", phaseName_)),
94  time_.name(),
95  mesh_,
98  ),
99  mesh_,
100  dimensionedScalar(phi.dimensions()/dimLength, Zero),
101  PhiPatchFieldTypes
102  )
103  );
104 
105  mesh_.schemes().setFluxRequired(PhiPtr_->name());
106  }
107 
108  return PhiPtr_();
109 }
110 
111 
113 Foam::functionObjects::phaseScalarTransport::alphaPhi()
114 {
115  if (!solveAlphaPhi_)
116  {
117  return mesh_.lookupObject<surfaceScalarField>(alphaPhiName_);
118  }
119 
120  const volScalarField& alpha =
121  mesh_.lookupObject<volScalarField>(alphaName_);
122  const surfaceScalarField& phi =
123  mesh_.lookupObject<surfaceScalarField>(phiName_);
124 
125  // Make a crude guess of the phase flux using default interpolation
126  tmp<surfaceScalarField> tAlphaPhi
127  (
129  (
130  alphaPhiName_,
132  )
133  );
134  surfaceScalarField& alphaPhi = tAlphaPhi.ref();
135 
136  // Get the potential field
137  volScalarField& Phi(this->Phi());
138 
139  // Construct the scheme names
140  const word laplacianScheme = "laplacian(" + pName_ + ")";
141 
142  // Debug writing. Write the material derivative of alpha, before and after
143  // the solution of the potential and the correction of alphaPhi. Before
144  // correction the field should be non-zero, and after it should be
145  // comparable to the solution tolerance.
146  auto writeDDt = [&](const label i)
147  {
148  const volScalarField DDtAlpha
149  (
150  "DDt("
152  (
153  IOobject::member(alpha.name()) + Foam::name(i),
154  IOobject::group(alpha.name())
155  )
156  + ")",
157  fvc::ddt(alpha) + fvc::div(alphaPhi)
158  );
159  Info<< type() << ": Writing " << DDtAlpha.name() << endl;
160  DDtAlpha.write();
161  };
162  if (debug && time_.writeTime())
163  {
164  writeDDt(0);
165  }
166 
167  // Lookup the non-orthogonal solution control
168  nonOrthogonalSolutionControl& control =
169  mesh_.lookupObjectRef<nonOrthogonalSolutionControl>
170  (
171  solutionControl::typeName
172  );
173 
174  // Solve for the potential and correct alphaPhi with the resulting flux
175  if (phi.dimensions() == dimVolume/dimTime)
176  {
177  while (control.correctNonOrthogonal())
178  {
179  fvScalarMatrix PhiEqn
180  (
181  fvm::laplacian(Phi, laplacianScheme)
182  + fvc::ddt(alpha)
183  + fvc::div(alphaPhi)
184  );
185 
186  PhiEqn.solve(pName_);
187 
188  if (control.finalNonOrthogonalIter())
189  {
190  alphaPhi += PhiEqn.flux();
191  }
192  }
193  }
194  else if (phi.dimensions() == dimMass/dimTime)
195  {
196  const volScalarField& rho =
197  mesh_.lookupObject<volScalarField>(rhoName_);
198 
199  while (control.correctNonOrthogonal())
200  {
201  fvScalarMatrix PhiEqn
202  (
203  fvm::laplacian(Phi, laplacianScheme)
204  + fvc::ddt(rho, alpha)
205  + fvc::div(alphaPhi)
206  );
207 
208  PhiEqn.solve(pName_);
209 
210  if (control.finalNonOrthogonalIter())
211  {
212  alphaPhi += PhiEqn.flux();
213  }
214  }
215  }
216  else
217  {
219  }
220 
221  // Debug writing
222  if (debug && time_.writeTime())
223  {
224  writeDDt(1);
225  }
226 
227  return tAlphaPhi;
228 }
229 
230 
232 Foam::functionObjects::phaseScalarTransport::D
233 (
234  const surfaceScalarField& alphaPhi
235 ) const
236 {
237  const word Dname("D" + s_.name());
238 
239  if (diffusivity_ == scalarTransport::diffusivityType::constant)
240  {
241  return volScalarField::New
242  (
243  Dname,
244  mesh_,
246  );
247  }
248  else
249  {
250  const word& nameNoPhase = momentumTransportModel::typeName;
251  const word namePhase = IOobject::groupName(nameNoPhase, phaseName_);
252 
253  // Try looking up the phase transport model, then try the mixture
254  // transport model, then fail with an error relating to the phase
255  // transport model
257  mesh_.foundObject<momentumTransportModel>(namePhase)
258  ? mesh_.lookupObject<momentumTransportModel>(namePhase)
259  : mesh_.foundObject<momentumTransportModel>(nameNoPhase)
260  ? mesh_.lookupObject<momentumTransportModel>(nameNoPhase)
261  : mesh_.lookupObject<momentumTransportModel>(namePhase);
262 
263  return volScalarField::New
264  (
265  Dname,
266  alphal_*turbulence.nu() + alphat_*turbulence.nut()
267  );
268  }
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
273 
275 (
276  const word& name,
277  const Time& runTime,
278  const dictionary& dict
279 )
280 :
281  fvMeshFunctionObject(name, runTime, dict),
282  fieldName_(dict.lookup("field")),
283  phaseName_(IOobject::group(fieldName_)),
284  s_
285  (
286  IOobject
287  (
288  fieldName_,
289  time_.name(),
290  mesh_,
291  IOobject::MUST_READ,
292  IOobject::NO_WRITE
293  ),
294  mesh_
295  ),
296  PhiPtr_(nullptr)
297 {
298  if (phaseName_ == word::null)
299  {
301  << "Field \"" << fieldName_ << "\" does not have a phase extension "
302  << "in its name. If it is associated with \"phaseA\" then it "
303  << "should be named \"" << fieldName_ << ".phaseA\"."
304  << exit(FatalError);
305  }
306 
307  read(dict);
308 }
309 
310 
311 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
312 
314 {}
315 
316 
317 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
318 
320 {
322 
323  solveAlphaPhi_ = dict.lookupOrDefault<bool>("solveAlphaPhi", false);
324 
325  alphaName_ =
326  dict.lookupOrDefault<word>
327  (
328  "alpha",
329  IOobject::groupName("alpha", phaseName_)
330  );
331  const word defaultAlphaPhiName =
332  IOobject::groupName("alphaPhi", phaseName_);
333  alphaPhiName_ =
334  solveAlphaPhi_
335  ? typedName(defaultAlphaPhiName)
336  : dict.lookupOrDefault<word>("alphaPhi", defaultAlphaPhiName);
337  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
338  rhoName_ =
339  dict.lookupOrDefault<word>
340  (
341  "rho",
342  IOobject::groupName("rho", phaseName_)
343  );
344  pName_ = dict.lookupOrDefault<word>("p", "p");
345  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
346 
347  diffusivity_ =
348  scalarTransport::diffusivityTypeNames_.read(dict.lookup("diffusivity"));
349 
350  switch(diffusivity_)
351  {
353  break;
354 
356  dict.lookup("D") >> D_;
357  break;
358 
360  alphal_ =
361  dict.lookupBackwardsCompatible<scalar>({"alphal", "alphaD"});
362  alphat_ =
363  dict.lookupBackwardsCompatible<scalar>({"alphat", "alphaDt"});
364  break;
365  }
366 
367  nCorr_ = dict.lookupOrDefaultBackwardsCompatible<label>
368  (
369  {"nCorrectors", "nCorr"},
370  0
371  );
372  residualAlpha_ = dict.lookupOrDefault<scalar>("residualAlpha", rootSmall);
373  writeAlphaField_ = dict.lookupOrDefault<bool>("writeAlphaField", true);
374 
375  return true;
376 }
377 
378 
380 {
381  return wordList{alphaName_, alphaPhiName_, phiName_, pName_};
382 }
383 
384 
386 {
387  Info<< type() << ": Executing" << endl;
388 
389  const volScalarField& alpha =
390  mesh_.lookupObject<volScalarField>(alphaName_);
391 
392  // Get the phase flux
393  tmp<surfaceScalarField> tAlphaPhi(this->alphaPhi());
394  const surfaceScalarField& alphaPhi = tAlphaPhi();
395 
396  // Get the relaxation coefficient
397  const scalar relaxCoeff =
398  mesh_.solution().relaxEquation(schemesField_)
399  ? mesh_.solution().equationRelaxationFactor(schemesField_)
400  : 0;
401 
402  // Models and constraints
405 
406  // Solve
407  if (alphaPhi.dimensions() == dimVolume/dimTime)
408  {
409  for (int i=0; i<=nCorr_; i++)
410  {
411  fvScalarMatrix fieldEqn
412  (
413  fvm::ddt(alpha, s_)
414  + fvm::div
415  (
416  alphaPhi,
417  s_,
418  "div(" + alphaPhi.name() + "," + schemesField_ + ")"
419  )
420  ==
421  fvModels.source(alpha, s_)
422  - fvm::ddt(residualAlpha_, s_)
423  + fvc::ddt(residualAlpha_, s_)
424  );
425 
426  if (diffusivity_ != scalarTransport::diffusivityType::none)
427  {
428  const volScalarField D(this->D(alphaPhi));
429 
430  fieldEqn -=
432  (
434  s_,
435  "laplacian(" + D.name() + "," + schemesField_ + ")"
436  );
437  }
438 
439  fieldEqn.relax(relaxCoeff);
440  fvConstraints.constrain(fieldEqn);
441  fieldEqn.solve(schemesField_);
443  }
444  }
445  else if (alphaPhi.dimensions() == dimMass/dimTime)
446  {
447  const volScalarField& rho =
448  mesh_.lookupObject<volScalarField>(rhoName_);
449 
450  for (int i=0; i<=nCorr_; i++)
451  {
452  fvScalarMatrix fieldEqn
453  (
454  fvm::ddt(alpha, rho, s_)
455  + fvm::div
456  (
457  alphaPhi,
458  s_,
459  "div(" + alphaPhi.name() + "," + schemesField_ + ")"
460  )
461  ==
462  fvModels.source(alpha, rho, s_)
463  - fvm::ddt(residualAlpha_*rho, s_)
464  + fvc::ddt(residualAlpha_*rho, s_)
465  );
466 
467  if (diffusivity_ != scalarTransport::diffusivityType::none)
468  {
469  const volScalarField D(this->D(alphaPhi));
470 
471  fieldEqn -=
473  (
475  s_,
476  "laplacian(" + D.name() + "," + schemesField_ + ")"
477  );
478  }
479 
480  fieldEqn.relax(relaxCoeff);
481  fvConstraints.constrain(fieldEqn);
482  fieldEqn.solve(schemesField_);
484  }
485  }
486  else
487  {
488  PhiDimensionErrorInFunction(alphaPhi);
489  }
490 
491  Info<< endl;
492 
493  return true;
494 }
495 
496 
498 {
499  s_.write();
500 
501  if (writeAlphaField_)
502  {
503  const volScalarField& alpha =
504  mesh_.lookupObject<volScalarField>(alphaName_);
505 
506  volScalarField alphaS
507  (
508  IOobject
509  (
510  "alpha" + fieldName_.capitalise(),
511  time_.name(),
512  mesh_
513  ),
514  alpha*s_
515  );
516 
517  alphaS.write();
518  }
519 
520  if (PhiPtr_.valid())
521  {
522  PhiPtr_->write();
523  }
524 
525  return true;
526 }
527 
528 
529 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
static const char *const typeName
Definition: Field.H:106
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
word group() const
Return group (extension part of name)
Definition: IOobject.C:321
word member() const
Return member (name without the extension)
Definition: IOobject.C:327
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:54
virtual Ostream & write(const char)=0
Write character.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
const word & name() const
Return const reference to name.
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
virtual bool read(const dictionary &)
Read optional controls.
Evolves a passive scalar transport equation within one phase of a multiphase simulation....
virtual wordList fields() const
Return the list of fields required.
virtual bool execute()
Solve for the evolution of the field.
virtual bool write()
Do nothing. The field is registered and written automatically.
virtual bool read(const dictionary &)
Read the settings from the given dictionary.
phaseScalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
static const NamedEnum< diffusivityType, 3 > diffusivityTypeNames_
Diffusivity type names.
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:447
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the first temporal derivative.
Calculate the divergence of the given field.
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.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
compressibleMomentumTransportModel momentumTransportModel
const char *const group
Group name for atomic constants.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
static const coefficient D("D", dimTemperature, 257.14)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const dimensionSet dimKinematicViscosity
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define PhiDimensionErrorInFunction(phi)
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
volScalarField & p
Foam::surfaceFields.