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-2024 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  if (constantD_)
238  {
239  return volScalarField::New
240  (
241  "D" + s_.name(),
242  mesh_,
244  );
245  }
246 
247  const word& nameNoPhase = momentumTransportModel::typeName;
248  const word namePhase = IOobject::groupName(nameNoPhase, phaseName_);
249 
250  const word& name =
251  mesh_.foundObject<momentumTransportModel>(namePhase)
252  ? namePhase
253  : mesh_.foundObject<momentumTransportModel>(nameNoPhase)
254  ? nameNoPhase
255  : word::null;
256 
257  if (name == word::null)
258  {
259  return volScalarField::New
260  (
261  "D" + s_.name(),
262  mesh_,
263  dimensionedScalar(alphaPhi.dimensions()/dimLength, 0)
264  );
265  }
266 
268  mesh_.lookupObject<momentumTransportModel>(name);
269 
270  return alphaD_*turbulence.nu() + alphaDt_*turbulence.nut();
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
275 
277 (
278  const word& name,
279  const Time& runTime,
280  const dictionary& dict
281 )
282 :
283  fvMeshFunctionObject(name, runTime, dict),
284  fieldName_(dict.lookup("field")),
285  phaseName_(IOobject::group(fieldName_)),
286  s_
287  (
288  IOobject
289  (
290  fieldName_,
291  time_.name(),
292  mesh_,
293  IOobject::MUST_READ,
294  IOobject::AUTO_WRITE
295  ),
296  mesh_
297  ),
298  alphaSPtr_(nullptr),
299  PhiPtr_(nullptr)
300 {
301  if (phaseName_ == word::null)
302  {
304  << "Field \"" << fieldName_ << "\" does not have a phase extension "
305  << "in its name. If it is associated with \"phaseA\" then it "
306  << "should be named \"" << fieldName_ << ".phaseA\"."
307  << exit(FatalError);
308  }
309 
310  read(dict);
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
315 
317 {}
318 
319 
320 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
321 
323 {
325 
326  solveAlphaPhi_ = dict.lookupOrDefault<bool>("solveAlphaPhi", false);
327 
328  alphaName_ =
329  dict.lookupOrDefault<word>
330  (
331  "alpha",
332  IOobject::groupName("alpha", phaseName_)
333  );
334  const word defaultAlphaPhiName =
335  IOobject::groupName("alphaPhi", phaseName_);
336  alphaPhiName_ =
337  solveAlphaPhi_
338  ? typedName(defaultAlphaPhiName)
339  : dict.lookupOrDefault<word>("alphaPhi", defaultAlphaPhiName);
340  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
341  rhoName_ =
342  dict.lookupOrDefault<word>
343  (
344  "rho",
345  IOobject::groupName("rho", phaseName_)
346  );
347  pName_ = dict.lookupOrDefault<word>("p", "p");
348  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
349 
350  constantD_ = dict.readIfPresent("D", D_);
351  alphaD_ = dict.lookupOrDefault<scalar>("alphaD", 1);
352  alphaDt_ = dict.lookupOrDefault<scalar>("alphaDt", 1);
353 
354  nCorr_ = dict.lookupOrDefault<label>("nCorr", 0);
355  residualAlpha_ = dict.lookupOrDefault<scalar>("residualAlpha", rootSmall);
356  writeAlphaField_ = dict.lookupOrDefault<bool>("writeAlphaField", true);
357 
358  return true;
359 }
360 
361 
363 {
364  return wordList{alphaName_, alphaPhiName_, phiName_, pName_};
365 }
366 
367 
369 {
370  Info<< type() << ": Executing" << endl;
371 
372  const volScalarField& alpha =
373  mesh_.lookupObject<volScalarField>(alphaName_);
374 
375  // Get the phase flux
376  tmp<surfaceScalarField> tAlphaPhi(this->alphaPhi());
377  const surfaceScalarField& alphaPhi = tAlphaPhi();
378 
379  // Get the diffusivity
380  const volScalarField D(this->D(alphaPhi));
381 
382  // Construct the scheme names
383  const word divScheme =
384  "div(" + alphaPhi.name() + "," + schemesField_ + ")";
385  const word laplacianScheme =
386  "laplacian(" + D.name() + "," + schemesField_ + ")";
387 
388  // Get the relaxation coefficient
389  const scalar relaxCoeff =
390  mesh_.solution().relaxEquation(schemesField_)
391  ? mesh_.solution().equationRelaxationFactor(schemesField_)
392  : 0;
393 
396  (
398  );
399 
400  // Solve
401  if (alphaPhi.dimensions() == dimVolume/dimTime)
402  {
403  for (int i=0; i<=nCorr_; i++)
404  {
405  fvScalarMatrix fieldEqn
406  (
407  fvm::ddt(alpha, s_)
408  + fvm::div(alphaPhi, s_, divScheme)
410  (
412  s_,
413  laplacianScheme
414  )
415  ==
416  fvModels.source(alpha, s_)
417  - fvm::ddt(residualAlpha_, s_)
418  + fvc::ddt(residualAlpha_, s_)
419  );
420 
421  fieldEqn.relax(relaxCoeff);
422  fvConstraints.constrain(fieldEqn);
423  fieldEqn.solve(schemesField_);
425  }
426  }
427  else if (alphaPhi.dimensions() == dimMass/dimTime)
428  {
429  const volScalarField& rho =
430  mesh_.lookupObject<volScalarField>(rhoName_);
431 
432  for (int i=0; i<=nCorr_; i++)
433  {
434  fvScalarMatrix fieldEqn
435  (
436  fvm::ddt(alpha, rho, s_)
437  + fvm::div(alphaPhi, s_, divScheme)
439  (
441  s_,
442  laplacianScheme
443  )
444  ==
445  fvModels.source(alpha, rho, s_)
446  - fvm::ddt(residualAlpha_*rho, s_)
447  + fvc::ddt(residualAlpha_*rho, s_)
448  );
449 
450  fieldEqn.relax(relaxCoeff);
451  fvConstraints.constrain(fieldEqn);
452  fieldEqn.solve(schemesField_);
454  }
455  }
456  else
457  {
458  PhiDimensionErrorInFunction(alphaPhi);
459  }
460 
461  // Update
462  if (writeAlphaField_)
463  {
464  if (!alphaSPtr_.valid())
465  {
466  alphaSPtr_.set
467  (
468  new volScalarField
469  (
470  IOobject
471  (
472  "alpha"
473  + word(toupper(fieldName_[0]))
474  + fieldName_(1, fieldName_.size() - 1),
475  time_.name(),
476  mesh_,
479  ),
480  mesh_,
481  dimensionedScalar(s_.dimensions(), Zero)
482  )
483  );
484  }
485 
486  alphaSPtr_() = alpha*s_;
487  }
488  else
489  {
490  if (alphaSPtr_.valid())
491  {
492  alphaSPtr_().clear();
493  }
494  }
495 
496  Info<< endl;
497 
498  return true;
499 }
500 
501 
503 {
504  return true;
505 }
506 
507 
508 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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:324
word member() const
Return member (name without the extension)
Definition: IOobject.C:330
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
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.
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.
virtual bool read(const dictionary &)
Read optional controls.
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:1751
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:503
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:404
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
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
#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
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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:64
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
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.