1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "solidWallMixedTemperatureCoupledFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "directMappedPatchBase.H"
32 #include "regionProperties.H"
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 bool Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::interfaceOwner
38 const polyMesh& nbrRegion
41 const fvMesh& myRegion = patch().boundaryMesh().mesh();
43 const regionProperties& props =
44 myRegion.objectRegistry::parent().lookupObject<regionProperties>
49 label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
52 label i = findIndex(props.solidRegionNames(), myRegion.name());
58 "solidWallMixedTemperatureCoupledFvPatchScalarField"
59 "::interfaceOwner(const polyMesh&) const"
60 ) << "Cannot find region " << myRegion.name()
61 << " neither in fluids " << props.fluidRegionNames()
62 << " nor in solids " << props.solidRegionNames()
65 myIndex = props.fluidRegionNames().size() + i;
67 label nbrIndex = findIndex(props.fluidRegionNames(), nbrRegion.name());
70 label i = findIndex(props.solidRegionNames(), nbrRegion.name());
74 FatalErrorIn("coupleManager::interfaceOwner(const polyMesh&) const")
75 << "Cannot find region " << nbrRegion.name()
76 << " neither in fluids " << props.fluidRegionNames()
77 << " nor in solids " << props.solidRegionNames()
80 nbrIndex = props.fluidRegionNames().size() + i;
83 return myIndex < nbrIndex;
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
90 solidWallMixedTemperatureCoupledFvPatchScalarField
93 const DimensionedField<scalar, volMesh>& iF
96 mixedFvPatchScalarField(p, iF),
97 neighbourFieldName_("undefined-neighbourFieldName"),
100 this->refValue() = 0.0;
101 this->refGrad() = 0.0;
102 this->valueFraction() = 1.0;
103 this->fixesValue_ = true;
107 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
108 solidWallMixedTemperatureCoupledFvPatchScalarField
110 const solidWallMixedTemperatureCoupledFvPatchScalarField& ptf,
112 const DimensionedField<scalar, volMesh>& iF,
113 const fvPatchFieldMapper& mapper
116 mixedFvPatchScalarField(ptf, p, iF, mapper),
117 neighbourFieldName_(ptf.neighbourFieldName_),
119 fixesValue_(ptf.fixesValue_)
123 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
124 solidWallMixedTemperatureCoupledFvPatchScalarField
127 const DimensionedField<scalar, volMesh>& iF,
128 const dictionary& dict
131 mixedFvPatchScalarField(p, iF),
132 neighbourFieldName_(dict.lookup("neighbourFieldName")),
133 KName_(dict.lookup("K"))
135 if (!isA<directMappedPatchBase>(this->patch().patch()))
139 "solidWallMixedTemperatureCoupledFvPatchScalarField::"
140 "solidWallMixedTemperatureCoupledFvPatchScalarField\n"
142 " const fvPatch& p,\n"
143 " const DimensionedField<scalar, volMesh>& iF,\n"
144 " const dictionary& dict\n"
146 ) << "\n patch type '" << p.type()
147 << "' not type '" << directMappedPatchBase::typeName << "'"
148 << "\n for patch " << p.name()
149 << " of field " << dimensionedInternalField().name()
150 << " in file " << dimensionedInternalField().objectPath()
154 fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
156 if (dict.found("refValue"))
159 refValue() = scalarField("refValue", dict, p.size());
160 refGrad() = scalarField("refGradient", dict, p.size());
161 valueFraction() = scalarField("valueFraction", dict, p.size());
162 fixesValue_ = readBool(dict.lookup("fixesValue"));
166 // Start from user entered data. Assume fixedValue.
169 valueFraction() = 1.0;
175 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
176 solidWallMixedTemperatureCoupledFvPatchScalarField
178 const solidWallMixedTemperatureCoupledFvPatchScalarField& wtcsf,
179 const DimensionedField<scalar, volMesh>& iF
182 mixedFvPatchScalarField(wtcsf, iF),
183 neighbourFieldName_(wtcsf.neighbourFieldName_),
184 KName_(wtcsf.KName_),
185 fixesValue_(wtcsf.fixesValue_)
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 const Foam::fvPatchScalarField&
192 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::K() const
194 return this->patch().lookupPatchField<volScalarField, scalar>(KName_);
198 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::updateCoeffs()
205 // Get the coupling information from the directMappedPatchBase
206 const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
210 const polyMesh& nbrMesh = mpp.sampleMesh();
212 tmp<scalarField> intFld = patchInternalField();
214 if (interfaceOwner(nbrMesh))
216 // Note: other side information could be cached - it only needs
217 // to be updated the first time round the iteration (i.e. when
218 // switching regions) but unfortunately we don't have this information.
220 const mapDistribute& distMap = mpp.map();
221 const fvPatch& nbrPatch = refCast<const fvMesh>
224 ).boundary()[mpp.samplePolyPatch().index()];
227 // Calculate the temperature by harmonic averaging
228 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
230 const solidWallMixedTemperatureCoupledFvPatchScalarField& nbrField =
231 refCast<const solidWallMixedTemperatureCoupledFvPatchScalarField>
233 nbrPatch.lookupPatchField<volScalarField, scalar>
239 // Swap to obtain full local values of neighbour internal field
240 scalarField nbrIntFld = nbrField.patchInternalField();
241 mapDistribute::distribute
243 Pstream::defaultCommsType,
245 distMap.constructSize(),
246 distMap.subMap(), // what to send
247 distMap.constructMap(), // what to receive
251 // Swap to obtain full local values of neighbour K*delta
252 scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
253 mapDistribute::distribute
255 Pstream::defaultCommsType,
257 distMap.constructSize(),
258 distMap.subMap(), // what to send
259 distMap.constructMap(), // what to receive
264 tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();
266 // Calculate common wall temperature. Reuse *this to store common value.
269 (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
270 / (myKDelta() + nbrKDelta)
273 fvPatchScalarField::operator=(Twall);
274 // Distribute back and assign to neighbour
275 mapDistribute::distribute
277 Pstream::defaultCommsType,
280 distMap.constructMap(), // reverse : what to send
284 const_cast<solidWallMixedTemperatureCoupledFvPatchScalarField&>
287 ).fvPatchScalarField::operator=(Twall);
291 // Switch between fixed value (of harmonic avg) or gradient
292 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
296 // Like snGrad but bypass switching on refValue/refGrad.
297 tmp<scalarField> normalGradient = (*this-intFld())*patch().deltaCoeffs();
301 scalar Q = gSum(K()*patch().magSf()*normalGradient());
303 Info<< "solidWallMixedTemperatureCoupledFvPatchScalarField::"
304 << "updateCoeffs() :"
305 << " patch:" << patch().name()
307 << " walltemperature "
308 << " min:" << gMin(*this)
309 << " max:" << gMax(*this)
310 << " avg:" << gAverage(*this)
316 // if outgoing flux use fixed value.
317 if (normalGradient()[i] < 0.0)
319 this->refValue()[i] = operator[](i);
320 this->refGrad()[i] = 0.0; // not used
321 this->valueFraction()[i] = 1.0;
326 this->refValue()[i] = 0.0; // not used
327 this->refGrad()[i] = normalGradient()[i];
328 this->valueFraction()[i] = 0.0;
332 reduce(nFixed, sumOp<label>());
334 fixesValue_ = (nFixed > 0);
338 label nTotSize = returnReduce(this->size(), sumOp<label>());
340 Info<< "solidWallMixedTemperatureCoupledFvPatchScalarField::"
341 << "updateCoeffs() :"
342 << " patch:" << patch().name()
343 << " out of:" << nTotSize
344 << " fixedBC:" << nFixed
345 << " gradient:" << nTotSize-nFixed << endl;
348 mixedFvPatchScalarField::updateCoeffs();
352 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::write
357 mixedFvPatchScalarField::write(os);
358 os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
359 << token::END_STATEMENT << nl;
360 os.writeKeyword("K") << KName_ << token::END_STATEMENT << nl;
361 os.writeKeyword("fixesValue") << fixesValue_ << token::END_STATEMENT << nl;
365 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 solidWallMixedTemperatureCoupledFvPatchScalarField
376 } // End namespace Foam
378 // ************************************************************************* //