MovingPhaseModel.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 | Copyright (C) 2015 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 
26 #include "MovingPhaseModel.H"
27 #include "phaseSystem.H"
30 #include "slipFvPatchFields.H"
32 
33 #include "fvmDdt.H"
34 #include "fvmDiv.H"
35 #include "fvmSup.H"
36 #include "fvcDdt.H"
37 #include "fvcDiv.H"
38 #include "surfaceInterpolate.H"
39 #include "fvMatrix.H"
40 
41 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
42 
43 template<class BasePhaseModel>
46 {
47  word phiName(IOobject::groupName("phi", this->name()));
48 
49  IOobject phiHeader
50  (
51  phiName,
52  U.mesh().time().timeName(),
53  U.mesh(),
54  IOobject::NO_READ
55  );
56 
57  if (phiHeader.headerOk())
58  {
59  Info<< "Reading face flux field " << phiName << endl;
60 
61  return tmp<surfaceScalarField>
62  (
64  (
65  IOobject
66  (
67  phiName,
68  U.mesh().time().timeName(),
69  U.mesh(),
70  IOobject::MUST_READ,
71  IOobject::AUTO_WRITE
72  ),
73  U.mesh()
74  )
75  );
76  }
77  else
78  {
79  Info<< "Calculating face flux field " << phiName << endl;
80 
81  wordList phiTypes
82  (
83  U.boundaryField().size(),
84  calculatedFvPatchScalarField::typeName
85  );
86 
87  forAll(U.boundaryField(), i)
88  {
89  if
90  (
91  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
92  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
93  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
94  )
95  {
96  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
97  }
98  }
99 
100  return tmp<surfaceScalarField>
101  (
103  (
104  IOobject
105  (
106  phiName,
107  U.mesh().time().timeName(),
108  U.mesh(),
109  IOobject::NO_READ,
110  IOobject::AUTO_WRITE
111  ),
112  fvc::interpolate(U) & U.mesh().Sf(),
113  phiTypes
114  )
115  );
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 template<class BasePhaseModel>
124 (
125  const phaseSystem& fluid,
126  const word& phaseName,
127  const label index
128 )
129 :
130  BasePhaseModel(fluid, phaseName, index),
131  U_
132  (
133  IOobject
134  (
135  IOobject::groupName("U", this->name()),
136  fluid.mesh().time().timeName(),
137  fluid.mesh(),
138  IOobject::MUST_READ,
139  IOobject::AUTO_WRITE
140  ),
141  fluid.mesh()
142  ),
143  phi_(phi(U_)),
144  alphaPhi_
145  (
146  IOobject
147  (
148  IOobject::groupName("alphaPhi", this->name()),
149  fluid.mesh().time().timeName(),
150  fluid.mesh()
151  ),
152  fluid.mesh(),
153  dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
154  ),
155  alphaRhoPhi_
156  (
157  IOobject
158  (
159  IOobject::groupName("alphaRhoPhi", this->name()),
160  fluid.mesh().time().timeName(),
161  fluid.mesh()
162  ),
163  fluid.mesh(),
164  dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0)
165  ),
166  DUDt_
167  (
168  IOobject
169  (
170  IOobject::groupName("DUDt", this->name()),
171  fluid.mesh().time().timeName(),
172  fluid.mesh()
173  ),
174  fluid.mesh(),
176  ),
177  divU_(NULL),
178  turbulence_
179  (
181  (
182  *this,
183  this->thermo().rho(),
184  U_,
185  alphaRhoPhi_,
186  phi_,
187  *this
188  )
189  ),
190  continuityError_
191  (
192  IOobject
193  (
194  IOobject::groupName("continuityError", this->name()),
195  fluid.mesh().time().timeName(),
196  fluid.mesh()
197  ),
198  fluid.mesh(),
200  )
201 {
202  phi_.writeOpt() = IOobject::AUTO_WRITE;
203  correctKinematics();
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
208 
209 template<class BasePhaseModel>
211 {}
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
216 template<class BasePhaseModel>
218 {
220 
221  this->fluid().MRF().correctBoundaryVelocity(U_);
222 
223  volScalarField& rho = this->thermo().rho();
224 
225  continuityError_ =
226  fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_)
227  - (this->fluid().fvOptions()(*this, rho)&rho);
228 }
229 
230 
231 template<class BasePhaseModel>
233 {
234  BasePhaseModel::correctKinematics();
235 
236  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
237 }
238 
239 
240 template<class BasePhaseModel>
242 {
243  BasePhaseModel::correctTurbulence();
244 
245  turbulence_->correct();
246 }
247 
248 
249 template<class BasePhaseModel>
251 {
252  BasePhaseModel::correctEnergyTransport();
253  turbulence_->correctEnergyTransport();
254 }
255 
256 
257 template<class BasePhaseModel>
260 {
261  return
262  (
263  fvm::ddt(*this, this->thermo().rho(), U_)
264  + fvm::div(alphaRhoPhi_, U_)
265  - fvm::Sp(continuityError_, U_)
266  + this->fluid().MRF().DDt(*this*this->thermo().rho(), U_)
267  + turbulence_->divDevRhoReff(U_)
268  );
269 }
270 
271 
272 template<class BasePhaseModel>
275 {
276  if (DbyA_.valid())
277  {
278  return DbyA_;
279  }
280  else
281  {
282  return surfaceScalarField::null();
283  }
284 }
285 
286 
287 template<class BasePhaseModel>
289 (
290  const tmp<surfaceScalarField>& DbyA
291 )
292 {
293  DbyA_ = DbyA;
294 }
295 
296 
297 template<class BasePhaseModel>
300 {
301  return U_;
302 }
303 
304 
305 template<class BasePhaseModel>
308 {
309  return U_;
310 }
311 
312 
313 template<class BasePhaseModel>
316 {
317  return DUDt_;
318 }
319 
320 
321 template<class BasePhaseModel>
324 {
325  return divU_;
326 }
327 
328 
329 template<class BasePhaseModel>
330 void
332 (
333  const tmp<volScalarField>& divU
334 )
335 {
336  divU_ = divU;
337 }
338 
339 
340 template<class BasePhaseModel>
343 {
344  return continuityError_;
345 }
346 
347 
348 template<class BasePhaseModel>
351 {
352  return phi_;
353 }
354 
355 
356 template<class BasePhaseModel>
359 {
360  return phi_;
361 }
362 
363 
364 template<class BasePhaseModel>
367 {
368  return alphaPhi_;
369 }
370 
371 
372 template<class BasePhaseModel>
375 {
376  return alphaPhi_;
377 }
378 
379 
380 template<class BasePhaseModel>
383 {
384  return alphaRhoPhi_;
385 }
386 
387 
388 template<class BasePhaseModel>
391 {
392  return alphaRhoPhi_;
393 }
394 
395 
396 template<class BasePhaseModel>
399 {
400  return turbulence_;
401 }
402 
403 
404 // ************************************************************************* //
partialSlipFvPatchFields.H
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::MovingPhaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Foam::MovingPhaseModel::U
virtual tmp< volVectorField > U() const
Constant access the velocity.
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
slipFvPatchFields.H
Foam::dimDensity
const dimensionSet dimDensity
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::phaseCompressibleTurbulenceModel
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
Definition: phaseCompressibleTurbulenceModel.H:45
fvcDiv.H
Calculate the divergence of the given field.
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
phaseSystem.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::MovingPhaseModel::~MovingPhaseModel
virtual ~MovingPhaseModel()
Destructor.
Foam::MovingPhaseModel::phi
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::MovingPhaseModel::DbyA
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::MovingPhaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::MovingPhaseModel::DUDt
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
Foam::MovingPhaseModel::MovingPhaseModel
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Foam::label
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
Foam::Info
messageStream Info
Foam::MovingPhaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
correct
fvOptions correct(rho)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::MovingPhaseModel::UEqn
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:15
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::MovingPhaseModel::turbulence
virtual const phaseCompressibleTurbulenceModel & turbulence() const
Return the turbulence model.
surfaceInterpolate.H
Surface Interpolation.
phaseCompressibleTurbulenceModel.H
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
fvmSup.H
Calculate the matrix for implicit and explicit sources.
rho
rho
Definition: pEqn.H:3
Foam::dimAcceleration
const dimensionSet dimAcceleration
Foam::MovingPhaseModel::alphaRhoPhi
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Constant access the mass flux of the phase.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::MovingPhaseModel::divU
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::rhoThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: rhoThermo.C:158
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
fixedValueFvPatchFields.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::MovingPhaseModel::correct
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
Foam::MovingPhaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
fvcDdt.H
Calculate the first temporal derivative.
Foam::MovingPhaseModel::continuityError
virtual tmp< volScalarField > continuityError() const
Constant access the continuity error.
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Foam::GeometricField::null
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:30
MovingPhaseModel.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fvmDdt.H
Calulate the matrix for the first temporal derivative.
Foam::fvc::DDt
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45