thermalBaffleModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "thermalBaffleModel.H"
30 #include "fvMesh.H"
32 #include "wedgePolyPatch.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 namespace thermalBaffleModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(thermalBaffleModel, 0);
46 defineRunTimeSelectionTable(thermalBaffleModel, mesh);
47 defineRunTimeSelectionTable(thermalBaffleModel, dictionary);
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 {
55  return true;
56 }
57 
58 
60 {
62  return true;
63 }
64 
65 
66 void thermalBaffleModel::init()
67 {
68  if (active_)
69  {
70  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
71 
72  // Check if region mesh in 1-D
73  label nTotalEdges = 0;
74  const label patchi = intCoupledPatchIDs_[0];
75  nTotalEdges = 2*nLayers_*rbm[patchi].nInternalEdges();
76  nTotalEdges +=
77  nLayers_*(rbm[patchi].nEdges() - rbm[patchi].nInternalEdges());
78 
79  reduce(nTotalEdges, sumOp<label>());
80 
81  label nFaces = 0;
82  forAll(rbm, patchi)
83  {
84  if (
85  rbm[patchi].size()
86  &&
87  (
88  isA<wedgePolyPatch>(rbm[patchi])
89  || isA<emptyPolyPatch>(rbm[patchi])
90  )
91  )
92  {
93  nFaces += rbm[patchi].size();
94  }
95  }
96  reduce(nFaces, sumOp<label>());
97 
98  if (nTotalEdges == nFaces)
99  {
100  oneD_ = true;
101  Info << "\nThe thermal baffle is 1D \n" << endl;
102  }
103  else
104  {
105  Info << "\nThe thermal baffle is 3D \n" << endl;
106  }
107 
109  {
110  const label patchi = intCoupledPatchIDs_[i];
111  const polyPatch& pp = rbm[patchi];
112 
113  if
114  (
115  !isA<mappedVariableThicknessWallPolyPatch>(pp)
116  && oneD_
118  )
119  {
121  << "' not type '"
122  << mappedVariableThicknessWallPolyPatch::typeName
123  << "'. This is necessary for 1D solution "
124  << " and variable thickness"
125  << "\n for patch. " << pp.name()
126  << exit(FatalError);
127  }
128  else if (!isA<mappedWallPolyPatch>(pp))
129  {
131  << "' not type '"
132  << mappedWallPolyPatch::typeName
133  << "'. This is necessary for 3D solution"
134  << "\n for patch. " << pp.name()
135  << exit(FatalError);
136  }
137  }
138 
139  if (oneD_ && !constantThickness_)
140  {
141  const label patchi = intCoupledPatchIDs_[0];
142  const polyPatch& pp = rbm[patchi];
143  const mappedVariableThicknessWallPolyPatch& ppCoupled =
144  refCast
145  <
146  const mappedVariableThicknessWallPolyPatch
147  >(pp);
148 
149  thickness_ = ppCoupled.thickness();
150 
151  // Check that thickness has the right size
152  if (thickness_.size() != pp.size())
153  {
155  << " coupled patches in thermalBaffle are " << nl
156  << " different sizes from list thickness" << nl
157  << exit(FatalError);
158  }
159 
160  // Calculate thickness of the baffle on the first face only.
161  if (delta_.value() == 0.0)
162  {
163  forAll(ppCoupled, localFacei)
164  {
165  label facei = ppCoupled.start() + localFacei;
166 
167  label faceO =
168  boundaryFaceOppositeFace_[localFacei];
169 
170  delta_.value() = mag
171  (
172  regionMesh().faceCentres()[facei]
173  - regionMesh().faceCentres()[faceO]
174  );
175  break;
176  }
177  }
178  }
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
184 
185 thermalBaffleModel::thermalBaffleModel(const fvMesh& mesh)
186 :
187  regionModel1D(mesh, "thermalBaffle"),
188  thickness_(),
189  delta_("delta", dimLength, Zero),
190  oneD_(false),
191  constantThickness_(true)
192 {}
193 
194 
195 thermalBaffleModel::thermalBaffleModel
196 (
197  const word& modelType,
198  const fvMesh& mesh,
199  const dictionary& dict
200 
201 )
202 :
203  regionModel1D(mesh, "thermalBaffle", modelType, dict, true),
204  thickness_(),
205  delta_("delta", dimLength, Zero),
206  oneD_(false),
207  constantThickness_(dict.getOrDefault("constantThickness", true))
208 {
209  init();
210 }
211 
212 
213 thermalBaffleModel::thermalBaffleModel
214 (
215  const word& modelType,
216  const fvMesh& mesh
217 )
218 :
219  regionModel1D(mesh, "thermalBaffle", modelType),
220  thickness_(),
221  delta_("delta", dimLength, Zero),
222  oneD_(false),
223  constantThickness_(getOrDefault("constantThickness", true))
224 {
225  init();
226 }
227 
228 
229 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
230 
232 {}
233 
234 
235 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
236 
238 {}
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 } // End namespace thermalBaffleModels
244 } // End namespace regionModels
245 } // End namespace Foam
246 
247 // ************************************************************************* //
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::oneD_
bool oneD_
Is it one dimension.
Definition: thermalBaffleModel.H:89
Foam::regionModels::regionModel1D::read
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel1D.C:170
Foam::regionModels::regionModel1D::boundaryFaceOppositeFace_
labelList boundaryFaceOppositeFace_
Global boundary face IDs oppositte coupled patch.
Definition: regionModel1D.H:86
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::delta_
dimensionedScalar delta_
Baffle mesh thickness.
Definition: thermalBaffleModel.H:86
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::constantThickness_
bool constantThickness_
Is thickness constant.
Definition: thermalBaffleModel.H:92
wedgePolyPatch.H
thermalBaffleModel.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::sumOp
Definition: ops.H:213
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::thickness_
scalarField thickness_
Baffle physical thickness.
Definition: thermalBaffleModel.H:83
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
init
mesh init(true)
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::~thermalBaffleModel
virtual ~thermalBaffleModel()
Destructor.
Definition: thermalBaffleModel.C:231
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::regionModels::regionModel::active_
Switch active_
Active flag.
Definition: regionModel.H:93
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
mappedVariableThicknessWallPolyPatch.H
Foam::regionModels::regionModel1D::nLayers_
label nLayers_
Number of layers in the region.
Definition: regionModel1D.H:89
Foam::regionModels::thermalBaffleModels::defineTypeNameAndDebug
defineTypeNameAndDebug(noThermo, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::regionModel::intCoupledPatchIDs_
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:114
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: thermalBaffleModel.C:237
Foam::regionModels::regionModel1D
Base class for 1-D region models.
Definition: regionModel1D.H:54
Foam::regionModels::thermalBaffleModels::defineRunTimeSelectionTable
defineRunTimeSelectionTable(thermalBaffleModel, mesh)
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::read
virtual bool read()
Read control parameters from IO dictionary.
Definition: thermalBaffleModel.C:52