1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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"
30 #include "addToRunTimeSelectionTable.H"
31 #include "zeroGradientFvPatchFields.H"
32 #include "fvMatrices.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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()
65 Info<< "thermoBaffle2D::solveEnergy()" << endl;
68 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
70 tmp<volScalarField> tQ
77 regionMesh().time().timeName(),
84 dimensionedScalar("zero", dimEnergy/dimVolume/dimTime, 0.0)
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_)
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)
103 const labelList& cells = boundaryFaceCells_[localFaceI];
106 const label cellId = cells[i];
109 Qs_.boundaryField()[patchI][localFaceI]
110 /thickness_[localFaceI];
112 rhoCp[cellId] *= delta_.value()/thickness_[localFaceI];
114 K[cellId] *= delta_.value()/thickness_[localFaceI];
126 - fvm::laplacian(K, T_)
138 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 thermoBaffle2D::thermoBaffle2D
143 const word& modelType,
145 const dictionary& dict
148 thermoBaffleModel(modelType, mesh, dict),
149 nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
150 thermo_(basicSolidThermo::New(regionMesh(), dict)),
157 regionMesh().time().timeName(),
159 IOobject::READ_IF_PRESENT,
166 dimEnergy/dimArea/dimTime,
167 pTraits<scalar>::zero
175 regionMesh().time().timeName(),
177 IOobject::READ_IF_PRESENT,
184 dimEnergy/dimVolume/dimTime,
185 pTraits<scalar>::zero
194 thermoBaffle2D::thermoBaffle2D
196 const word& modelType,
200 thermoBaffleModel(modelType, mesh),
201 nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
202 thermo_(basicSolidThermo::New(regionMesh())),
209 regionMesh().time().timeName(),
211 IOobject::READ_IF_PRESENT,
218 dimEnergy/dimArea/dimTime,
219 pTraits<scalar>::zero
227 regionMesh().time().timeName(),
229 IOobject::READ_IF_PRESENT,
236 dimEnergy/dimVolume/dimTime,
237 pTraits<scalar>::zero
246 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
248 thermoBaffle2D::~thermoBaffle2D()
252 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
254 void thermoBaffle2D::init()
256 if (oneD_ && !constantThickness_)
258 label patchI = intCoupledPatchIDs_[0];
259 const label Qsb = Qs_.boundaryField()[patchI].size();
260 if (Qsb!= thickness_.size())
264 "thermoBaffle2D::thermoBaffle2D"
266 " const word& modelType,"
267 " const fvMesh& mesh,"
268 " const dictionary& dict"
270 ) << "the boundary field of Qs is "
271 << Qsb << " and " << nl
272 << "the field 'thickness' is " << thickness_.size() << nl
279 void thermoBaffle2D::preEvolveRegion()
283 void thermoBaffle2D::evolveRegion()
285 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
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
316 const volScalarField& thermoBaffle2D::T() const
322 const basicSolidThermo& thermoBaffle2D::thermo() const
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)
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 <<
342 mag(regionMesh().Sf().boundaryField()[patchI])
344 * thermo_->K().boundaryField()[patchI]
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 } // end namespace thermoBaffleModels
353 } // end namespace regionModels
354 } // end namespace Foam
356 // ************************************************************************* //