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-2022 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 {
55  defineTypeNameAndDebug(phaseScalarTransport, 0);
56 
58  (
59  functionObject,
60  phaseScalarTransport,
61  dictionary
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 =
74  mesh_.lookupObject<surfaceScalarField>(phiName_);
75  const volScalarField& p =
76  mesh_.lookupObject<volScalarField>(pName_);
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  "Phi" + s_.name(),
94  mesh_.time().timeName(),
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 alphaPhi exists then return it
116  if (mesh_.foundObject<surfaceScalarField>(alphaPhiName_))
117  {
118  return mesh_.lookupObject<surfaceScalarField>(alphaPhiName_);
119  }
120 
121  // Otherwise generate it ...
122  Info<< type() << ": " << surfaceScalarField::typeName << " "
123  << alphaPhiName_ << " was not found, so generating it" << endl;
124 
125  const volScalarField& alpha =
126  mesh_.lookupObject<volScalarField>(alphaName_);
127  const surfaceScalarField& phi =
128  mesh_.lookupObject<surfaceScalarField>(phiName_);
129 
130  // Make a crude guess of the phase flux using default interpolation
131  tmp<surfaceScalarField> tAlphaPhi
132  (
134  (
135  alphaPhiName_,
136  phi*fvc::interpolate(alpha)
137  )
138  );
139  surfaceScalarField& alphaPhi = tAlphaPhi.ref();
140 
141  // Get the potential field
142  volScalarField& Phi(this->Phi());
143 
144  // Construct the scheme names
145  const word laplacianScheme = "laplacian(" + pName_ + ")";
146 
147  // Debug writing. Write the material derivative of alpha, before and after
148  // the solution of the potential and the correction of alphaPhi. Before
149  // correction the field should be non-zero, and after it should be
150  // comparable to the solution tolerance.
151  auto writeDDt = [&](const label i)
152  {
153  const volScalarField DDtAlpha
154  (
155  "DDt("
157  (
158  IOobject::member(alpha.name()) + Foam::name(i),
159  IOobject::group(alpha.name())
160  )
161  + ")",
162  fvc::ddt(alpha) + fvc::div(alphaPhi)
163  );
164  Info<< type() << ": Writing " << DDtAlpha.name() << endl;
165  DDtAlpha.write();
166  };
167  if (debug && mesh_.time().writeTime())
168  {
169  writeDDt(0);
170  }
171 
172  // Lookup the non-orthogonal solution control
173  nonOrthogonalSolutionControl& control =
174  mesh_.lookupObjectRef<nonOrthogonalSolutionControl>
175  (
176  solutionControl::typeName
177  );
178 
179  // Solve for the potential and correct alphaPhi with the resulting flux
180  if (phi.dimensions() == dimVolume/dimTime)
181  {
182  while (control.correctNonOrthogonal())
183  {
184  fvScalarMatrix PhiEqn
185  (
186  fvm::laplacian(Phi, laplacianScheme)
187  + fvc::ddt(alpha)
188  + fvc::div(alphaPhi)
189  );
190 
191  PhiEqn.solve(pName_);
192 
193  if (control.finalNonOrthogonalIter())
194  {
195  alphaPhi += PhiEqn.flux();
196  }
197  }
198  }
199  else if (phi.dimensions() == dimMass/dimTime)
200  {
201  const volScalarField& rho =
202  mesh_.lookupObject<volScalarField>(rhoName_);
203 
204  while (control.correctNonOrthogonal())
205  {
206  fvScalarMatrix PhiEqn
207  (
208  fvm::laplacian(Phi, laplacianScheme)
209  + fvc::ddt(rho, alpha)
210  + fvc::div(alphaPhi)
211  );
212 
213  PhiEqn.solve(pName_);
214 
215  if (control.finalNonOrthogonalIter())
216  {
217  alphaPhi += PhiEqn.flux();
218  }
219  }
220  }
221  else
222  {
224  }
225 
226  // Debug writing
227  if (debug && mesh_.time().writeTime())
228  {
229  writeDDt(1);
230  }
231 
232  return tAlphaPhi;
233 }
234 
235 
237 Foam::functionObjects::phaseScalarTransport::D
238 (
240 ) const
241 {
242  if (constantD_)
243  {
244  return volScalarField::New
245  (
246  "D" + s_.name(),
247  mesh_,
249  );
250  }
251 
252  const word& nameNoPhase = momentumTransportModel::typeName;
253  const word namePhase = IOobject::groupName(nameNoPhase, phaseName_);
254 
255  const word& name =
256  mesh_.foundObject<momentumTransportModel>(namePhase)
257  ? namePhase
258  : mesh_.foundObject<momentumTransportModel>(nameNoPhase)
259  ? nameNoPhase
260  : word::null;
261 
262  if (name == word::null)
263  {
264  return volScalarField::New
265  (
266  "D" + s_.name(),
267  mesh_,
268  dimensionedScalar(alphaPhi.dimensions()/dimLength, 0)
269  );
270  }
271 
273  mesh_.lookupObject<momentumTransportModel>(name);
274 
275  return alphaD_*turbulence.nu() + alphaDt_*turbulence.nut();
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
280 
282 (
283  const word& name,
284  const Time& runTime,
285  const dictionary& dict
286 )
287 :
288  fvMeshFunctionObject(name, runTime, dict),
289  fieldName_(dict.lookup("field")),
290  phaseName_(IOobject::group(fieldName_)),
291  nCorr_(0),
292  residualAlpha_(rootSmall),
293  s_
294  (
295  IOobject
296  (
297  fieldName_,
298  mesh_.time().timeName(),
299  mesh_,
302  ),
303  mesh_
304  ),
305  alphaSPtr_(nullptr),
306  PhiPtr_(nullptr)
307 {
308  if (phaseName_ == word::null)
309  {
311  << "Field \"" << fieldName_ << "\" does not have a phase extension "
312  << "in its name. If it is associated with \"phaseA\" then it "
313  << "should be named \"" << fieldName_ << ".phaseA\"."
314  << exit(FatalError);
315  }
316 
317  read(dict);
318 }
319 
320 
321 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
322 
324 {}
325 
326 
327 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
328 
330 {
332 
333  alphaName_ =
334  dict.lookupOrDefault<word>
335  (
336  "alpha",
337  IOobject::groupName("alpha", phaseName_)
338  );
339  alphaPhiName_ =
340  dict.lookupOrDefault<word>
341  (
342  "alphaPhi",
343  IOobject::groupName("alphaPhi", phaseName_)
344  );
345  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
346  rhoName_ =
347  dict.lookupOrDefault<word>
348  (
349  "rho",
350  IOobject::groupName("rho", phaseName_)
351  );
352  pName_ = dict.lookupOrDefault<word>("p", "p");
353  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
354 
355  constantD_ = dict.readIfPresent("D", D_);
356  alphaD_ = dict.lookupOrDefault<scalar>("alphaD", 1);
357  alphaDt_ = dict.lookupOrDefault<scalar>("alphaDt", 1);
358 
359  dict.readIfPresent("nCorr", nCorr_);
360  dict.readIfPresent("residualAlpha", residualAlpha_);
361  writeAlphaField_ = dict.lookupOrDefault<bool>("writeAlphaField", true);
362 
363  return true;
364 }
365 
366 
368 {
369  return wordList{alphaName_, alphaPhiName_, phiName_, pName_};
370 }
371 
372 
374 {
375  Info<< type() << ": Executing" << endl;
376 
377  const volScalarField& alpha =
378  mesh_.lookupObject<volScalarField>(alphaName_);
379 
380  // Get the phase flux
381  tmp<surfaceScalarField> tAlphaPhi(this->alphaPhi());
382  const surfaceScalarField& alphaPhi = tAlphaPhi();
383 
384  // Get the diffusivity
385  const volScalarField D(this->D(alphaPhi));
386 
387  // Construct the scheme names
388  const word divScheme =
389  "div(" + alphaPhi.name() + "," + schemesField_ + ")";
390  const word laplacianScheme =
391  "laplacian(" + D.name() + "," + schemesField_ + ")";
392 
393  // Get the relaxation coefficient
394  const scalar relaxCoeff =
395  mesh_.solution().relaxEquation(schemesField_)
396  ? mesh_.solution().equationRelaxationFactor(schemesField_)
397  : 0;
398 
401  (
403  );
404 
405  // Solve
406  if (alphaPhi.dimensions() == dimVolume/dimTime)
407  {
408  for (int i=0; i<=nCorr_; i++)
409  {
410  fvScalarMatrix fieldEqn
411  (
412  fvm::ddt(alpha, s_)
413  + fvm::div(alphaPhi, s_, divScheme)
415  (
417  s_,
418  laplacianScheme
419  )
420  ==
421  fvModels.source(alpha, s_)
422  - fvm::ddt(residualAlpha_, s_)
423  + fvc::ddt(residualAlpha_, s_)
424  );
425 
426  fieldEqn.relax(relaxCoeff);
427  fvConstraints.constrain(fieldEqn);
428  fieldEqn.solve(schemesField_);
429  fvConstraints.constrain(s_);
430  }
431  }
432  else if (alphaPhi.dimensions() == dimMass/dimTime)
433  {
434  const volScalarField& rho =
435  mesh_.lookupObject<volScalarField>(rhoName_);
436 
437  for (int i=0; i<=nCorr_; i++)
438  {
439  fvScalarMatrix fieldEqn
440  (
441  fvm::ddt(alpha, rho, s_)
442  + fvm::div(alphaPhi, s_, divScheme)
444  (
445  fvc::interpolate(alpha)*fvc::interpolate(rho*D),
446  s_,
447  laplacianScheme
448  )
449  ==
450  fvModels.source(alpha, rho, s_)
451  - fvm::ddt(residualAlpha_*rho, s_)
452  + fvc::ddt(residualAlpha_*rho, s_)
453  );
454 
455  fieldEqn.relax(relaxCoeff);
456  fvConstraints.constrain(fieldEqn);
457  fieldEqn.solve(schemesField_);
458  fvConstraints.constrain(s_);
459  }
460  }
461  else
462  {
463  PhiDimensionErrorInFunction(alphaPhi);
464  }
465 
466  // Update
467  if (writeAlphaField_)
468  {
469  if (!alphaSPtr_.valid())
470  {
471  alphaSPtr_.set
472  (
473  new volScalarField
474  (
475  IOobject
476  (
477  "alpha"
478  + word(toupper(fieldName_[0]))
479  + fieldName_(1, fieldName_.size() - 1),
480  mesh_.time().timeName(),
481  mesh_,
484  ),
485  mesh_,
486  dimensionedScalar(s_.dimensions(), Zero)
487  )
488  );
489  }
490 
491  alphaSPtr_() = alpha*s_;
492  }
493  else
494  {
495  if (alphaSPtr_.valid())
496  {
497  alphaSPtr_().clear();
498  }
499  }
500 
501  Info<< endl;
502 
503  return true;
504 }
505 
506 
508 {
509  return true;
510 }
511 
512 
513 // ************************************************************************* //
Foam::surfaceFields.
const dimensionSet dimViscosity
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
static const char *const typeName
Definition: Field.H:105
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Calculate the matrix for the laplacian of the field.
word member() const
Return member (name without the extension)
Definition: IOobject.C:330
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
virtual bool read(const dictionary &)
Read the settings from the given dictionary.
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.
const dimensionSet dimLength
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the first temporal derivative.
compressibleMomentumTransportModel momentumTransportModel
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool write()
Do nothing. The field is registered and written automatically.
Foam::fvConstraints & fvConstraints
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual wordList fields() const
Return the list of fields required.
#define PhiDimensionErrorInFunction(phi)
Finite volume models.
Definition: fvModels.H:60
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
static const word null
An empty word.
Definition: word.H:77
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
phi
Definition: correctPhi.H:3
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:72
static const zero Zero
Definition: zero.H:97
Calculate the divergence of the given field.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Finite volume constraints.
Definition: fvConstraints.H:57
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:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimVolume
messageStream Info
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual bool execute()
Solve for the evolution of the field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
phaseScalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864