Bugfix: Pstream default comms format
[foam-extend-3.2.git] / src / turbulenceModels / compressible / RAS / derivedFvPatchFields / turbulentTemperatureCoupledBaffle / turbulentTemperatureCoupledBaffleFvPatchScalarField.C
blob4418e1aa16d37f7faf31eaea9494b7e6fdc50a3d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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"
33 #include "RASModel.H"
35 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
37 bool Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::interfaceOwner
39     const polyMesh& nbrRegion,
40     const polyPatch& nbrPatch
41 ) const
43     const fvMesh& myRegion = patch().boundaryMesh().mesh();
45     if (nbrRegion.name() == myRegion.name())
46     {
47         return patch().index() < nbrPatch.index();
48     }
49     else
50     {
51         const regionProperties& props =
52             myRegion.objectRegistry::parent().lookupObject<regionProperties>
53             (
54                 "regionProperties"
55             );
57         label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
58         if (myIndex == -1)
59         {
60             label i = findIndex(props.solidRegionNames(), myRegion.name());
62             if (i == -1)
63             {
64                 FatalErrorIn
65                 (
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()
72                     << exit(FatalError);
73             }
74             myIndex = props.fluidRegionNames().size() + i;
75         }
76         label nbrIndex = findIndex
77         (
78             props.fluidRegionNames(),
79             nbrRegion.name()
80         );
81         if (nbrIndex == -1)
82         {
83             label i = findIndex(props.solidRegionNames(), nbrRegion.name());
85             if (i == -1)
86             {
87                 FatalErrorIn
88                 (
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()
94                     << exit(FatalError);
95             }
96             nbrIndex = props.fluidRegionNames().size() + i;
97         }
99         return myIndex < nbrIndex;
100     }
104 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
106 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
107 turbulentTemperatureCoupledBaffleFvPatchScalarField
109     const fvPatch& p,
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,
123     const fvPatch& p,
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
137     const fvPatch& p,
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()))
147     {
148         FatalErrorIn
149         (
150             "turbulentTemperatureCoupledBaffleFvPatchScalarField::"
151             "turbulentTemperatureCoupledBaffleFvPatchScalarField\n"
152             "(\n"
153             "    const fvPatch& p,\n"
154             "    const DimensionedField<scalar, volMesh>& iF,\n"
155             "    const dictionary& dict\n"
156             ")\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()
162             << exit(FatalError);
163     }
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")
188     {
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");
197         return
198             talpha().boundaryField()[patch().index()]
199            *thermo.Cp()().boundaryField()[patch().index()];
200     }
201     else if (mesh.objectRegistry::foundObject<volScalarField>(KappaName_))
202     {
203         return lookupPatchField<volScalarField, scalar>(KappaName_);
204     }
205     else if (mesh.objectRegistry::foundObject<volSymmTensorField>(KappaName_))
206     {
207         const symmTensorField& KappaWall =
208             lookupPatchField<volSymmTensorField, scalar>(KappaName_);
210         vectorField n = patch().nf();
212         return n & KappaWall & n;
213     }
214     else
215     {
216         FatalErrorIn
217         (
218             "turbulentTemperatureCoupledBaffleFvPatchScalarField::Kappa() const"
219         )   << "Did not find field " << KappaName_
220             << " on mesh " << mesh.name() << " patch " << patch().name()
221             << endl
222             << "Please set 'Kappa' to 'none', a valid volScalarField"
223             << " or a valid volSymmTensorField." << exit(FatalError);
225         return scalarField(0);
226     }
230 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::updateCoeffs()
232     if (updated())
233     {
234         return;
235     }
237     // Get the coupling information from the directMappedPatchBase
238     const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
239     (
240         patch().patch()
241     );
242     const polyMesh& nbrMesh = mpp.sampleMesh();
243     const fvPatch& nbrPatch = refCast<const fvMesh>
244     (
245         nbrMesh
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()))
255     {
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>
266         (
267             nbrPatch.lookupPatchField<volScalarField, scalar>
268             (
269                 neighbourFieldName_
270             )
271         );
273         // Swap to obtain full local values of neighbour internal field
274         scalarField nbrIntFld = nbrField.patchInternalField();
275         mapDistribute::distribute
276         (
277             Pstream::defaultComms(),
278             distMap.schedule(),
279             distMap.constructSize(),
280             distMap.subMap(),           // what to send
281             distMap.constructMap(),     // what to receive
282             nbrIntFld
283         );
285         // Swap to obtain full local values of neighbour Kappa*delta
286         scalarField nbrKappaDelta = nbrField.Kappa()*nbrPatch.deltaCoeffs();
287         mapDistribute::distribute
288         (
289             Pstream::defaultComms(),
290             distMap.schedule(),
291             distMap.constructSize(),
292             distMap.subMap(),           // what to send
293             distMap.constructMap(),     // what to receive
294             nbrKappaDelta
295         );
297         tmp<scalarField> myKappaDelta = Kappa()*patch().deltaCoeffs();
299         // Calculate common wall temperature.
300         // Reuse *this to store common value.
301         scalarField Twall
302         (
303             (myKappaDelta()*intFld() + nbrKappaDelta*nbrIntFld)
304           / (myKappaDelta() + nbrKappaDelta)
305         );
306         // Assign to me
307         fvPatchScalarField::operator=(Twall);
308         // Distribute back and assign to neighbour
309         mapDistribute::distribute
310         (
311             Pstream::defaultComms(),
312             distMap.schedule(),
313             nbrField.size(),
314             distMap.constructMap(),     // reverse : what to send
315             distMap.subMap(),
316             Twall
317         );
318         const_cast<turbulentTemperatureCoupledBaffleFvPatchScalarField&>
319         (
320             nbrField
321         ).fvPatchScalarField::operator=(Twall);
322     }
324     if (debug)
325     {
326         //tmp<scalarField> normalGradient =
327         //    (*this-intFld())
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() << " :"
338             << " heatFlux:" << Q
339             << " walltemperature "
340             << " min:" << gMin(*this)
341             << " max:" << gMax(*this)
342             << " avg:" << gAverage(*this)
343             << endl;
344     }
346     fixedValueFvPatchScalarField::updateCoeffs();
350 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::write
352     Ostream& os
353 ) const
355     fixedValueFvPatchScalarField::write(os);
356     os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
357         << token::END_STATEMENT << nl;
358     os.writeKeyword("Kappa") << KappaName_ << token::END_STATEMENT << nl;
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 namespace Foam
367 makePatchTypeField
369     fvPatchScalarField,
370     turbulentTemperatureCoupledBaffleFvPatchScalarField
373 } // End namespace Foam
375 // ************************************************************************* //