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 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "directMappedPatchBase.H"
32 #include "regionProperties.H"
33 #include "basicThermo.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 namespace compressible
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 bool turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::interfaceOwner
47 const polyMesh& nbrRegion,
48 const polyPatch& nbrPatch
51 const fvMesh& myRegion = patch().boundaryMesh().mesh();
53 if (nbrRegion.name() == myRegion.name())
55 return patch().index() < nbrPatch.index();
59 const regionProperties& props =
60 myRegion.objectRegistry::parent().lookupObject<regionProperties>
65 label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
68 label i = findIndex(props.solidRegionNames(), myRegion.name());
74 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField"
75 "::interfaceOwner(const polyMesh&"
76 ", const polyPatch&)const"
77 ) << "Cannot find region " << myRegion.name()
78 << " neither in fluids " << props.fluidRegionNames()
79 << " nor in solids " << props.solidRegionNames()
82 myIndex = props.fluidRegionNames().size() + i;
84 label nbrIndex = findIndex
86 props.fluidRegionNames(),
91 label i = findIndex(props.solidRegionNames(), nbrRegion.name());
97 "coupleManager::interfaceOwner"
98 "(const polyMesh&, const polyPatch&) const"
99 ) << "Cannot find region " << nbrRegion.name()
100 << " neither in fluids " << props.fluidRegionNames()
101 << " nor in solids " << props.solidRegionNames()
104 nbrIndex = props.fluidRegionNames().size() + i;
106 return myIndex < nbrIndex;
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
114 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
117 const DimensionedField<scalar, volMesh>& iF
120 mixedFvPatchScalarField(p, iF),
121 neighbourFieldName_("undefined-neighbourFieldName"),
122 KName_("undefined-K")
124 this->refValue() = 0.0;
125 this->refGrad() = 0.0;
126 this->valueFraction() = 1.0;
127 this->fixesValue_ = true;
131 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
132 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
134 const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& ptf,
136 const DimensionedField<scalar, volMesh>& iF,
137 const fvPatchFieldMapper& mapper
140 mixedFvPatchScalarField(ptf, p, iF, mapper),
141 neighbourFieldName_(ptf.neighbourFieldName_),
143 fixesValue_(ptf.fixesValue_)
147 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
148 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
151 const DimensionedField<scalar, volMesh>& iF,
152 const dictionary& dict
155 mixedFvPatchScalarField(p, iF),
156 neighbourFieldName_(dict.lookup("neighbourFieldName")),
157 KName_(dict.lookup("K"))
159 if (!isA<directMappedPatchBase>(this->patch().patch()))
163 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::"
164 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField\n"
166 " const fvPatch& p,\n"
167 " const DimensionedField<scalar, volMesh>& iF,\n"
168 " const dictionary& dict\n"
170 ) << "\n patch type '" << p.type()
171 << "' not type '" << directMappedPatchBase::typeName << "'"
172 << "\n for patch " << p.name()
173 << " of field " << dimensionedInternalField().name()
174 << " in file " << dimensionedInternalField().objectPath()
178 fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
180 if (dict.found("refValue"))
183 refValue() = scalarField("refValue", dict, p.size());
184 refGrad() = scalarField("refGradient", dict, p.size());
185 valueFraction() = scalarField("valueFraction", dict, p.size());
186 fixesValue_ = readBool(dict.lookup("fixesValue"));
190 // Start from user entered data. Assume fixedValue.
193 valueFraction() = 1.0;
199 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
200 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
202 const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
203 const DimensionedField<scalar, volMesh>& iF
206 mixedFvPatchScalarField(wtcsf, iF),
207 neighbourFieldName_(wtcsf.neighbourFieldName_),
208 KName_(wtcsf.KName_),
209 fixesValue_(wtcsf.fixesValue_)
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::K() const
218 const fvMesh& mesh = patch().boundaryMesh().mesh();
220 if (KName_ == "none")
222 const compressible::RASModel& model =
223 db().lookupObject<compressible::RASModel>("RASProperties");
225 tmp<volScalarField> talpha = model.alphaEff();
227 const basicThermo& thermo =
228 db().lookupObject<basicThermo>("thermophysicalProperties");
231 talpha().boundaryField()[patch().index()]
232 *thermo.Cp()().boundaryField()[patch().index()];
234 else if (mesh.objectRegistry::foundObject<volScalarField>(KName_))
236 return patch().lookupPatchField<volScalarField, scalar>(KName_);
238 else if (mesh.objectRegistry::foundObject<volSymmTensorField>(KName_))
240 const symmTensorField& KWall =
241 patch().lookupPatchField<volSymmTensorField, scalar>(KName_);
243 vectorField n = patch().nf();
245 return n & KWall & n;
251 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::K()"
253 ) << "Did not find field " << KName_
254 << " on mesh " << mesh.name() << " patch " << patch().name()
256 << "Please set 'K' to 'none', a valid volScalarField"
257 << " or a valid volSymmTensorField." << exit(FatalError);
259 return scalarField(0);
264 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs()
271 // Get the coupling information from the directMappedPatchBase
272 const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
276 const polyMesh& nbrMesh = mpp.sampleMesh();
277 const fvPatch& nbrPatch = refCast<const fvMesh>
280 ).boundary()[mpp.samplePolyPatch().index()];
282 // Force recalculation of mapping and schedule
283 const mapDistribute& distMap = mpp.map();
284 (void)distMap.schedule();
286 tmp<scalarField> intFld = patchInternalField();
288 if (interfaceOwner(nbrMesh, nbrPatch.patch()))
290 // Note: other side information could be cached - it only needs
291 // to be updated the first time round the iteration (i.e. when
292 // switching regions) but unfortunately we don't have this information.
295 // Calculate the temperature by harmonic averaging
296 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
298 const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField&
302 const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
305 nbrPatch.lookupPatchField<volScalarField, scalar>
311 // Swap to obtain full local values of neighbour internal field
312 scalarField nbrIntFld = nbrField.patchInternalField();
313 mapDistribute::distribute
315 Pstream::defaultCommsType,
317 distMap.constructSize(),
318 distMap.subMap(), // what to send
319 distMap.constructMap(), // what to receive
323 // Swap to obtain full local values of neighbour K*delta
324 scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
325 mapDistribute::distribute
327 Pstream::defaultCommsType,
329 distMap.constructSize(),
330 distMap.subMap(), // what to send
331 distMap.constructMap(), // what to receive
335 tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();
337 // Calculate common wall temperature. Reuse *this to store common value.
340 (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
341 / (myKDelta() + nbrKDelta)
344 fvPatchScalarField::operator=(Twall);
345 // Distribute back and assign to neighbour
346 mapDistribute::distribute
348 Pstream::defaultCommsType,
351 distMap.constructMap(), // reverse : what to send
355 const_cast<turbulentTemperatureCoupledBaffleMixedFvPatchScalarField&>
358 ).fvPatchScalarField::operator=(Twall);
362 // Switch between fixed value (of harmonic avg) or gradient
363 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
367 // Like snGrad but bypass switching on refValue/refGrad.
368 tmp<scalarField> normalGradient = (*this-intFld())*patch().deltaCoeffs();
372 scalar Q = gSum(K()*patch().magSf()*normalGradient());
374 Info<< patch().boundaryMesh().mesh().name() << ':'
375 << patch().name() << ':'
376 << this->dimensionedInternalField().name() << " -> "
377 << nbrMesh.name() << ':'
378 << nbrPatch.name() << ':'
379 << this->dimensionedInternalField().name() << " :"
381 << " walltemperature "
382 << " min:" << gMin(*this)
383 << " max:" << gMax(*this)
384 << " avg:" << gAverage(*this)
390 // if outgoing flux use fixed value.
391 if (normalGradient()[i] < 0.0)
393 this->refValue()[i] = operator[](i);
394 this->refGrad()[i] = 0.0; // not used by me
395 this->valueFraction()[i] = 1.0;
400 // Fixed gradient. Make sure to have valid refValue (even though
401 // I am not using it - other boundary conditions might)
402 this->refValue()[i] = operator[](i);
403 this->refGrad()[i] = normalGradient()[i];
404 this->valueFraction()[i] = 0.0;
408 reduce(nFixed, sumOp<label>());
410 fixesValue_ = (nFixed > 0);
414 label nTotSize = returnReduce(this->size(), sumOp<label>());
416 Info<< patch().boundaryMesh().mesh().name() << ':'
417 << patch().name() << ':'
418 << this->dimensionedInternalField().name() << " -> "
419 << nbrMesh.name() << ':'
420 << nbrPatch.name() << ':'
421 << this->dimensionedInternalField().name() << " :"
422 << " patch:" << patch().name()
423 << " out of:" << nTotSize
424 << " fixedBC:" << nFixed
425 << " gradient:" << nTotSize-nFixed << endl;
428 mixedFvPatchScalarField::updateCoeffs();
432 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
437 mixedFvPatchScalarField::write(os);
438 os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
439 << token::END_STATEMENT << nl;
440 os.writeKeyword("K") << KName_ << token::END_STATEMENT << nl;
441 os.writeKeyword("fixesValue") << fixesValue_ << token::END_STATEMENT << nl;
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
450 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
454 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
456 } // End namespace compressible
457 } // End namespace Foam
460 // ************************************************************************* //