/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 4.0 \\ / A nd | Web: http://www.foam-extend.org \\/ M anipulation | For copyright notice see file Copyright ------------------------------------------------------------------------------- License This file is part of foam-extend. 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 Free Software Foundation, either version 3 of the License, or (at your option) any later version. 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 along with foam-extend. If not, see . \*---------------------------------------------------------------------------*/ #include "solidWallMixedTemperatureCoupledFvPatchScalarField.H" #include "addToRunTimeSelectionTable.H" #include "fvPatchFieldMapper.H" #include "volFields.H" #include "directMappedPatchBase.H" #include "regionProperties.H" // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::solidWallMixedTemperatureCoupledFvPatchScalarField:: solidWallMixedTemperatureCoupledFvPatchScalarField ( const fvPatch& p, const DimensionedField& iF ) : mixedFvPatchScalarField(p, iF), neighbourFieldName_("undefined-neighbourFieldName"), 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& iF, const fvPatchFieldMapper& mapper ) : mixedFvPatchScalarField(ptf, p, iF, mapper), neighbourFieldName_(ptf.neighbourFieldName_), KappaName_(ptf.KappaName_) {} Foam::solidWallMixedTemperatureCoupledFvPatchScalarField:: solidWallMixedTemperatureCoupledFvPatchScalarField ( const fvPatch& p, const DimensionedField& iF, const dictionary& dict ) : mixedFvPatchScalarField(p, iF), neighbourFieldName_(dict.lookup("neighbourFieldName")), KappaName_(dict.lookup("Kappa")) { if (!isA(this->patch().patch())) { FatalErrorIn ( "solidWallMixedTemperatureCoupledFvPatchScalarField::" "solidWallMixedTemperatureCoupledFvPatchScalarField\n" "(\n" " const fvPatch& p,\n" " const DimensionedField& 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())); 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& iF ) : mixedFvPatchScalarField(wtcsf, iF), neighbourFieldName_(wtcsf.neighbourFieldName_), KappaName_(wtcsf.KappaName_) {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // const Foam::fvPatchScalarField& Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::Kappa() const { return this->lookupPatchField(KappaName_); } void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::updateCoeffs() { if (updated()) { return; } // Get the coupling information from the directMappedPatchBase const directMappedPatchBase& mpp = refCast ( patch().patch() ); const polyMesh& nbrMesh = mpp.sampleMesh(); const fvPatch& nbrPatch = refCast ( nbrMesh ).boundary()[mpp.samplePolyPatch().index()]; // Force recalculation of mapping and schedule const mapDistribute& distMap = mpp.map(); tmp intFld = patchInternalField(); const solidWallMixedTemperatureCoupledFvPatchScalarField& nbrField = refCast ( nbrPatch.lookupPatchField ( neighbourFieldName_ ) ); // Swap to obtain full local values of neighbour internal field scalarField nbrIntFld = nbrField.patchInternalField(); mapDistribute::distribute ( static_cast(Pstream::defaultCommsType()), distMap.schedule(), distMap.constructSize(), distMap.subMap(), // what to send distMap.constructMap(), // what to receive nbrIntFld ); // Swap to obtain full local values of neighbour Kappa*delta scalarField nbrKappaDelta = nbrField.Kappa()*nbrPatch.deltaCoeffs(); mapDistribute::distribute ( static_cast(Pstream::defaultCommsType()), distMap.schedule(), distMap.constructSize(), distMap.subMap(), // what to send distMap.constructMap(), // what to receive nbrKappaDelta ); tmp myKappaDelta = Kappa()*patch().deltaCoeffs(); // Both sides agree on // - temperature : (myKappaDelta*fld + nbrKappaDelta*nbrFld)/(myKappaDelta+nbrKappaDelta) // - 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 = nbrKappaDelta / (nbrKappaDelta + myKappaDelta()) this->refValue() = nbrIntFld; this->refGrad() = 0.0; this->valueFraction() = nbrKappaDelta / (nbrKappaDelta + myKappaDelta()); mixedFvPatchScalarField::updateCoeffs(); if (debug) { scalar Q = gSum(Kappa()*patch().magSf()*snGrad()); 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); os.writeKeyword("neighbourFieldName")<< neighbourFieldName_ << token::END_STATEMENT << nl; os.writeKeyword("Kappa") << KappaName_ << token::END_STATEMENT << nl; } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { makePatchTypeField ( fvPatchScalarField, solidWallMixedTemperatureCoupledFvPatchScalarField ); } // End namespace Foam // ************************************************************************* //