BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / regionModels / thermoBaffleModels / thermoBaffle2D / thermoBaffle2D.C
blobd46fb80a62cc6b8d79379fb1858d205c1e33fa41
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 the
13     Free Software Foundation; either version 3 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 "thermoBaffle2D.H"
29 #include "fvm.H"
30 #include "addToRunTimeSelectionTable.H"
31 #include "zeroGradientFvPatchFields.H"
32 #include "fvMatrices.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
39 namespace regionModels
41 namespace thermoBaffleModels
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 defineTypeNameAndDebug(thermoBaffle2D, 0);
48 addToRunTimeSelectionTable(thermoBaffleModel, thermoBaffle2D, mesh);
49 addToRunTimeSelectionTable(thermoBaffleModel, thermoBaffle2D, dictionary);
51 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
54 bool thermoBaffle2D::read()
56     this->solution().lookup("nNonOrthCorr") >> nNonOrthCorr_;
57     return regionModel1D::read();
61 void thermoBaffle2D::solveEnergy()
63     if (debug)
64     {
65         Info<< "thermoBaffle2D::solveEnergy()" << endl;
66     }
68     const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
70     tmp<volScalarField> tQ
71     (
72         new volScalarField
73         (
74             IOobject
75             (
76                 "tQ",
77                 regionMesh().time().timeName(),
78                 regionMesh(),
79                 IOobject::NO_READ,
80                 IOobject::NO_WRITE,
81                 false
82             ),
83             regionMesh(),
84             dimensionedScalar("zero", dimEnergy/dimVolume/dimTime, 0.0)
85         )
86     );
88     volScalarField& Q = tQ();
90     volScalarField rhoCp("rhoCp", thermo_->rho()*thermo_->Cp()());
91     volScalarField K("K", thermo_->K());
94     //If region is one-dimension variable thickness
95     if (oneD_ && !constantThickness_)
96     {
97         // Scale K and rhoCp and fill Q in the internal baffle region.
98         const label patchI = intCoupledPatchIDs_[0];
99         const polyPatch& ppCoupled = rbm[patchI];
101         forAll(ppCoupled, localFaceI)
102         {
103             const labelList& cells = boundaryFaceCells_[localFaceI];
104             forAll (cells, i)
105             {
106                 const label cellId = cells[i];
108                 Q[cellId] =
109                     Qs_.boundaryField()[patchI][localFaceI]
110                     /thickness_[localFaceI];
112                 rhoCp[cellId] *= delta_.value()/thickness_[localFaceI];
114                 K[cellId] *= delta_.value()/thickness_[localFaceI];
115             }
116         }
117     }
118     else
119     {
120         Q = Q_;
121     }
123     fvScalarMatrix TEqn
124     (
125         fvm::ddt(rhoCp, T_)
126       - fvm::laplacian(K, T_)
127      ==
128         Q
129     );
131     TEqn.relax();
132     TEqn.solve();
134     thermo_->correct();
138 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
141 thermoBaffle2D::thermoBaffle2D
143     const word& modelType,
144     const fvMesh& mesh,
145     const dictionary& dict
148     thermoBaffleModel(modelType, mesh, dict),
149     nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
150     thermo_(basicSolidThermo::New(regionMesh(), dict)),
151     T_(thermo_->T()),
152     Qs_
153     (
154         IOobject
155         (
156             "Qs",
157             regionMesh().time().timeName(),
158             regionMesh(),
159             IOobject::READ_IF_PRESENT,
160             IOobject::AUTO_WRITE
161         ),
162         regionMesh(),
163         dimensionedScalar
164         (
165             "zero",
166             dimEnergy/dimArea/dimTime,
167             pTraits<scalar>::zero
168         )
169     ),
170     Q_
171     (
172         IOobject
173         (
174             "Q",
175             regionMesh().time().timeName(),
176             regionMesh(),
177             IOobject::READ_IF_PRESENT,
178             IOobject::AUTO_WRITE
179         ),
180         regionMesh(),
181         dimensionedScalar
182         (
183             "zero",
184             dimEnergy/dimVolume/dimTime,
185             pTraits<scalar>::zero
186         )
187     )
189     init();
190     thermo_->correct();
194 thermoBaffle2D::thermoBaffle2D
196     const word& modelType,
197     const fvMesh& mesh
200     thermoBaffleModel(modelType, mesh),
201     nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
202     thermo_(basicSolidThermo::New(regionMesh())),
203     T_(thermo_->T()),
204     Qs_
205     (
206         IOobject
207         (
208             "Qs",
209             regionMesh().time().timeName(),
210             regionMesh(),
211             IOobject::READ_IF_PRESENT,
212             IOobject::AUTO_WRITE
213         ),
214         regionMesh(),
215         dimensionedScalar
216         (
217             "zero",
218             dimEnergy/dimArea/dimTime,
219             pTraits<scalar>::zero
220         )
221     ),
222     Q_
223     (
224         IOobject
225         (
226             "Q",
227             regionMesh().time().timeName(),
228             regionMesh(),
229             IOobject::READ_IF_PRESENT,
230             IOobject::AUTO_WRITE
231         ),
232         regionMesh(),
233         dimensionedScalar
234         (
235             "zero",
236             dimEnergy/dimVolume/dimTime,
237             pTraits<scalar>::zero
238         )
239     )
241     init();
242     thermo_->correct();
246 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
248 thermoBaffle2D::~thermoBaffle2D()
252 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
254 void thermoBaffle2D::init()
256     if (oneD_ && !constantThickness_)
257     {
258         label patchI = intCoupledPatchIDs_[0];
259         const label Qsb = Qs_.boundaryField()[patchI].size();
260         if (Qsb!= thickness_.size())
261         {
262             FatalErrorIn
263             (
264                 "thermoBaffle2D::thermoBaffle2D"
265                 "("
266                 "   const word& modelType,"
267                 "   const fvMesh& mesh,"
268                 "   const dictionary& dict"
269                 ")"
270             )   << "the boundary field of Qs is "
271                 << Qsb << " and " << nl
272                 << "the field 'thickness' is " << thickness_.size() << nl
273                 << exit(FatalError);
274         }
275     }
279 void thermoBaffle2D::preEvolveRegion()
283 void thermoBaffle2D::evolveRegion()
285     for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
286     {
287         solveEnergy();
288     }
292 const tmp<volScalarField> thermoBaffle2D::Cp() const
294     return thermo_->Cp();
298 const volScalarField& thermoBaffle2D::kappa() const
300     return thermo_->kappa();
304 const volScalarField& thermoBaffle2D::rho() const
306     return thermo_->rho();
310 const volScalarField& thermoBaffle2D::K() const
312     return thermo_->K();
316 const volScalarField& thermoBaffle2D::T() const
318     return T_;
322 const basicSolidThermo& thermoBaffle2D::thermo() const
324     return thermo_;
328 void thermoBaffle2D::info() const
330     Info<< indent << "min/max(T) = " << min(T_).value() << ", "
331         << max(T_).value() << nl;
333     const labelList& coupledPatches = intCoupledPatchIDs();
334     forAll (coupledPatches, i)
335     {
336         const label patchI = coupledPatches[i];
337         const fvPatchScalarField& pT = T_.boundaryField()[patchI];
338         const word patchName = regionMesh().boundary()[patchI].name();
339         Info << indent << "Q : " << patchName << indent <<
340             gSum
341             (
342                 mag(regionMesh().Sf().boundaryField()[patchI])
343               * pT.snGrad()
344               * thermo_->K().boundaryField()[patchI]
345             ) << endl;
346     }
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 } // end namespace thermoBaffleModels
353 } // end namespace regionModels
354 } // end namespace Foam
356 // ************************************************************************* //