1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011-2011 OpenCFD Ltd.
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);
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 bool thermoBaffle2D::read()
55 this->solution().lookup("nNonOrthCorr") >> nNonOrthCorr_;
56 return regionModel1D::read();
60 void thermoBaffle2D::solveEnergy()
64 Info<< "thermoBaffle2D::solveEnergy()" << endl;
67 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
69 tmp<volScalarField> tQ
76 regionMesh().time().timeName(),
83 dimensionedScalar("zero", dimEnergy/dimVolume/dimTime, 0.0)
87 volScalarField& Q = tQ();
89 volScalarField rhoCp("rhoCp", thermo_->rho()*thermo_->Cp()());
90 volScalarField K("K", thermo_->K());
93 //If region is one-dimension variable thickness can be used.
96 // Scale K and rhoCp and fill Q in the internal baffle region.
97 const label patchI = intCoupledPatchIDs_[0];
98 const polyPatch& ppCoupled = rbm[patchI];
100 forAll(ppCoupled, localFaceI)
102 const labelList& cells = boundaryFaceCells_[localFaceI];
105 const label cellId = cells[i];
108 Qs_.boundaryField()[patchI][localFaceI]
109 /thickness_[localFaceI];
111 rhoCp[cellId] *= delta_.value()/thickness_[localFaceI];
113 K[cellId] *= delta_.value()/thickness_[localFaceI];
125 - fvm::laplacian(K, T_)
137 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 thermoBaffle2D::thermoBaffle2D
141 const word& modelType,
145 thermoBaffleModel(modelType, mesh),
146 nNonOrthCorr_(readLabel(solution().lookup("nNonOrthCorr"))),
147 thermo_(basicSolidThermo::New(regionMesh())),
154 regionMesh().time().timeName(),
156 IOobject::READ_IF_PRESENT,
163 dimEnergy/dimArea/dimTime,
164 pTraits<scalar>::zero
172 regionMesh().time().timeName(),
174 IOobject::READ_IF_PRESENT,
181 dimEnergy/dimVolume/dimTime,
182 pTraits<scalar>::zero
188 label patchI = intCoupledPatchIDs_[0];
189 const label Qsb = Qs_.boundaryField()[patchI].size();
190 if (Qsb!= thickness_.size())
194 "thermoBaffle2D::thermoBaffle2D"
196 " const word& modelType,"
197 " const fvMesh& mesh"
199 ) << "the boundary field of Qs is "
200 << Qsb << " and " << nl
201 << "the field 'thickness' is " << thickness_.size() << nl
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212 thermoBaffle2D::~thermoBaffle2D()
216 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
218 void thermoBaffle2D::preEvolveRegion()
222 void thermoBaffle2D::evolveRegion()
224 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
231 const tmp<volScalarField> thermoBaffle2D::Cp() const
233 return thermo_->Cp();
237 const volScalarField& thermoBaffle2D::kappa() const
239 return thermo_->kappa();
243 const volScalarField& thermoBaffle2D::rho() const
245 return thermo_->rho();
249 const volScalarField& thermoBaffle2D::K() const
255 const volScalarField& thermoBaffle2D::T() const
261 void thermoBaffle2D::info() const
263 Info<< indent << "min/max(T) = " << min(T_).value() << ", "
264 << max(T_).value() << nl;
266 const labelList& coupledPatches = intCoupledPatchIDs();
267 forAll (coupledPatches, i)
269 const label patchI = coupledPatches[i];
270 const fvPatchScalarField& pT = T_.boundaryField()[patchI];
271 const word patchName = regionMesh().boundary()[patchI].name();
272 Info << indent << "Q : " << patchName << indent <<
275 mag(regionMesh().Sf().boundaryField()[patchI])
277 * thermo_->K().boundaryField()[patchI]
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 } // end namespace thermoBaffleModels
286 } // end namespace regionModels
287 } // end namespace Foam
289 // ************************************************************************* //