1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. 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"
32 #include "basicThermo.H"
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 bool Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::interfaceOwner
39 const polyMesh& nbrRegion,
40 const polyPatch& nbrPatch
43 const fvMesh& myRegion = patch().boundaryMesh().mesh();
45 if (nbrRegion.name() == myRegion.name())
47 return patch().index() < nbrPatch.index();
51 const regionProperties& props =
52 myRegion.objectRegistry::parent().lookupObject<regionProperties>
57 label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
60 label i = findIndex(props.solidRegionNames(), myRegion.name());
66 "turbulentTemperatureCoupledBaffleFvPatchScalarField"
67 "::interfaceOwner(const polyMesh&"
68 ", const polyPatch&)const"
69 ) << "Cannot find region " << myRegion.name()
70 << " neither in fluids " << props.fluidRegionNames()
71 << " nor in solids " << props.solidRegionNames()
74 myIndex = props.fluidRegionNames().size() + i;
76 label nbrIndex = findIndex
78 props.fluidRegionNames(),
83 label i = findIndex(props.solidRegionNames(), nbrRegion.name());
89 "coupleManager::interfaceOwner"
90 "(const polyMesh&, const polyPatch&) const"
91 ) << "Cannot find region " << nbrRegion.name()
92 << " neither in fluids " << props.fluidRegionNames()
93 << " nor in solids " << props.solidRegionNames()
96 nbrIndex = props.fluidRegionNames().size() + i;
99 return myIndex < nbrIndex;
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
106 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
107 turbulentTemperatureCoupledBaffleFvPatchScalarField
110 const DimensionedField<scalar, volMesh>& iF
113 fixedValueFvPatchScalarField(p, iF),
114 neighbourFieldName_("undefined-neighbourFieldName"),
115 KappaName_("undefined-Kappa")
119 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
120 turbulentTemperatureCoupledBaffleFvPatchScalarField
122 const turbulentTemperatureCoupledBaffleFvPatchScalarField& ptf,
124 const DimensionedField<scalar, volMesh>& iF,
125 const fvPatchFieldMapper& mapper
128 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
129 neighbourFieldName_(ptf.neighbourFieldName_),
130 KappaName_(ptf.KappaName_)
134 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
135 turbulentTemperatureCoupledBaffleFvPatchScalarField
138 const DimensionedField<scalar, volMesh>& iF,
139 const dictionary& dict
142 fixedValueFvPatchScalarField(p, iF, dict),
143 neighbourFieldName_(dict.lookup("neighbourFieldName")),
144 KappaName_(dict.lookup("Kappa"))
146 if (!isA<directMappedPatchBase>(this->patch().patch()))
150 "turbulentTemperatureCoupledBaffleFvPatchScalarField::"
151 "turbulentTemperatureCoupledBaffleFvPatchScalarField\n"
153 " const fvPatch& p,\n"
154 " const DimensionedField<scalar, volMesh>& iF,\n"
155 " const dictionary& dict\n"
157 ) << "\n patch type '" << p.type()
158 << "' not type '" << directMappedPatchBase::typeName << "'"
159 << "\n for patch " << p.name()
160 << " of field " << dimensionedInternalField().name()
161 << " in file " << dimensionedInternalField().objectPath()
167 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
168 turbulentTemperatureCoupledBaffleFvPatchScalarField
170 const turbulentTemperatureCoupledBaffleFvPatchScalarField& wtcsf,
171 const DimensionedField<scalar, volMesh>& iF
174 fixedValueFvPatchScalarField(wtcsf, iF),
175 neighbourFieldName_(wtcsf.neighbourFieldName_),
176 KappaName_(wtcsf.KappaName_)
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 Foam::tmp<Foam::scalarField>
183 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::Kappa() const
185 const fvMesh& mesh = patch().boundaryMesh().mesh();
187 if (KappaName_ == "none")
189 const compressible::RASModel& model =
190 db().lookupObject<compressible::RASModel>("RASProperties");
192 tmp<volScalarField> talpha = model.alphaEff();
194 const basicThermo& thermo =
195 db().lookupObject<basicThermo>("thermophysicalProperties");
198 talpha().boundaryField()[patch().index()]
199 *thermo.Cp()().boundaryField()[patch().index()];
201 else if (mesh.objectRegistry::foundObject<volScalarField>(KappaName_))
203 return lookupPatchField<volScalarField, scalar>(KappaName_);
205 else if (mesh.objectRegistry::foundObject<volSymmTensorField>(KappaName_))
207 const symmTensorField& KappaWall =
208 lookupPatchField<volSymmTensorField, scalar>(KappaName_);
210 vectorField n = patch().nf();
212 return n & KappaWall & n;
218 "turbulentTemperatureCoupledBaffleFvPatchScalarField::Kappa() const"
219 ) << "Did not find field " << KappaName_
220 << " on mesh " << mesh.name() << " patch " << patch().name()
222 << "Please set 'Kappa' to 'none', a valid volScalarField"
223 << " or a valid volSymmTensorField." << exit(FatalError);
225 return scalarField(0);
230 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::updateCoeffs()
237 // Get the coupling information from the directMappedPatchBase
238 const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
242 const polyMesh& nbrMesh = mpp.sampleMesh();
243 const fvPatch& nbrPatch = refCast<const fvMesh>
246 ).boundary()[mpp.samplePolyPatch().index()];
248 // Force recalculation of mapping and schedule
249 const mapDistribute& distMap = mpp.map();
250 (void)distMap.schedule();
252 tmp<scalarField> intFld = patchInternalField();
254 if (interfaceOwner(nbrMesh, nbrPatch.patch()))
256 // Note: other side information could be cached - it only needs
257 // to be updated the first time round the iteration (i.e. when
258 // switching regions) but unfortunately we don't have this information.
261 // Calculate the temperature by harmonic averaging
262 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
264 const turbulentTemperatureCoupledBaffleFvPatchScalarField& nbrField =
265 refCast<const turbulentTemperatureCoupledBaffleFvPatchScalarField>
267 nbrPatch.lookupPatchField<volScalarField, scalar>
273 // Swap to obtain full local values of neighbour internal field
274 scalarField nbrIntFld = nbrField.patchInternalField();
275 mapDistribute::distribute
277 Pstream::defaultComms(),
279 distMap.constructSize(),
280 distMap.subMap(), // what to send
281 distMap.constructMap(), // what to receive
285 // Swap to obtain full local values of neighbour Kappa*delta
286 scalarField nbrKappaDelta = nbrField.Kappa()*nbrPatch.deltaCoeffs();
287 mapDistribute::distribute
289 Pstream::defaultComms(),
291 distMap.constructSize(),
292 distMap.subMap(), // what to send
293 distMap.constructMap(), // what to receive
297 tmp<scalarField> myKappaDelta = Kappa()*patch().deltaCoeffs();
299 // Calculate common wall temperature.
300 // Reuse *this to store common value.
303 (myKappaDelta()*intFld() + nbrKappaDelta*nbrIntFld)
304 / (myKappaDelta() + nbrKappaDelta)
307 fvPatchScalarField::operator=(Twall);
308 // Distribute back and assign to neighbour
309 mapDistribute::distribute
311 Pstream::defaultComms(),
314 distMap.constructMap(), // reverse : what to send
318 const_cast<turbulentTemperatureCoupledBaffleFvPatchScalarField&>
321 ).fvPatchScalarField::operator=(Twall);
326 //tmp<scalarField> normalGradient =
328 // * patch().deltaCoeffs();
330 scalar Q = gSum(Kappa()*patch().magSf()*snGrad());
332 Info<< patch().boundaryMesh().mesh().name() << ':'
333 << patch().name() << ':'
334 << this->dimensionedInternalField().name() << " -> "
335 << nbrMesh.name() << ':'
336 << nbrPatch.name() << ':'
337 << this->dimensionedInternalField().name() << " :"
339 << " walltemperature "
340 << " min:" << gMin(*this)
341 << " max:" << gMax(*this)
342 << " avg:" << gAverage(*this)
346 fixedValueFvPatchScalarField::updateCoeffs();
350 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::write
355 fixedValueFvPatchScalarField::write(os);
356 os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
357 << token::END_STATEMENT << nl;
358 os.writeKeyword("Kappa") << KappaName_ << token::END_STATEMENT << nl;
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
370 turbulentTemperatureCoupledBaffleFvPatchScalarField
373 } // End namespace Foam
375 // ************************************************************************* //