BUG: Have sensible value for refValue.
[OpenFOAM-1.6.x.git] / applications / solvers / heatTransfer / chtMultiRegionFoam / derivedFvPatchFields / solidWallMixedTemperatureCoupled / solidWallMixedTemperatureCoupledFvPatchScalarField.C
blob5d5c8d999dcb39b1a226ef5fcd0f280a89d53192
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 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
19     for more details.
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 "solidWallMixedTemperatureCoupledFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "directMappedPatchBase.H"
32 #include "regionProperties.H"
34 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
36 bool Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::interfaceOwner
38     const polyMesh& nbrRegion
39 ) const
41     const fvMesh& myRegion = patch().boundaryMesh().mesh();
43     const regionProperties& props =
44         myRegion.objectRegistry::parent().lookupObject<regionProperties>
45         (
46             "regionProperties"
47         );
49     label myIndex = findIndex(props.fluidRegionNames(), myRegion.name());
50     if (myIndex == -1)
51     {
52         label i = findIndex(props.solidRegionNames(), myRegion.name());
54         if (i == -1)
55         {
56             FatalErrorIn
57             (
58                 "solidWallMixedTemperatureCoupledFvPatchScalarField"
59                 "::interfaceOwner(const polyMesh&) const"
60             )   << "Cannot find region " << myRegion.name()
61                 << " neither in fluids " << props.fluidRegionNames()
62                 << " nor in solids " << props.solidRegionNames()
63                 << exit(FatalError);
64         }
65         myIndex = props.fluidRegionNames().size() + i;
66     }
67     label nbrIndex = findIndex(props.fluidRegionNames(), nbrRegion.name());
68     if (nbrIndex == -1)
69     {
70         label i = findIndex(props.solidRegionNames(), nbrRegion.name());
72         if (i == -1)
73         {
74             FatalErrorIn("coupleManager::interfaceOwner(const polyMesh&) const")
75                 << "Cannot find region " << nbrRegion.name()
76                 << " neither in fluids " << props.fluidRegionNames()
77                 << " nor in solids " << props.solidRegionNames()
78                 << exit(FatalError);
79         }
80         nbrIndex = props.fluidRegionNames().size() + i;
81     }
83     return myIndex < nbrIndex;
87 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
89 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
90 solidWallMixedTemperatureCoupledFvPatchScalarField
92     const fvPatch& p,
93     const DimensionedField<scalar, volMesh>& iF
96     mixedFvPatchScalarField(p, iF),
97     neighbourFieldName_("undefined-neighbourFieldName"),
98     KName_("undefined-K")
100     this->refValue() = 0.0;
101     this->refGrad() = 0.0;
102     this->valueFraction() = 1.0;
103     this->fixesValue_ = true;
107 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
108 solidWallMixedTemperatureCoupledFvPatchScalarField
110     const solidWallMixedTemperatureCoupledFvPatchScalarField& ptf,
111     const fvPatch& p,
112     const DimensionedField<scalar, volMesh>& iF,
113     const fvPatchFieldMapper& mapper
116     mixedFvPatchScalarField(ptf, p, iF, mapper),
117     neighbourFieldName_(ptf.neighbourFieldName_),
118     KName_(ptf.KName_),
119     fixesValue_(ptf.fixesValue_)
123 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
124 solidWallMixedTemperatureCoupledFvPatchScalarField
126     const fvPatch& p,
127     const DimensionedField<scalar, volMesh>& iF,
128     const dictionary& dict
131     mixedFvPatchScalarField(p, iF),
132     neighbourFieldName_(dict.lookup("neighbourFieldName")),
133     KName_(dict.lookup("K"))
135     if (!isA<directMappedPatchBase>(this->patch().patch()))
136     {
137         FatalErrorIn
138         (
139             "solidWallMixedTemperatureCoupledFvPatchScalarField::"
140             "solidWallMixedTemperatureCoupledFvPatchScalarField\n"
141             "(\n"
142             "    const fvPatch& p,\n"
143             "    const DimensionedField<scalar, volMesh>& iF,\n"
144             "    const dictionary& dict\n"
145             ")\n"
146         )   << "\n    patch type '" << p.type()
147             << "' not type '" << directMappedPatchBase::typeName << "'"
148             << "\n    for patch " << p.name()
149             << " of field " << dimensionedInternalField().name()
150             << " in file " << dimensionedInternalField().objectPath()
151             << exit(FatalError);
152     }
154     fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
156     if (dict.found("refValue"))
157     {
158         // Full restart
159         refValue() = scalarField("refValue", dict, p.size());
160         refGrad() = scalarField("refGradient", dict, p.size());
161         valueFraction() = scalarField("valueFraction", dict, p.size());
162         fixesValue_ = readBool(dict.lookup("fixesValue"));
163     }
164     else
165     {
166         // Start from user entered data. Assume fixedValue.
167         refValue() = *this;
168         refGrad() = 0.0;
169         valueFraction() = 1.0;
170         fixesValue_ = true;
171     }
175 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
176 solidWallMixedTemperatureCoupledFvPatchScalarField
178     const solidWallMixedTemperatureCoupledFvPatchScalarField& wtcsf,
179     const DimensionedField<scalar, volMesh>& iF
182     mixedFvPatchScalarField(wtcsf, iF),
183     neighbourFieldName_(wtcsf.neighbourFieldName_),
184     KName_(wtcsf.KName_),
185     fixesValue_(wtcsf.fixesValue_)
189 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
191 const Foam::fvPatchScalarField&
192 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::K() const
194     return this->patch().lookupPatchField<volScalarField, scalar>(KName_);
198 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::updateCoeffs()
200     if (updated())
201     {
202         return;
203     }
205     // Get the coupling information from the directMappedPatchBase
206     const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
207     (
208         patch().patch()
209     );
210     const polyMesh& nbrMesh = mpp.sampleMesh();
211     // Force recalculation of mapping and schedule
212     const mapDistribute& distMap = mpp.map();
213     (void)distMap.schedule();
215     tmp<scalarField> intFld = patchInternalField();
217     if (interfaceOwner(nbrMesh))
218     {
219         // Note: other side information could be cached - it only needs
220         // to be updated the first time round the iteration (i.e. when
221         // switching regions) but unfortunately we don't have this information.
223         const fvPatch& nbrPatch = refCast<const fvMesh>
224         (
225             nbrMesh
226         ).boundary()[mpp.samplePolyPatch().index()];
229         // Calculate the temperature by harmonic averaging
230         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
232         const solidWallMixedTemperatureCoupledFvPatchScalarField& nbrField =
233         refCast<const solidWallMixedTemperatureCoupledFvPatchScalarField>
234         (
235             nbrPatch.lookupPatchField<volScalarField, scalar>
236             (
237                 neighbourFieldName_
238             )
239         );
241         // Swap to obtain full local values of neighbour internal field
242         scalarField nbrIntFld = nbrField.patchInternalField();
243         mapDistribute::distribute
244         (
245             Pstream::defaultCommsType,
246             distMap.schedule(),
247             distMap.constructSize(),
248             distMap.subMap(),           // what to send
249             distMap.constructMap(),     // what to receive
250             nbrIntFld
251         );
253         // Swap to obtain full local values of neighbour K*delta
254         scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
255         mapDistribute::distribute
256         (
257             Pstream::defaultCommsType,
258             distMap.schedule(),
259             distMap.constructSize(),
260             distMap.subMap(),           // what to send
261             distMap.constructMap(),     // what to receive
262             nbrKDelta
263         );
266         tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();
268         // Calculate common wall temperature. Reuse *this to store common value.
269         scalarField Twall
270         (
271             (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
272           / (myKDelta() + nbrKDelta)
273         );
274         // Assign to me
275         fvPatchScalarField::operator=(Twall);
276         // Distribute back and assign to neighbour
277         mapDistribute::distribute
278         (
279             Pstream::defaultCommsType,
280             distMap.schedule(),
281             nbrField.size(),
282             distMap.constructMap(),     // reverse : what to send
283             distMap.subMap(),
284             Twall
285         );
286         const_cast<solidWallMixedTemperatureCoupledFvPatchScalarField&>
287         (
288             nbrField
289         ).fvPatchScalarField::operator=(Twall);
290     }
293     // Switch between fixed value (of harmonic avg) or gradient
294     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
296     label nFixed = 0;
298     // Like snGrad but bypass switching on refValue/refGrad.
299     tmp<scalarField> normalGradient = (*this-intFld())*patch().deltaCoeffs();
301     if (debug)
302     {
303         scalar Q = gSum(K()*patch().magSf()*normalGradient());
305         Info<< "solidWallMixedTemperatureCoupledFvPatchScalarField::"
306             << "updateCoeffs() :"
307             << " patch:" << patch().name()
308             << " heatFlux:" << Q
309             << " walltemperature "
310             << " min:" << gMin(*this)
311             << " max:" << gMax(*this)
312             << " avg:" << gAverage(*this)
313             << endl;
314     }
316     forAll(*this, i)
317     {
318         // if outgoing flux use fixed value.
319         if (normalGradient()[i] < 0.0)
320         {
321             this->refValue()[i] = operator[](i);
322             this->refGrad()[i] = 0.0;   // not used by me
323             this->valueFraction()[i] = 1.0;
324             nFixed++;
325         }
326         else
327         {
328             // Fixed gradient. Make sure to have valid refValue (even though
329             // I am not using it - other boundary conditions might)
330             this->refValue()[i] = operator[](i);
331             this->refGrad()[i] = normalGradient()[i];
332             this->valueFraction()[i] = 0.0;
333         }
334     }
336     reduce(nFixed, sumOp<label>());
338     fixesValue_ = (nFixed > 0);
340     if (debug)
341     {
342         label nTotSize = returnReduce(this->size(), sumOp<label>());
344         Info<< "solidWallMixedTemperatureCoupledFvPatchScalarField::"
345             << "updateCoeffs() :"
346             << " patch:" << patch().name()
347             << " out of:" << nTotSize
348             << " fixedBC:" << nFixed
349             << " gradient:" << nTotSize-nFixed << endl;
350     }
352     mixedFvPatchScalarField::updateCoeffs();
356 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::write
358     Ostream& os
359 ) const
361     mixedFvPatchScalarField::write(os);
362     os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
363         << token::END_STATEMENT << nl;
364     os.writeKeyword("K") << KName_ << token::END_STATEMENT << nl;
365     os.writeKeyword("fixesValue") << fixesValue_ << token::END_STATEMENT << nl;
369 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 namespace Foam
374 makePatchTypeField
376     fvPatchScalarField,
377     solidWallMixedTemperatureCoupledFvPatchScalarField
380 } // End namespace Foam
382 // ************************************************************************* //