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-2016 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 "fvcFlux.H"
39 #include "surfaceInterpolate.H"
40 #include "fvMatrix.H"
41 
42 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
43 
44 template<class BasePhaseModel>
47 {
48  word phiName(IOobject::groupName("phi", this->name()));
49 
50  IOobject phiHeader
51  (
52  phiName,
53  U.mesh().time().timeName(),
54  U.mesh(),
55  IOobject::NO_READ
56  );
57 
58  if (phiHeader.headerOk())
59  {
60  Info<< "Reading face flux field " << phiName << endl;
61 
62  return tmp<surfaceScalarField>
63  (
65  (
66  IOobject
67  (
68  phiName,
69  U.mesh().time().timeName(),
70  U.mesh(),
71  IOobject::MUST_READ,
72  IOobject::AUTO_WRITE
73  ),
74  U.mesh()
75  )
76  );
77  }
78  else
79  {
80  Info<< "Calculating face flux field " << phiName << endl;
81 
82  wordList phiTypes
83  (
84  U.boundaryField().size(),
85  calculatedFvPatchScalarField::typeName
86  );
87 
88  forAll(U.boundaryField(), i)
89  {
90  if
91  (
92  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
93  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
94  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
95  )
96  {
97  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
98  }
99  }
100 
101  return tmp<surfaceScalarField>
102  (
104  (
105  IOobject
106  (
107  phiName,
108  U.mesh().time().timeName(),
109  U.mesh(),
110  IOobject::NO_READ,
111  IOobject::AUTO_WRITE
112  ),
113  fvc::flux(U),
114  phiTypes
115  )
116  );
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
123 template<class BasePhaseModel>
125 (
126  const phaseSystem& fluid,
127  const word& phaseName,
128  const label index
129 )
130 :
131  BasePhaseModel(fluid, phaseName, index),
132  U_
133  (
134  IOobject
135  (
136  IOobject::groupName("U", this->name()),
137  fluid.mesh().time().timeName(),
138  fluid.mesh(),
139  IOobject::MUST_READ,
140  IOobject::AUTO_WRITE
141  ),
142  fluid.mesh()
143  ),
144  phi_(phi(U_)),
145  alphaPhi_
146  (
147  IOobject
148  (
149  IOobject::groupName("alphaPhi", this->name()),
150  fluid.mesh().time().timeName(),
151  fluid.mesh()
152  ),
153  fluid.mesh(),
154  dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
155  ),
156  alphaRhoPhi_
157  (
158  IOobject
159  (
160  IOobject::groupName("alphaRhoPhi", this->name()),
161  fluid.mesh().time().timeName(),
162  fluid.mesh()
163  ),
164  fluid.mesh(),
165  dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0)
166  ),
167  DUDt_
168  (
169  IOobject
170  (
171  IOobject::groupName("DUDt", this->name()),
172  fluid.mesh().time().timeName(),
173  fluid.mesh()
174  ),
175  fluid.mesh(),
177  ),
178  divU_(NULL),
179  turbulence_
180  (
182  (
183  *this,
184  this->thermo().rho(),
185  U_,
186  alphaRhoPhi_,
187  phi_,
188  *this
189  )
190  ),
191  continuityError_
192  (
193  IOobject
194  (
195  IOobject::groupName("continuityError", this->name()),
196  fluid.mesh().time().timeName(),
197  fluid.mesh()
198  ),
199  fluid.mesh(),
201  )
202 {
203  phi_.writeOpt() = IOobject::AUTO_WRITE;
204  correctKinematics();
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
209 
210 template<class BasePhaseModel>
212 {}
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 template<class BasePhaseModel>
219 {
221 
222  this->fluid().MRF().correctBoundaryVelocity(U_);
223 
224  volScalarField& rho = this->thermo().rho();
225 
226  continuityError_ =
227  fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_)
228  - (this->fluid().fvOptions()(*this, rho)&rho);
229 }
230 
231 
232 template<class BasePhaseModel>
234 {
235  BasePhaseModel::correctKinematics();
236 
237  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
238 }
239 
240 
241 template<class BasePhaseModel>
243 {
244  BasePhaseModel::correctTurbulence();
245 
246  turbulence_->correct();
247 }
248 
249 
250 template<class BasePhaseModel>
252 {
253  BasePhaseModel::correctEnergyTransport();
254  turbulence_->correctEnergyTransport();
255 }
256 
257 
258 template<class BasePhaseModel>
261 {
262  return
263  (
264  fvm::ddt(*this, this->thermo().rho(), U_)
265  + fvm::div(alphaRhoPhi_, U_)
266  - fvm::Sp(continuityError_, U_)
267  + this->fluid().MRF().DDt(*this*this->thermo().rho(), U_)
268  + turbulence_->divDevRhoReff(U_)
269  );
270 }
271 
272 
273 template<class BasePhaseModel>
276 {
277  if (DbyA_.valid())
278  {
279  return DbyA_;
280  }
281  else
282  {
283  return surfaceScalarField::null();
284  }
285 }
286 
287 
288 template<class BasePhaseModel>
290 (
291  const tmp<surfaceScalarField>& DbyA
292 )
293 {
294  DbyA_ = DbyA;
295 }
296 
297 
298 template<class BasePhaseModel>
301 {
302  return U_;
303 }
304 
305 
306 template<class BasePhaseModel>
309 {
310  return U_;
311 }
312 
313 
314 template<class BasePhaseModel>
317 {
318  return DUDt_;
319 }
320 
321 
322 template<class BasePhaseModel>
325 {
326  return divU_;
327 }
328 
329 
330 template<class BasePhaseModel>
331 void
333 (
334  const tmp<volScalarField>& divU
335 )
336 {
337  divU_ = divU;
338 }
339 
340 
341 template<class BasePhaseModel>
344 {
345  return continuityError_;
346 }
347 
348 
349 template<class BasePhaseModel>
352 {
353  return phi_;
354 }
355 
356 
357 template<class BasePhaseModel>
360 {
361  return phi_;
362 }
363 
364 
365 template<class BasePhaseModel>
368 {
369  return alphaPhi_;
370 }
371 
372 
373 template<class BasePhaseModel>
376 {
377  return alphaPhi_;
378 }
379 
380 
381 template<class BasePhaseModel>
384 {
385  return alphaRhoPhi_;
386 }
387 
388 
389 template<class BasePhaseModel>
392 {
393  return alphaRhoPhi_;
394 }
395 
396 
397 template<class BasePhaseModel>
400 {
401  return turbulence_;
402 }
403 
404 
405 // ************************************************************************* //
surfaceScalarField & phi
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
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
multiphaseSystem & fluid
Definition: createFields.H:10
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual void correctTurbulence()
Correct the turbulence.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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:52
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:31
const dimensionSet dimAcceleration
dynamicFvMesh & mesh
IOMRFZoneList & MRF
virtual ~MovingPhaseModel()
Destructor.
virtual const phaseCompressibleTurbulenceModel & turbulence() const
Return the turbulence model.
Calculate the face-flux of the given field.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calulate the matrix for the first temporal derivative.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/ubuntu/OpenFOAM-4.1/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H"1{alphav=max(min((rho-rholSat)/(rhovSat-rholSat), scalar(1)), scalar(0));alphal=1.0-alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
const dimensionSet dimDensity
static const GeometricField< scalar, fvsPatchField, surfaceMesh > & null()
Return a null geometric field.
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
messageStream Info
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volVectorField > U() const
Constant access the velocity.
virtual tmp< volScalarField > continuityError() const
Constant access the continuity error.
Calculate the matrix for implicit and explicit sources.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Constant access the mass flux of the phase.