This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/heatTransfer/chtMultiRegionFoam/derivedFvPatchFields/solidWallMixedTemperatureCoupled/solidWallMixedTemperatureCoupledFvPatchScalarField.C

270 lines
8.5 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2016-06-20 15:00:40 +00:00
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend 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
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "solidWallMixedTemperatureCoupledFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
2010-08-26 14:22:03 +00:00
#include "directMappedPatchBase.H"
#include "regionProperties.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
solidWallMixedTemperatureCoupledFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
mixedFvPatchScalarField(p, iF),
2010-08-26 14:22:03 +00:00
neighbourFieldName_("undefined-neighbourFieldName"),
2015-08-05 11:34:57 +00:00
KappaName_("undefined-Kappa")
{
this->refValue() = 0.0;
this->refGrad() = 0.0;
this->valueFraction() = 1.0;
}
Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
solidWallMixedTemperatureCoupledFvPatchScalarField
(
const solidWallMixedTemperatureCoupledFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
mixedFvPatchScalarField(ptf, p, iF, mapper),
2010-08-26 14:22:03 +00:00
neighbourFieldName_(ptf.neighbourFieldName_),
2015-08-05 11:34:57 +00:00
KappaName_(ptf.KappaName_)
{}
Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
solidWallMixedTemperatureCoupledFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
mixedFvPatchScalarField(p, iF),
2010-08-26 14:22:03 +00:00
neighbourFieldName_(dict.lookup("neighbourFieldName")),
2015-08-05 11:34:57 +00:00
KappaName_(dict.lookup("Kappa"))
{
2010-08-26 14:22:03 +00:00
if (!isA<directMappedPatchBase>(this->patch().patch()))
{
FatalErrorIn
(
"solidWallMixedTemperatureCoupledFvPatchScalarField::"
"solidWallMixedTemperatureCoupledFvPatchScalarField\n"
"(\n"
" const fvPatch& p,\n"
" const DimensionedField<scalar, volMesh>& iF,\n"
" const dictionary& dict\n"
")\n"
) << "\n patch type '" << p.type()
<< "' not type '" << directMappedPatchBase::typeName << "'"
<< "\n for patch " << p.name()
<< " of field " << dimensionedInternalField().name()
<< " in file " << dimensionedInternalField().objectPath()
<< exit(FatalError);
}
fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
2010-08-26 14:22:03 +00:00
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;
}
}
Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
solidWallMixedTemperatureCoupledFvPatchScalarField
(
const solidWallMixedTemperatureCoupledFvPatchScalarField& wtcsf,
const DimensionedField<scalar, volMesh>& iF
)
:
mixedFvPatchScalarField(wtcsf, iF),
2010-08-26 14:22:03 +00:00
neighbourFieldName_(wtcsf.neighbourFieldName_),
2015-08-05 11:34:57 +00:00
KappaName_(wtcsf.KappaName_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
const Foam::fvPatchScalarField&
2015-08-05 11:34:57 +00:00
Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::Kappa() const
{
2015-08-05 11:34:57 +00:00
return this->lookupPatchField<volScalarField, scalar>(KappaName_);
}
void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
2010-08-26 14:22:03 +00:00
// Get the coupling information from the directMappedPatchBase
const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
(
patch().patch()
);
const polyMesh& nbrMesh = mpp.sampleMesh();
const fvPatch& nbrPatch = refCast<const fvMesh>
(
nbrMesh
).boundary()[mpp.samplePolyPatch().index()];
2010-08-26 14:22:03 +00:00
// Force recalculation of mapping and schedule
const mapDistribute& distMap = mpp.map();
2010-08-26 14:22:03 +00:00
tmp<scalarField> intFld = patchInternalField();
2010-08-26 14:22:03 +00:00
const solidWallMixedTemperatureCoupledFvPatchScalarField& nbrField =
refCast<const solidWallMixedTemperatureCoupledFvPatchScalarField>
(
nbrPatch.lookupPatchField<volScalarField, scalar>
(
neighbourFieldName_
)
);
// Swap to obtain full local values of neighbour internal field
scalarField nbrIntFld = nbrField.patchInternalField();
mapDistribute::distribute
(
static_cast<Pstream::commsTypes>(Pstream::defaultCommsType()),
2010-08-26 14:22:03 +00:00
distMap.schedule(),
distMap.constructSize(),
distMap.subMap(), // what to send
distMap.constructMap(), // what to receive
nbrIntFld
);
2015-08-05 11:34:57 +00:00
// Swap to obtain full local values of neighbour Kappa*delta
scalarField nbrKappaDelta = nbrField.Kappa()*nbrPatch.deltaCoeffs();
2010-08-26 14:22:03 +00:00
mapDistribute::distribute
(
static_cast<Pstream::commsTypes>(Pstream::defaultCommsType()),
2010-08-26 14:22:03 +00:00
distMap.schedule(),
distMap.constructSize(),
distMap.subMap(), // what to send
distMap.constructMap(), // what to receive
2015-08-05 11:34:57 +00:00
nbrKappaDelta
2010-08-26 14:22:03 +00:00
);
2015-08-05 11:34:57 +00:00
tmp<scalarField> myKappaDelta = Kappa()*patch().deltaCoeffs();
2010-08-26 14:22:03 +00:00
// Both sides agree on
2015-08-05 11:34:57 +00:00
// - temperature : (myKappaDelta*fld + nbrKappaDelta*nbrFld)/(myKappaDelta+nbrKappaDelta)
2010-08-26 14:22:03 +00:00
// - 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
2015-08-05 11:34:57 +00:00
// - mixFraction = nbrKappaDelta / (nbrKappaDelta + myKappaDelta())
2010-08-26 14:22:03 +00:00
this->refValue() = nbrIntFld;
2010-08-26 14:22:03 +00:00
this->refGrad() = 0.0;
2015-08-05 11:34:57 +00:00
this->valueFraction() = nbrKappaDelta / (nbrKappaDelta + myKappaDelta());
mixedFvPatchScalarField::updateCoeffs();
2010-08-26 14:22:03 +00:00
if (debug)
{
2015-08-05 11:34:57 +00:00
scalar Q = gSum(Kappa()*patch().magSf()*snGrad());
2010-08-26 14:22:03 +00:00
Info<< patch().boundaryMesh().mesh().name() << ':'
<< patch().name() << ':'
<< this->dimensionedInternalField().name() << " -> "
<< nbrMesh.name() << ':'
<< nbrPatch.name() << ':'
<< this->dimensionedInternalField().name() << " :"
<< " heatFlux:" << Q
<< " walltemperature "
<< " min:" << gMin(*this)
<< " max:" << gMax(*this)
<< " avg:" << gAverage(*this)
<< endl;
}
}
void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::write
(
Ostream& os
) const
{
mixedFvPatchScalarField::write(os);
2010-08-26 14:22:03 +00:00
os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
<< token::END_STATEMENT << nl;
2015-08-05 11:34:57 +00:00
os.writeKeyword("Kappa") << KappaName_ << token::END_STATEMENT << nl;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
solidWallMixedTemperatureCoupledFvPatchScalarField
);
} // End namespace Foam
// ************************************************************************* //