41 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::regionType
44turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::regionTypeNames_
46 { regionType::solid,
"solid" },
47 { regionType::fluid,
"fluid" },
54 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::KMethodType
57turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::KMethodTypeNames_
59 { KMethodType::mtSolidThermo,
"solidThermo" },
60 { KMethodType::mtLookup,
"lookup" },
61 { KMethodType::mtPhaseSystem,
"phaseSystem" }
69turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
75 const polyMesh&
mesh =
patch().boundaryMesh().mesh();
76 const label patchi =
patch().index();
82 const solidThermo&
thermo =
83 mesh.lookupObject<solidThermo>(basicThermo::dictName);
85 return thermo.kappa(patchi);
97 return patch().patchField(*ptr);
110 return n & KWall &
n;
116 <<
"Did not find field " << kappaName_
117 <<
" on mesh " <<
mesh.name() <<
" patch " <<
patch().name()
119 <<
" Please set 'kappa' to the name of a volScalarField"
120 <<
" or volSymmTensorField."
129 const phaseSystem&
fluid =
131 mesh.lookupObject<phaseSystem>(
"phaseProperties")
134 auto tkappaEff = tmp<scalarField>::New(
patch().size(), Zero);
154 <<
"Unimplemented method " << KMethodTypeNames_[method_] <<
nl
155 <<
"Please set 'kappaMethod' to one of "
157 <<
"and 'kappa' to the name of the volScalar"
163 return tmp<scalarField>::New();
181 const DimensionedField<scalar, volMesh>& iF
184 mixedFvPatchScalarField(
p, iF),
188 otherPhaseName_(
"vapor"),
189 TnbrName_(
"undefined-Tnbr"),
190 qrNbrName_(
"undefined-qrNbr"),
191 qrName_(
"undefined-qr")
193 this->refValue() = 0.0;
194 this->refGrad() = 0.0;
195 this->valueFraction() = 1.0;
199turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
200turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
202 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField& psf,
204 const DimensionedField<scalar, volMesh>& iF,
205 const fvPatchFieldMapper& mapper
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_),
219turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
220turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
223 const DimensionedField<scalar, volMesh>& iF,
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"))
236 if (!isA<mappedPatchBase>(this->
patch().
patch()))
239 <<
"' not type '" << mappedPatchBase::typeName <<
"'"
240 <<
"\n for patch " <<
p.name()
241 <<
" of field " << internalField().name()
242 <<
" in file " << internalField().objectPath()
247 this->readValueEntry(
dict, IOobjectOption::MUST_READ);
249 if (this->readMixedEntries(
dict))
258 valueFraction() = 1.0;
263turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
264turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
266 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField& psf,
267 const DimensionedField<scalar, volMesh>& iF
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_),
283void turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::
291 const polyMesh&
mesh =
patch().boundaryMesh().mesh();
295 const int oldTag = UPstream::incrMsgType();
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];
308 const auto& nbrField =
311 const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
316 mpp.distribute(TcNbr);
321 KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
322 mpp.distribute(KDeltaNbr);
327 if (qrName_ !=
"none")
333 if (qrNbrName_ !=
"none")
336 mpp.distribute(qrNbr);
340 if (regionType_ == solid)
343 const phaseSystem&
fluid =
345 nbrMesh.lookupObject<phaseSystem>(
"phaseProperties")
349 const phaseModel& liquid
351 fluid.phases()[nbrField.internalField().group()]
354 const phaseModel& vapor(
fluid.phases()[otherPhaseName_]);
360 alphal*(liquid.kappaEff(samplePatchi))*nbrPatch.deltaCoeffs();
361 mpp.distribute(KDeltaLiqNbr);
366 alphav*(vapor.kappaEff(samplePatchi))*nbrPatch.deltaCoeffs();
367 mpp.distribute(KDeltaVapNbr);
371 vapor.thermo().
T().boundaryField()[samplePatchi];
372 TvNbr = Tv.patchInternalField();
373 mpp.distribute(TvNbr);
381 scalarField KDeltaLiqVapNbr(KDeltaLiqNbr + KDeltaVapNbr);
382 valueFraction() = KDeltaLiqVapNbr/(KDeltaLiqVapNbr + KDelta);
383 refValue() =
c/KDeltaLiqVapNbr;
384 refGrad() = (qr + qrNbr)/
kappa(Tp);
396 <<
" heat transfer rate from solid:" << Q
397 <<
" walltemperature "
398 <<
" min:" <<
limits.min()
399 <<
" max:" <<
limits.max()
400 <<
" avg:" << avg <<
nl
404 else if (regionType_ ==
fluid)
406 const phaseSystem&
fluid =
408 mesh.lookupObject<phaseSystem>(
"phaseProperties")
411 const phaseModel& liquid
416 const phaseModel& vapor(
fluid.phases()[otherPhaseName_]);
419 vapor.thermo().
T().boundaryField()[patchi];
425 alphav*(vapor.kappaEff(patchi))*
patch().deltaCoeffs()
432 alphal*(liquid.kappaEff(patchi))*
patch().deltaCoeffs()
437 const scalarField c(TcNbr*KDeltaNbr + Tv.patchInternalField()*KdeltaVap);
441 valueFraction() = a/(a + KdeltaLiq);
443 refGrad() = (qr + qrNbr)/
kappa(Tp);
449 scalarField qVap((Tp - Tv.patchInternalField())*KdeltaVap);
459 Info<<
" qLiq: " << infoLiq.min() <<
" - " << infoLiq.max() <<
nl
460 <<
" qVap: " << infoVap.min() <<
" - " << infoVap.max() <<
nl;
462 scalar QLiq =
gSum(qLiq*
patch().magSf());
463 scalar QVap =
gSum(qVap*
patch().magSf());
465 Info<<
" Heat transfer to Liq: " << QLiq <<
endl;
466 Info<<
" Heat transfer to Vap: " << QVap <<
endl;
468 Info<<
" walltemperature "
469 <<
" min:" << infoT.min()
470 <<
" max:" << infoT.max()
478 <<
"Unknown phase type. Valid types are: "
479 << regionTypeNames_ <<
nl <<
exit(FatalError);
482 UPstream::msgType(oldTag);
484 mixedFvPatchScalarField::updateCoeffs();
488void turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField::write
493 mixedFvPatchField<scalar>::write(
os);
494 os.writeEntry(
"kappaMethod", KMethodTypeNames_[method_]);
495 os.writeEntryIfDifferent<word>(
"kappa",
"none", kappaName_);
497 os.writeEntry(
"Tnbr", TnbrName_);
499 os.writeEntryIfDifferent<word>(
"qrNbr",
"none", qrNbrName_);
500 os.writeEntryIfDifferent<word>(
"qr",
"none", qrName_);
502 os.writeEntry(
"region", regionTypeNames_[regionType_]);
503 os.writeEntry(
"otherPhase", otherPhaseName_);
512 turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
Macros for easy insertion into run-time selection tables.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors).
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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)
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.
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.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
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)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.