alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.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 
27 #include "fvPatchFieldMapper.H"
29 
30 #include "twoPhaseSystem.H"
34 #include "ThermalDiffusivity.H"
36 #include "wallFvPatch.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace compressible
43 {
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
47 scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxExp_
48  = 50.0;
49 scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::tolerance_
50  = 0.01;
51 label alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxIters_
52  = 10;
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
57 {
58  if (!isA<wallFvPatch>(patch()))
59  {
61  << "Patch type for patch " << patch().name() << " must be wall\n"
62  << "Current patch type is " << patch().type() << nl
63  << exit(FatalError);
64  }
65 }
66 
67 
68 tmp<scalarField>
70 (
71  const scalarField& Prat
72 ) const
73 {
74  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
75 }
76 
77 
78 tmp<scalarField>
80 (
81  const scalarField& P,
82  const scalarField& Prat
83 ) const
84 {
85  tmp<scalarField> typsf(new scalarField(this->size()));
86  scalarField& ypsf = typsf.ref();
87 
88  forAll(ypsf, facei)
89  {
90  scalar ypt = 11.0;
91 
92  for (int i=0; i<maxIters_; i++)
93  {
94  scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
95  scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
96  scalar yptNew = ypt - f/df;
97 
98  if (yptNew < VSMALL)
99  {
100  ypsf[facei] = 0;
101  }
102  else if (mag(yptNew - ypt) < tolerance_)
103  {
104  ypsf[facei] = yptNew;
105  }
106  else
107  {
108  ypt = yptNew;
109  }
110  }
111 
112  ypsf[facei] = ypt;
113  }
114 
115  return typsf;
116 }
117 
118 tmp<scalarField>
120 (
121  const scalarField& prevAlphat
122 ) const
123 {
124  // Lookup the fluid model
125  const ThermalPhaseChangePhaseSystem
126  <
127  MomentumTransferPhaseSystem<twoPhaseSystem>
128  >& fluid =
129  refCast
130  <
131  const ThermalPhaseChangePhaseSystem
132  <
133  MomentumTransferPhaseSystem<twoPhaseSystem>
134  >
135  >
136  (
137  db().lookupObject<phaseSystem>("phaseProperties")
138  );
139 
140  const phaseModel& liquid
141  (
142  fluid.phase1().name() == internalField().group()
143  ? fluid.phase1()
144  : fluid.phase2()
145  );
146 
147  const label patchi = patch().index();
148 
149  // Retrieve turbulence properties from model
150  const phaseCompressibleTurbulenceModel& turbModel = liquid.turbulence();
151 
152  const scalar Cmu25 = pow025(Cmu_);
153 
154  const scalarField& y = turbModel.y()[patchi];
155 
156  const tmp<scalarField> tmuw = turbModel.mu(patchi);
157  const scalarField& muw = tmuw();
158 
159  const tmp<scalarField> talphaw = liquid.thermo().alpha(patchi);
160  const scalarField& alphaw = talphaw();
161 
162  const tmp<volScalarField> tk = turbModel.k();
163  const volScalarField& k = tk();
164  const fvPatchScalarField& kw = k.boundaryField()[patchi];
165 
166  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
167  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
168  const scalarField magGradUw(mag(Uw.snGrad()));
169 
170  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
171  const fvPatchScalarField& hew =
172  liquid.thermo().he().boundaryField()[patchi];
173 
174  const fvPatchScalarField& Tw =
175  liquid.thermo().T().boundaryField()[patchi];
176 
177  scalarField Tp(Tw.patchInternalField());
178 
179  // Heat flux [W/m2] - lagging alphatw
180  const scalarField qDot
181  (
182  (prevAlphat + alphaw)*hew.snGrad()
183  );
184 
185  scalarField uTau(Cmu25*sqrt(kw));
186 
187  scalarField yPlus(uTau*y/(muw/rhow));
188 
189  scalarField Pr(muw/alphaw);
190 
191  // Molecular-to-turbulent Prandtl number ratio
192  scalarField Prat(Pr/Prt_);
193 
194  // Thermal sublayer thickness
195  scalarField P(this->Psmooth(Prat));
196 
197  scalarField yPlusTherm(this->yPlusTherm(P, Prat));
198 
199  tmp<scalarField> talphatConv(new scalarField(this->size()));
200  scalarField& alphatConv = talphatConv.ref();
201 
202  // Populate boundary values
203  forAll(alphatConv, facei)
204  {
205  // Evaluate new effective thermal diffusivity
206  scalar alphaEff = 0.0;
207  if (yPlus[facei] < yPlusTherm[facei])
208  {
209  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
210  scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
211  scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
212  alphaEff = A/(B + C + VSMALL);
213  }
214  else
215  {
216  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
217  scalar B =
218  qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
219  scalar magUc =
220  uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
221  scalar C =
222  0.5*rhow[facei]*uTau[facei]
223  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
224  alphaEff = A/(B + C + VSMALL);
225  }
226 
227  // Update convective heat transfer turbulent thermal diffusivity
228  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
229  }
230 
231  return talphatConv;
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
236 
239 (
240  const fvPatch& p,
241  const DimensionedField<scalar, volMesh>& iF
242 )
243 :
245  Prt_(0.85),
246  Cmu_(0.09),
247  kappa_(0.41),
248  E_(9.8)
249 {
250  checkType();
251 }
252 
253 
256 (
257  const fvPatch& p,
258  const DimensionedField<scalar, volMesh>& iF,
259  const dictionary& dict
260 )
261 :
263  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
264  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
265  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
266  E_(dict.lookupOrDefault<scalar>("E", 9.8))
267 {}
268 
269 
272 (
274  const fvPatch& p,
275  const DimensionedField<scalar, volMesh>& iF,
276  const fvPatchFieldMapper& mapper
277 )
278 :
280  Prt_(ptf.Prt_),
281  Cmu_(ptf.Cmu_),
282  kappa_(ptf.kappa_),
283  E_(ptf.E_)
284 {}
285 
286 
289 (
291 )
292 :
294  Prt_(awfpsf.Prt_),
295  Cmu_(awfpsf.Cmu_),
296  kappa_(awfpsf.kappa_),
297  E_(awfpsf.E_)
298 {}
299 
300 
303 (
305  const DimensionedField<scalar, volMesh>& iF
306 )
307 :
309  Prt_(awfpsf.Prt_),
310  Cmu_(awfpsf.Cmu_),
311  kappa_(awfpsf.kappa_),
312  E_(awfpsf.E_)
313 {}
314 
315 
316 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
317 
319 {
320  if (updated())
321  {
322  return;
323  }
324 
325  operator==(calcAlphat(*this));
326 
327  fixedValueFvPatchScalarField::updateCoeffs();
328 }
329 
330 
332 (
333  Ostream& os
334 ) const
335 {
337  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
338  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
339  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
340  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
341  dmdt_.writeEntry("dmdt", os);
342  writeEntry("value", os);
343 }
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
349 (
352 );
353 
354 
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 
357 } // End namespace compressible
358 } // End namespace Foam
359 
360 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
fvPatchField< vector > fvPatchVectorField
dimensionedScalar log(const dimensionedScalar &ds)
multiphaseSystem & fluid
Definition: createFields.H:10
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Macros for easy insertion into run-time selection tables.
scalar magUp
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar uTau
dimensionedScalar exp(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const char nl
Definition: Ostream.H:262
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
alphatPhaseChangeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
volScalarField & C
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
scalar yPlus
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
bool compressible
Definition: pEqn.H:30