BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / turbulenceModel / derivedFvPatchFields / turbulentTemperatureCoupledBaffle / turbulentTemperatureCoupledBaffleFvPatchScalarField.C
blob012d0130121d8706eb2bc1bb28f8ce01f3fea7fb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
39 ) const
41     const fvMesh& myRegion = patch().boundaryMesh().mesh();
43     if (nbrRegion.name() == myRegion.name())
44     {
45         return patch().index() < nbrPatch.index();
46     }
47     else
48     {
49         const regionProperties& props =
50             myRegion.objectRegistry::parent().lookupObject<regionProperties>
51             (
52                 "regionProperties"
53             );
55         label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
56         if (myIndex == -1)
57         {
58             label i = findIndex(props.solidRegionNames(), myRegion.name());
60             if (i == -1)
61             {
62                 FatalErrorIn
63                 (
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()
70                     << exit(FatalError);
71             }
72             myIndex = props.fluidRegionNames().size() + i;
73         }
74         label nbrIndex = findIndex
75         (
76             props.fluidRegionNames(),
77             nbrRegion.name()
78         );
79         if (nbrIndex == -1)
80         {
81             label i = findIndex(props.solidRegionNames(), nbrRegion.name());
83             if (i == -1)
84             {
85                 FatalErrorIn
86                 (
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()
92                     << exit(FatalError);
93             }
94             nbrIndex = props.fluidRegionNames().size() + i;
95         }
97         return myIndex < nbrIndex;
98     }
102 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
104 Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::
105 turbulentTemperatureCoupledBaffleFvPatchScalarField
107     const fvPatch& p,
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,
121     const fvPatch& p,
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
135     const fvPatch& p,
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()))
145     {
146         FatalErrorIn
147         (
148             "turbulentTemperatureCoupledBaffleFvPatchScalarField::"
149             "turbulentTemperatureCoupledBaffleFvPatchScalarField\n"
150             "(\n"
151             "    const fvPatch& p,\n"
152             "    const DimensionedField<scalar, volMesh>& iF,\n"
153             "    const dictionary& dict\n"
154             ")\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()
160             << exit(FatalError);
161     }
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()
182     if (updated())
183     {
184         return;
185     }
187     // Since we're inside initEvaluate/evaluate there might be processor
188     // comms underway. Change the tag we use.
189     int oldTag = UPstream::msgType();
190     UPstream::msgType() = oldTag+1;
192     // Get the coupling information from the directMappedPatchBase
193     const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
194     (
195         this->patch().patch()
196     );
197     const polyMesh& nbrMesh = mpp.sampleMesh();
198     const fvPatch& nbrPatch = refCast<const fvMesh>
199     (
200         nbrMesh
201     ).boundary()[mpp.samplePolyPatch().index()];
203     // Force recalculation of mapping and schedule
204     const mapDistribute& distMap = mpp.map();
205     (void)distMap.schedule();
207     tmp<scalarField> intFld = patchInternalField();
209     if (interfaceOwner(nbrMesh, nbrPatch.patch()))
210     {
211         // Note: other side information could be cached - it only needs
212         // to be updated the first time round the iteration (i.e. when
213         // switching regions) but unfortunately we don't have this information.
216         // Calculate the temperature by harmonic averaging
217         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
219         const turbulentTemperatureCoupledBaffleFvPatchScalarField& nbrField =
220         refCast<const turbulentTemperatureCoupledBaffleFvPatchScalarField>
221         (
222             nbrPatch.lookupPatchField<volScalarField, scalar>
223             (
224                 neighbourFieldName_
225             )
226         );
228         // Swap to obtain full local values of neighbour internal field
229         scalarField nbrIntFld(nbrField.patchInternalField());
230         distMap.distribute(nbrIntFld);
232         // Swap to obtain full local values of neighbour K*delta
233         scalarField nbrKDelta(nbrField.K(nbrField)*nbrPatch.deltaCoeffs());
234         distMap.distribute(nbrKDelta);
236         tmp<scalarField> myKDelta = K(*this)*patch().deltaCoeffs();
238         // Calculate common wall temperature. Reuse *this to store common value.
239         scalarField Twall
240         (
241             (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
242           / (myKDelta() + nbrKDelta)
243         );
244         // Assign to me
245         fvPatchScalarField::operator=(Twall);
246         // Distribute back and assign to neighbour
247         distMap.reverseDistribute(nbrField.size(), Twall);
248         const_cast<turbulentTemperatureCoupledBaffleFvPatchScalarField&>
249         (
250             nbrField
251         ).fvPatchScalarField::operator=(Twall);
252     }
254     if (debug)
255     {
256         //tmp<scalarField> normalGradient =
257         //    (*this-intFld())
258         //  * patch().deltaCoeffs();
260         scalar Q = gSum(K(*this)*patch().magSf()*snGrad());
262         Info<< patch().boundaryMesh().mesh().name() << ':'
263             << patch().name() << ':'
264             << this->dimensionedInternalField().name() << " <- "
265             << nbrMesh.name() << ':'
266             << nbrPatch.name() << ':'
267             << this->dimensionedInternalField().name() << " :"
268             << " heat[W]:" << Q
269             << " walltemperature "
270             << " min:" << gMin(*this)
271             << " max:" << gMax(*this)
272             << " avg:" << gAverage(*this)
273             << endl;
274     }
276     // Restore tag
277     UPstream::msgType() = oldTag;
279     fixedValueFvPatchScalarField::updateCoeffs();
283 void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::write
285     Ostream& os
286 ) const
288     fixedValueFvPatchScalarField::write(os);
289     os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
290         << token::END_STATEMENT << nl;
291     temperatureCoupledBase::write(os);
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 namespace Foam
300 makePatchTypeField
302     fvPatchScalarField,
303     turbulentTemperatureCoupledBaffleFvPatchScalarField
306 } // End namespace Foam
308 // ************************************************************************* //