Loading...
Searching...
No Matches
turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018-2020 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "phaseSystem.H"
33#include "mappedPatchBase.H"
34#include "solidThermo.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38const Foam::Enum
39<
40 Foam::compressible::
41 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::regionType
42>
43Foam::compressible::
44turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::regionTypeNames_
45{
46 { regionType::solid, "solid" },
47 { regionType::fluid, "fluid" },
48};
49
50
51const Foam::Enum
52<
53 Foam::compressible::
54 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::KMethodType
55>
56Foam::compressible::
57turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::KMethodTypeNames_
58{
59 { KMethodType::mtSolidThermo, "solidThermo" },
60 { KMethodType::mtLookup, "lookup" },
61 { KMethodType::mtPhaseSystem, "phaseSystem" }
62};
63
64
65// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66
67
68Foam::tmp<Foam::scalarField> Foam::compressible::
69turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
70kappa
71(
72 const scalarField& Tp
73) const
74{
75 const polyMesh& mesh = patch().boundaryMesh().mesh();
76 const label patchi = patch().index();
77
78 switch (method_)
79 {
80 case mtSolidThermo:
81 {
82 const solidThermo& thermo =
83 mesh.lookupObject<solidThermo>(basicThermo::dictName);
84
85 return thermo.kappa(patchi);
86 break;
87 }
88
89 case mtLookup:
90 {
91 {
92 const auto* ptr =
93 mesh.cfindObject<volScalarField>(kappaName_);
94
95 if (ptr)
96 {
97 return patch().patchField(*ptr);
98 }
99 }
100 {
101 const auto* ptr =
102 mesh.cfindObject<volSymmTensorField>(kappaName_);
103
104 if (ptr)
105 {
106 const symmTensorField& KWall = patch().patchField(*ptr);
107
108 const vectorField n(patch().nf());
109
110 return n & KWall & n;
111 }
112 }
113
114 {
116 << "Did not find field " << kappaName_
117 << " on mesh " << mesh.name() << " patch " << patch().name()
118 << nl
119 << " Please set 'kappa' to the name of a volScalarField"
120 << " or volSymmTensorField."
121 << exit(FatalError);
122 }
123 break;
124 }
125
126 case mtPhaseSystem:
127 {
128 // Lookup the fluid model
129 const phaseSystem& fluid =
130 (
131 mesh.lookupObject<phaseSystem>("phaseProperties")
132 );
133
134 auto tkappaEff = tmp<scalarField>::New(patch().size(), Zero);
135 auto& kappaEff = tkappaEff.ref();
136
137 forAll(fluid.phases(), phasei)
138 {
139 const phaseModel& phase = fluid.phases()[phasei];
140
141 const fvPatchScalarField& alpha = phase.boundaryField()[patchi];
142
143 kappaEff += alpha*phase.kappaEff(patchi)();
144 }
145
146 return tkappaEff;
147
148 break;
149 }
150
151 default:
152 {
154 << "Unimplemented method " << KMethodTypeNames_[method_] << nl
155 << "Please set 'kappaMethod' to one of "
156 << flatOutput(KMethodTypeNames_.sortedToc()) << nl
157 << "and 'kappa' to the name of the volScalar"
158 << exit(FatalError);
159 }
160 }
161
162 // Return zero-sized (not nullptr)
163 return tmp<scalarField>::New();
164}
165
166
167// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168
169
170namespace Foam
171{
172namespace compressible
173{
174
175// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
176
179(
180 const fvPatch& p,
181 const DimensionedField<scalar, volMesh>& iF
182)
183:
184 mixedFvPatchScalarField(p, iF),
185 regionType_(fluid),
186 method_(mtLookup),
187 kappaName_("none"),
188 otherPhaseName_("vapor"),
189 TnbrName_("undefined-Tnbr"),
190 qrNbrName_("undefined-qrNbr"),
191 qrName_("undefined-qr")
192{
193 this->refValue() = 0.0;
194 this->refGrad() = 0.0;
195 this->valueFraction() = 1.0;
196}
197
198
199turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
200turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
201(
202 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField& psf,
203 const fvPatch& p,
204 const DimensionedField<scalar, volMesh>& iF,
205 const fvPatchFieldMapper& mapper
206)
207:
208 mixedFvPatchScalarField(psf, p, iF, mapper),
209 regionType_(psf.regionType_),
210 method_(psf.method_),
211 kappaName_(psf.kappaName_),
212 otherPhaseName_(psf.otherPhaseName_),
213 TnbrName_(psf.TnbrName_),
214 qrNbrName_(psf.qrNbrName_),
215 qrName_(psf.qrName_)
216{}
217
218
219turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
220turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
221(
222 const fvPatch& p,
223 const DimensionedField<scalar, volMesh>& iF,
224 const dictionary& dict
225)
226:
227 mixedFvPatchScalarField(p, iF),
228 regionType_(regionTypeNames_.get("region", dict)),
229 method_(KMethodTypeNames_.get("kappaMethod", dict)),
230 kappaName_(dict.getOrDefault<word>("kappa", "none")),
231 otherPhaseName_(dict.get<word>("otherPhase")),
232 TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
233 qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
234 qrName_(dict.getOrDefault<word>("qr", "none"))
235{
236 if (!isA<mappedPatchBase>(this->patch().patch()))
237 {
239 << "' not type '" << mappedPatchBase::typeName << "'"
240 << "\n for patch " << p.name()
241 << " of field " << internalField().name()
242 << " in file " << internalField().objectPath()
243 << exit(FatalError);
244 }
245
246
247 this->readValueEntry(dict, IOobjectOption::MUST_READ);
248
249 if (this->readMixedEntries(dict))
250 {
251 // Full restart
252 }
253 else
254 {
255 // Start from user entered data. Assume fixedValue.
256 refValue() = *this;
257 refGrad() = 0.0;
258 valueFraction() = 1.0;
259 }
260}
261
262
263turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
264turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
265(
266 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField& psf,
267 const DimensionedField<scalar, volMesh>& iF
268)
269:
270 mixedFvPatchScalarField(psf, iF),
271 regionType_(psf.regionType_),
272 method_(psf.method_),
273 kappaName_(psf.kappaName_),
274 otherPhaseName_(psf.otherPhaseName_),
275 TnbrName_(psf.TnbrName_),
276 qrNbrName_(psf.qrNbrName_),
277 qrName_(psf.qrName_)
278{}
279
280
281// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
282
283void turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
284updateCoeffs()
285{
286 if (updated())
287 {
288 return;
289 }
290
291 const polyMesh& mesh = patch().boundaryMesh().mesh();
292
293 // Since we're inside initEvaluate/evaluate there might be processor
294 // comms underway. Change the tag we use.
295 const int oldTag = UPstream::incrMsgType();
296
297 // Get the coupling information from the mappedPatchBase
298 const label patchi = patch().index();
299 const mappedPatchBase& mpp =
300 refCast<const mappedPatchBase>(patch().patch());
301 const polyMesh& nbrMesh = mpp.sampleMesh();
302 const label samplePatchi = mpp.samplePolyPatch().index();
303 const fvPatch& nbrPatch =
304 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
305
306 scalarField& Tp = *this;
307
308 const auto& nbrField =
309 refCast
310 <
311 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
312 >(nbrPatch.lookupPatchField<volScalarField>(TnbrName_));
313
314 // Swap to obtain full local values of neighbour internal field
315 scalarField TcNbr(nbrField.patchInternalField());
316 mpp.distribute(TcNbr);
317
318
319 // Swap to obtain full local values of neighbour K*delta
320 scalarField KDeltaNbr;
321 KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
322 mpp.distribute(KDeltaNbr);
323
324 scalarField KDelta(kappa(Tp)*patch().deltaCoeffs());
325
326 scalarField qr(Tp.size(), 0.0);
327 if (qrName_ != "none")
328 {
329 qr = patch().lookupPatchField<volScalarField>(qrName_);
330 }
331
332 scalarField qrNbr(Tp.size(), 0.0);
333 if (qrNbrName_ != "none")
334 {
335 qrNbr = nbrPatch.lookupPatchField<volScalarField>(qrNbrName_);
336 mpp.distribute(qrNbr);
337 }
338
339
340 if (regionType_ == solid)
341 {
342 // Lookup the fluid model in the nbrFvMesh
343 const phaseSystem& fluid =
344 (
345 nbrMesh.lookupObject<phaseSystem>("phaseProperties")
346 );
347
348 // The BC is applied to the liquid phase of the fluid
349 const phaseModel& liquid
350 (
351 fluid.phases()[nbrField.internalField().group()]
352 );
353
354 const phaseModel& vapor(fluid.phases()[otherPhaseName_]);
355
356
357 scalarField KDeltaLiqNbr;
358 const fvPatchScalarField& alphal = liquid.boundaryField()[samplePatchi];
359 KDeltaLiqNbr =
360 alphal*(liquid.kappaEff(samplePatchi))*nbrPatch.deltaCoeffs();
361 mpp.distribute(KDeltaLiqNbr);
362
363 scalarField KDeltaVapNbr;
364 const fvPatchScalarField& alphav = vapor.boundaryField()[samplePatchi];
365 KDeltaVapNbr =
366 alphav*(vapor.kappaEff(samplePatchi))*nbrPatch.deltaCoeffs();
367 mpp.distribute(KDeltaVapNbr);
368
369 scalarField TvNbr;
370 const fvPatchScalarField& Tv =
371 vapor.thermo().T().boundaryField()[samplePatchi];
372 TvNbr = Tv.patchInternalField();
373 mpp.distribute(TvNbr);
374
375 // TcNbr: liquid Tp
376 // TvNbr: vapour Tp
377 scalarField c(TcNbr*KDeltaLiqNbr + TvNbr*KDeltaVapNbr);
378
379 //valueFraction() = KDeltaNbr/(KDeltaNbr + KDelta);
380 //refValue() = c/KDeltaNbr;
381 scalarField KDeltaLiqVapNbr(KDeltaLiqNbr + KDeltaVapNbr);
382 valueFraction() = KDeltaLiqVapNbr/(KDeltaLiqVapNbr + KDelta);
383 refValue() = c/KDeltaLiqVapNbr;
384 refGrad() = (qr + qrNbr)/kappa(Tp);
385
386 if (debug)
387 {
388 scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
389
390 auto limits = gMinMax(Tp);
391 auto avg = gAverage(Tp);
392
393 Info<< "T solid : " << nl << endl;
394
395 Info
396 << " heat transfer rate from solid:" << Q
397 << " walltemperature "
398 << " min:" << limits.min()
399 << " max:" << limits.max()
400 << " avg:" << avg << nl
401 << endl;
402 }
403 }
404 else if (regionType_ == fluid)
405 {
406 const phaseSystem& fluid =
407 (
408 mesh.lookupObject<phaseSystem>("phaseProperties")
409 );
410
411 const phaseModel& liquid
412 (
413 fluid.phases()[internalField().group()]
414 );
415
416 const phaseModel& vapor(fluid.phases()[otherPhaseName_]);
417
418 const fvPatchScalarField& Tv =
419 vapor.thermo().T().boundaryField()[patchi];
420
421 const fvPatchScalarField& alphav = vapor.boundaryField()[patchi];
422
423 const scalarField KdeltaVap
424 (
425 alphav*(vapor.kappaEff(patchi))*patch().deltaCoeffs()
426 );
427
428 const fvPatchScalarField& alphal = liquid.boundaryField()[patchi];
429
430 const scalarField KdeltaLiq
431 (
432 alphal*(liquid.kappaEff(patchi))*patch().deltaCoeffs()
433 );
434
435 // TcNbr: solid Tp
436 // Tv: vapour Tp
437 const scalarField c(TcNbr*KDeltaNbr + Tv.patchInternalField()*KdeltaVap);
438
439 const scalarField a(KdeltaVap + KDeltaNbr);
440
441 valueFraction() = a/(a + KdeltaLiq);
442 refValue() = c/a;
443 refGrad() = (qr + qrNbr)/kappa(Tp);
444
445 if (debug)
446 {
447 scalarField Tc(patchInternalField());
448 scalarField qLiq((Tp - Tc)*KdeltaLiq);
449 scalarField qVap((Tp - Tv.patchInternalField())*KdeltaVap);
450
451 auto infoT = gMinMax(Tp);
452 auto avgT = gAverage(Tp);
453
454 auto infoLiq = gMinMax(qLiq);
455 auto infoVap = gMinMax(qVap);
456
457 Info<< "T flow : " << nl << endl;
458
459 Info<< " qLiq: " << infoLiq.min() << " - " << infoLiq.max() << nl
460 << " qVap: " << infoVap.min() << " - " << infoVap.max() << nl;
461
462 scalar QLiq = gSum(qLiq*patch().magSf());
463 scalar QVap = gSum(qVap*patch().magSf());
464
465 Info<< " Heat transfer to Liq: " << QLiq << endl;
466 Info<< " Heat transfer to Vap: " << QVap << endl;
467
468 Info<< " walltemperature "
469 << " min:" << infoT.min()
470 << " max:" << infoT.max()
471 << " avg:" << avgT
472 << endl;
473 }
474 }
475 else
476 {
478 << "Unknown phase type. Valid types are: "
479 << regionTypeNames_ << nl << exit(FatalError);
480 }
481
482 UPstream::msgType(oldTag); // Restore tag
483
484 mixedFvPatchScalarField::updateCoeffs();
485}
486
487
488void turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::write
489(
490 Ostream& os
491) const
492{
493 mixedFvPatchField<scalar>::write(os);
494 os.writeEntry("kappaMethod", KMethodTypeNames_[method_]);
495 os.writeEntryIfDifferent<word>("kappa","none", kappaName_);
496
497 os.writeEntry("Tnbr", TnbrName_);
498
499 os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
500 os.writeEntryIfDifferent<word>("qr", "none", qrName_);
501
502 os.writeEntry("region", regionTypeNames_[regionType_]);
503 os.writeEntry("otherPhase", otherPhaseName_);
504}
505
506
507// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
508
510(
511 fvPatchScalarField,
512 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
513);
514
515
516// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
517
518} // End namespace compressible
519} // End namespace Foam
520
521
522// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
alphal
Definition alphavPsi.H:3
twoPhaseSystem & fluid
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors).
Definition Field.C:745
turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
volScalarField & p
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
label phasei
Definition pEqn.H:27
kappaEff
Definition TEqn.H:10
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
constexpr const char *const group
Group name for atomic constants.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
Field< vector > vectorField
Specialisation of Field<T> for vector.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299