BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / turbulenceModel / derivedFvPatchFields / turbulentTemperatureCoupledBaffleMixed / turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C
blob28ca25cda74b6564d18f358a6ca15d83d8a5df3d
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 "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "directMappedPatchBase.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
36 namespace compressible
39 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
42 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
44 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
45 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
47     const fvPatch& p,
48     const DimensionedField<scalar, volMesh>& iF
51     mixedFvPatchScalarField(p, iF),
52     temperatureCoupledBase(patch(), "undefined", "undefined-K"),
53     neighbourFieldName_("undefined-neighbourFieldName")
55     this->refValue() = 0.0;
56     this->refGrad() = 0.0;
57     this->valueFraction() = 1.0;
61 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
62 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
64     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& ptf,
65     const fvPatch& p,
66     const DimensionedField<scalar, volMesh>& iF,
67     const fvPatchFieldMapper& mapper
70     mixedFvPatchScalarField(ptf, p, iF, mapper),
71     temperatureCoupledBase(patch(), ptf.KMethod(), ptf.KName()),
72     neighbourFieldName_(ptf.neighbourFieldName_)
76 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
77 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
79     const fvPatch& p,
80     const DimensionedField<scalar, volMesh>& iF,
81     const dictionary& dict
84     mixedFvPatchScalarField(p, iF),
85     temperatureCoupledBase(patch(), dict),
86     neighbourFieldName_(dict.lookup("neighbourFieldName"))
88     if (!isA<directMappedPatchBase>(this->patch().patch()))
89     {
90         FatalErrorIn
91         (
92             "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::"
93             "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField\n"
94             "(\n"
95             "    const fvPatch& p,\n"
96             "    const DimensionedField<scalar, volMesh>& iF,\n"
97             "    const dictionary& dict\n"
98             ")\n"
99         )   << "\n    patch type '" << p.type()
100             << "' not type '" << directMappedPatchBase::typeName << "'"
101             << "\n    for patch " << p.name()
102             << " of field " << dimensionedInternalField().name()
103             << " in file " << dimensionedInternalField().objectPath()
104             << exit(FatalError);
105     }
107     fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
109     if (dict.found("refValue"))
110     {
111         // Full restart
112         refValue() = scalarField("refValue", dict, p.size());
113         refGrad() = scalarField("refGradient", dict, p.size());
114         valueFraction() = scalarField("valueFraction", dict, p.size());
115     }
116     else
117     {
118         // Start from user entered data. Assume fixedValue.
119         refValue() = *this;
120         refGrad() = 0.0;
121         valueFraction() = 1.0;
122     }
126 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
127 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
129     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
130     const DimensionedField<scalar, volMesh>& iF
133     mixedFvPatchScalarField(wtcsf, iF),
134     temperatureCoupledBase(patch(), wtcsf.KMethod(), wtcsf.KName()),
135     neighbourFieldName_(wtcsf.neighbourFieldName_)
139 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
141 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs()
143     if (updated())
144     {
145         return;
146     }
148     // Since we're inside initEvaluate/evaluate there might be processor
149     // comms underway. Change the tag we use.
150     int oldTag = UPstream::msgType();
151     UPstream::msgType() = oldTag+1;
153     // Get the coupling information from the directMappedPatchBase
154     const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
155     (
156         patch().patch()
157     );
158     const polyMesh& nbrMesh = mpp.sampleMesh();
159     const fvPatch& nbrPatch = refCast<const fvMesh>
160     (
161         nbrMesh
162     ).boundary()[mpp.samplePolyPatch().index()];
164     // Force recalculation of mapping and schedule
165     const mapDistribute& distMap = mpp.map();
167     tmp<scalarField> intFld = patchInternalField();
170     // Calculate the temperature by harmonic averaging
171     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
173     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& nbrField =
174     refCast
175     <
176         const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
177     >
178     (
179         nbrPatch.lookupPatchField<volScalarField, scalar>
180         (
181             neighbourFieldName_
182         )
183     );
185     // Swap to obtain full local values of neighbour internal field
186     scalarField nbrIntFld(nbrField.patchInternalField());
187     distMap.distribute(nbrIntFld);
189     // Swap to obtain full local values of neighbour K*delta
190     scalarField nbrKDelta(nbrField.K(nbrField)*nbrPatch.deltaCoeffs());
191     distMap.distribute(nbrKDelta);
193     tmp<scalarField> myKDelta = K(*this)*patch().deltaCoeffs();
196     // Both sides agree on
197     // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
198     // - gradient    : (temperature-fld)*delta
199     // We've got a degree of freedom in how to implement this in a mixed bc.
200     // (what gradient, what fixedValue and mixing coefficient)
201     // Two reasonable choices:
202     // 1. specify above temperature on one side (preferentially the high side)
203     //    and above gradient on the other. So this will switch between pure
204     //    fixedvalue and pure fixedgradient
205     // 2. specify gradient and temperature such that the equations are the
206     //    same on both sides. This leads to the choice of
207     //    - refGradient = zero gradient
208     //    - refValue = neighbour value
209     //    - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
212     this->refValue() = nbrIntFld;
214     this->refGrad() = 0.0;
216     this->valueFraction() = nbrKDelta / (nbrKDelta + myKDelta());
218     mixedFvPatchScalarField::updateCoeffs();
220     if (debug)
221     {
222         scalar Q = gSum(K(*this)*patch().magSf()*snGrad());
224         Info<< patch().boundaryMesh().mesh().name() << ':'
225             << patch().name() << ':'
226             << this->dimensionedInternalField().name() << " <- "
227             << nbrMesh.name() << ':'
228             << nbrPatch.name() << ':'
229             << this->dimensionedInternalField().name() << " :"
230             << " heat[W]:" << Q
231             << " walltemperature "
232             << " min:" << gMin(*this)
233             << " max:" << gMax(*this)
234             << " avg:" << gAverage(*this)
235             << endl;
236     }
238     // Restore tag
239     UPstream::msgType() = oldTag;
243 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
245     Ostream& os
246 ) const
248     mixedFvPatchScalarField::write(os);
249     os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
250         << token::END_STATEMENT << nl;
251     temperatureCoupledBase::write(os);
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 makePatchTypeField
259     fvPatchScalarField,
260     turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 } // End namespace compressible
267 } // End namespace Foam
270 // ************************************************************************* //