mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
ENH: use simpler lookupPatchField form
This commit is contained in:
@ -88,26 +88,29 @@ kappa
|
||||
|
||||
case mtLookup:
|
||||
{
|
||||
if (mesh.foundObject<volScalarField>(kappaName_))
|
||||
{
|
||||
return patch().lookupPatchField<volScalarField, scalar>
|
||||
(
|
||||
kappaName_
|
||||
);
|
||||
const auto* ptr =
|
||||
mesh.cfindObject<volScalarField>(kappaName_);
|
||||
|
||||
if (ptr)
|
||||
{
|
||||
return patch().patchField(*ptr);
|
||||
}
|
||||
}
|
||||
else if (mesh.foundObject<volSymmTensorField>(kappaName_))
|
||||
{
|
||||
const symmTensorField& KWall =
|
||||
patch().lookupPatchField<volSymmTensorField, scalar>
|
||||
(
|
||||
kappaName_
|
||||
);
|
||||
const auto* ptr =
|
||||
mesh.cfindObject<volSymmTensorField>(kappaName_);
|
||||
|
||||
const vectorField n(patch().nf());
|
||||
if (ptr)
|
||||
{
|
||||
const symmTensorField& KWall = patch().patchField(*ptr);
|
||||
|
||||
return n & KWall & n;
|
||||
const vectorField n(patch().nf());
|
||||
|
||||
return n & KWall & n;
|
||||
}
|
||||
}
|
||||
else
|
||||
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Did not find field " << kappaName_
|
||||
@ -117,9 +120,6 @@ kappa
|
||||
<< " or volSymmTensorField."
|
||||
<< exit(FatalError);
|
||||
}
|
||||
|
||||
|
||||
|
||||
break;
|
||||
}
|
||||
|
||||
@ -131,10 +131,8 @@ kappa
|
||||
mesh.lookupObject<phaseSystem>("phaseProperties")
|
||||
);
|
||||
|
||||
tmp<scalarField> kappaEff
|
||||
(
|
||||
new scalarField(patch().size(), 0.0)
|
||||
);
|
||||
auto tkappaEff = tmp<scalarField>::New(patch().size(), Zero);
|
||||
auto& kappaEff = tkappaEff.ref();
|
||||
|
||||
forAll(fluid.phases(), phasei)
|
||||
{
|
||||
@ -142,10 +140,10 @@ kappa
|
||||
|
||||
const fvPatchScalarField& alpha = phase.boundaryField()[patchi];
|
||||
|
||||
kappaEff.ref() += alpha*phase.kappaEff(patchi)();
|
||||
kappaEff += alpha*phase.kappaEff(patchi)();
|
||||
}
|
||||
|
||||
return kappaEff;
|
||||
return tkappaEff;
|
||||
|
||||
break;
|
||||
}
|
||||
@ -161,9 +159,11 @@ kappa
|
||||
}
|
||||
}
|
||||
|
||||
return scalarField(0);
|
||||
// Return zero-sized (not nullptr)
|
||||
return tmp<scalarField>::New();
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
|
||||
@ -308,12 +308,11 @@ updateCoeffs()
|
||||
|
||||
scalarField& Tp = *this;
|
||||
|
||||
const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField&
|
||||
nbrField = refCast
|
||||
<const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField>
|
||||
(
|
||||
nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
|
||||
);
|
||||
const auto& nbrField =
|
||||
refCast
|
||||
<
|
||||
const turbulentTemperatureTwoPhaseRadCoupledMixedFvPatchScalarField
|
||||
>(nbrPatch.lookupPatchField<volScalarField>(TnbrName_));
|
||||
|
||||
// Swap to obtain full local values of neighbour internal field
|
||||
scalarField TcNbr(nbrField.patchInternalField());
|
||||
@ -330,13 +329,13 @@ updateCoeffs()
|
||||
scalarField qr(Tp.size(), 0.0);
|
||||
if (qrName_ != "none")
|
||||
{
|
||||
qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
|
||||
qr = patch().lookupPatchField<volScalarField>(qrName_);
|
||||
}
|
||||
|
||||
scalarField qrNbr(Tp.size(), 0.0);
|
||||
if (qrNbrName_ != "none")
|
||||
{
|
||||
qrNbr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrNbrName_);
|
||||
qrNbr = nbrPatch.lookupPatchField<volScalarField>(qrNbrName_);
|
||||
mpp.distribute(qrNbr);
|
||||
}
|
||||
|
||||
|
||||
Reference in New Issue
Block a user