Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "mappedPatchBase.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace compressible
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
mixedFvPatchScalarField(p, iF),
temperatureCoupledBase(patch(), "undefined", "undefined-K"),
neighbourFieldName_("undefined-neighbourFieldName")
{
this->refValue() = 0.0;
this->refGrad() = 0.0;
this->valueFraction() = 1.0;
}
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
(
const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
mixedFvPatchScalarField(ptf, p, iF, mapper),
temperatureCoupledBase(patch(), ptf.KMethod(), ptf.kappaName()),
neighbourFieldName_(ptf.neighbourFieldName_)
{}
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
mixedFvPatchScalarField(p, iF),
temperatureCoupledBase(patch(), dict),
neighbourFieldName_(dict.lookup("neighbourFieldName"))
if (!isA<mappedPatchBase>(this->patch().patch()))
{
FatalErrorIn
(
"turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::"
"turbulentTemperatureCoupledBaffleMixedFvPatchScalarField\n"
"(\n"
" const fvPatch& p,\n"
" const DimensionedField<scalar, volMesh>& iF,\n"
" const dictionary& dict\n"
")\n"
) << "\n patch type '" << p.type()
<< "' not type '" << mappedPatchBase::typeName << "'"
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
<< "\n for patch " << p.name()
<< " of field " << dimensionedInternalField().name()
<< " in file " << dimensionedInternalField().objectPath()
<< exit(FatalError);
}
fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
if (dict.found("refValue"))
{
// Full restart
refValue() = scalarField("refValue", dict, p.size());
refGrad() = scalarField("refGradient", dict, p.size());
valueFraction() = scalarField("valueFraction", dict, p.size());
}
else
{
// Start from user entered data. Assume fixedValue.
refValue() = *this;
refGrad() = 0.0;
valueFraction() = 1.0;
}
}
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
(
const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
const DimensionedField<scalar, volMesh>& iF
)
:
mixedFvPatchScalarField(wtcsf, iF),
temperatureCoupledBase(patch(), wtcsf.KMethod(), wtcsf.kappaName()),
neighbourFieldName_(wtcsf.neighbourFieldName_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
// Since we're inside initEvaluate/evaluate there might be processor
// comms underway. Change the tag we use.
int oldTag = UPstream::msgType();
UPstream::msgType() = oldTag+1;
// Get the coupling information from the mappedPatchBase
const mappedPatchBase& mpp = refCast<const mappedPatchBase>
);
const polyMesh& nbrMesh = mpp.sampleMesh();
const fvPatch& nbrPatch = refCast<const fvMesh>
(
nbrMesh
).boundary()[mpp.samplePolyPatch().index()];
// Force recalculation of mapping and schedule
const mapDistribute& distMap = mpp.map();
tmp<scalarField> intFld = patchInternalField();
// Calculate the temperature by harmonic averaging
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& nbrField =
refCast
<
const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
>
(
nbrPatch.lookupPatchField<volScalarField, scalar>
neighbourFieldName_
)
);
// Swap to obtain full local values of neighbour internal field
scalarField nbrIntFld(nbrField.patchInternalField());
// Swap to obtain full local values of neighbour kappa*delta
scalarField nbrKDelta(nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs());
tmp<scalarField> myKDelta = kappa(*this)*patch().deltaCoeffs();
// Both sides agree on
// - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
// - gradient : (temperature-fld)*delta
// We've got a degree of freedom in how to implement this in a mixed bc.
// (what gradient, what fixedValue and mixing coefficient)
// Two reasonable choices:
// 1. specify above temperature on one side (preferentially the high side)
// and above gradient on the other. So this will switch between pure
// fixedvalue and pure fixedgradient
// 2. specify gradient and temperature such that the equations are the
// same on both sides. This leads to the choice of
// - refGradient = zero gradient
// - refValue = neighbour value
// - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
this->refValue() = nbrIntFld;
this->refGrad() = 0.0;
this->valueFraction() = nbrKDelta / (nbrKDelta + myKDelta());
mixedFvPatchScalarField::updateCoeffs();
if (debug)
{
scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
Info<< patch().boundaryMesh().mesh().name() << ':'
<< patch().name() << ':'
<< this->dimensionedInternalField().name() << " <- "
<< nbrMesh.name() << ':'
<< nbrPatch.name() << ':'
<< this->dimensionedInternalField().name() << " :"
<< " heat transfer rate:" << Q
<< " walltemperature "
<< " min:" << gMin(*this)
<< " max:" << gMax(*this)
<< " avg:" << gAverage(*this)
<< endl;
}
// Restore tag
UPstream::msgType() = oldTag;
}
void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
(
Ostream& os
) const
{
mixedFvPatchScalarField::write(os);
os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
<< token::END_STATEMENT << nl;
temperatureCoupledBase::write(os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField
(
fvPatchScalarField,
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace compressible
} // End namespace Foam
// ************************************************************************* //