1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "turbulentTemperatureCoupledBaffleFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "directMappedPatchBase.H"
31 #include "regionProperties.H"
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 bool Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::interfaceOwner
37 const polyMesh& nbrRegion,
38 const polyPatch& nbrPatch
41 const fvMesh& myRegion = patch().boundaryMesh().mesh();
43 if (nbrRegion.name() == myRegion.name())
45 return patch().index() < nbrPatch.index();
49 const regionProperties& props =
50 myRegion.objectRegistry::parent().lookupObject<regionProperties>
55 label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
58 label i = findIndex(props.solidRegionNames(), myRegion.name());
64 "turbulentTemperatureCoupledBaffleFvPatchScalarField"
65 "::interfaceOwner(const polyMesh&"
66 ", const polyPatch&)const"
67 ) << "Cannot find region " << myRegion.name()
68 << " neither in fluids " << props.fluidRegionNames()
69 << " nor in solids " << props.solidRegionNames()
72 myIndex = props.fluidRegionNames().size() + i;
74 label nbrIndex = findIndex
76 props.fluidRegionNames(),
81 label i = findIndex(props.solidRegionNames(), nbrRegion.name());
87 "coupleManager::interfaceOwner"
88 "(const polyMesh&, const polyPatch&) const"
89 ) << "Cannot find region " << nbrRegion.name()
90 << " neither in fluids " << props.fluidRegionNames()
91 << " nor in solids " << props.solidRegionNames()
94 nbrIndex = props.fluidRegionNames().size() + i;
97 return myIndex < nbrIndex;
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
104 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
105 turbulentTemperatureCoupledBaffleFvPatchScalarField
108 const DimensionedField<scalar, volMesh>& iF
111 fixedValueFvPatchScalarField(p, iF),
112 temperatureCoupledBase(patch(), "undefined", "undefined-K"),
113 neighbourFieldName_("undefined-neighbourFieldName")
117 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
118 turbulentTemperatureCoupledBaffleFvPatchScalarField
120 const turbulentTemperatureCoupledBaffleFvPatchScalarField& ptf,
122 const DimensionedField<scalar, volMesh>& iF,
123 const fvPatchFieldMapper& mapper
126 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
127 temperatureCoupledBase(patch(), ptf.KMethod(), ptf.KName()),
128 neighbourFieldName_(ptf.neighbourFieldName_)
132 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
133 turbulentTemperatureCoupledBaffleFvPatchScalarField
136 const DimensionedField<scalar, volMesh>& iF,
137 const dictionary& dict
140 fixedValueFvPatchScalarField(p, iF, dict),
141 temperatureCoupledBase(patch(), dict),
142 neighbourFieldName_(dict.lookup("neighbourFieldName"))
144 if (!isA<directMappedPatchBase>(this->patch().patch()))
148 "turbulentTemperatureCoupledBaffleFvPatchScalarField::"
149 "turbulentTemperatureCoupledBaffleFvPatchScalarField\n"
151 " const fvPatch& p,\n"
152 " const DimensionedField<scalar, volMesh>& iF,\n"
153 " const dictionary& dict\n"
155 ) << "\n patch type '" << p.type()
156 << "' not type '" << directMappedPatchBase::typeName << "'"
157 << "\n for patch " << p.name()
158 << " of field " << dimensionedInternalField().name()
159 << " in file " << dimensionedInternalField().objectPath()
165 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
166 turbulentTemperatureCoupledBaffleFvPatchScalarField
168 const turbulentTemperatureCoupledBaffleFvPatchScalarField& wtcsf,
169 const DimensionedField<scalar, volMesh>& iF
172 fixedValueFvPatchScalarField(wtcsf, iF),
173 temperatureCoupledBase(patch(), wtcsf.KMethod(), wtcsf.KName()),
174 neighbourFieldName_(wtcsf.neighbourFieldName_)
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::updateCoeffs()
187 // Get the coupling information from the directMappedPatchBase
188 const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
190 this->patch().patch()
192 const polyMesh& nbrMesh = mpp.sampleMesh();
193 const fvPatch& nbrPatch = refCast<const fvMesh>
196 ).boundary()[mpp.samplePolyPatch().index()];
198 // Force recalculation of mapping and schedule
199 const mapDistribute& distMap = mpp.map();
200 (void)distMap.schedule();
202 tmp<scalarField> intFld = patchInternalField();
204 if (interfaceOwner(nbrMesh, nbrPatch.patch()))
206 // Note: other side information could be cached - it only needs
207 // to be updated the first time round the iteration (i.e. when
208 // switching regions) but unfortunately we don't have this information.
211 // Calculate the temperature by harmonic averaging
212 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
214 const turbulentTemperatureCoupledBaffleFvPatchScalarField& nbrField =
215 refCast<const turbulentTemperatureCoupledBaffleFvPatchScalarField>
217 nbrPatch.lookupPatchField<volScalarField, scalar>
223 // Swap to obtain full local values of neighbour internal field
224 scalarField nbrIntFld(nbrField.patchInternalField());
225 distMap.distribute(nbrIntFld);
227 // Swap to obtain full local values of neighbour K*delta
228 scalarField nbrKDelta(nbrField.K(nbrField)*nbrPatch.deltaCoeffs());
229 distMap.distribute(nbrKDelta);
231 tmp<scalarField> myKDelta = K(*this)*patch().deltaCoeffs();
233 // Calculate common wall temperature. Reuse *this to store common value.
236 (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
237 / (myKDelta() + nbrKDelta)
240 fvPatchScalarField::operator=(Twall);
241 // Distribute back and assign to neighbour
242 distMap.reverseDistribute(nbrField.size(), Twall);
243 const_cast<turbulentTemperatureCoupledBaffleFvPatchScalarField&>
246 ).fvPatchScalarField::operator=(Twall);
251 //tmp<scalarField> normalGradient =
253 // * patch().deltaCoeffs();
255 scalar Q = gSum(K(*this)*patch().magSf()*snGrad());
257 Info<< patch().boundaryMesh().mesh().name() << ':'
258 << patch().name() << ':'
259 << this->dimensionedInternalField().name() << " <- "
260 << nbrMesh.name() << ':'
261 << nbrPatch.name() << ':'
262 << this->dimensionedInternalField().name() << " :"
264 << " walltemperature "
265 << " min:" << gMin(*this)
266 << " max:" << gMax(*this)
267 << " avg:" << gAverage(*this)
271 fixedValueFvPatchScalarField::updateCoeffs();
275 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::write
280 fixedValueFvPatchScalarField::write(os);
281 os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
282 << token::END_STATEMENT << nl;
283 temperatureCoupledBase::write(os);
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 turbulentTemperatureCoupledBaffleFvPatchScalarField
298 } // End namespace Foam
300 // ************************************************************************* //