regionModel1D.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "regionModel1D.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace regionModels
33 {
34  defineTypeNameAndDebug(regionModel1D, 0);
35 }
36 }
37 
38 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39 
40 void Foam::regionModels::regionModel1D::constructMeshObjects()
41 {
42  nMagSfPtr_.reset
43  (
45  (
46  IOobject
47  (
48  "nMagSf",
49  time().timeName(),
50  regionMesh(),
53  ),
54  regionMesh(),
56  )
57  );
58 }
59 
60 
61 void Foam::regionModels::regionModel1D::initialise()
62 {
64 
65  // Calculate boundaryFaceFaces and boundaryFaceCells
66 
67  DynamicList<label> faceIDs;
68  DynamicList<label> cellIDs;
69 
70  label localPyrolysisFacei = 0;
71 
72  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
73 
74  forAll(intCoupledPatchIDs_, i)
75  {
76  const label patchi = intCoupledPatchIDs_[i];
77  const polyPatch& ppCoupled = rbm[patchi];
78  forAll(ppCoupled, localFacei)
79  {
80  label facei = ppCoupled.start() + localFacei;
81  label celli = -1;
82  label nFaces = 0;
83  label nCells = 0;
84  do
85  {
86  label ownCelli = regionMesh().faceOwner()[facei];
87  if (ownCelli != celli)
88  {
89  celli = ownCelli;
90  }
91  else
92  {
93  celli = regionMesh().faceNeighbour()[facei];
94  }
95  nCells++;
96  cellIDs.append(celli);
97  const cell& cFaces = regionMesh().cells()[celli];
98  facei = cFaces.opposingFaceLabel(facei, regionMesh().faces());
99  faceIDs.append(facei);
100  nFaces++;
101  } while (regionMesh().isInternalFace(facei));
102 
103  boundaryFaceOppositeFace_[localPyrolysisFacei] = facei;
104  faceIDs.remove(); // remove boundary face.
105  nFaces--;
106 
107  boundaryFaceFaces_[localPyrolysisFacei].transfer(faceIDs);
108  boundaryFaceCells_[localPyrolysisFacei].transfer(cellIDs);
109 
110  localPyrolysisFacei++;
111  nLayers_ = nCells;
112  }
113  }
114 
115  boundaryFaceOppositeFace_.setSize(localPyrolysisFacei);
116  boundaryFaceFaces_.setSize(localPyrolysisFacei);
117  boundaryFaceCells_.setSize(localPyrolysisFacei);
118 
119  surfaceScalarField& nMagSf = nMagSfPtr_();
120 
121  surfaceScalarField::Boundary nMagSfBf =
122  nMagSf.boundaryFieldRef();
123 
124  localPyrolysisFacei = 0;
125 
126  forAll(intCoupledPatchIDs_, i)
127  {
128  const label patchi = intCoupledPatchIDs_[i];
129  const polyPatch& ppCoupled = rbm[patchi];
130  const vectorField& pNormals = ppCoupled.faceNormals();
131  nMagSfBf[patchi] = regionMesh().Sf().boundaryField()[patchi] & pNormals;
132  forAll(pNormals, localFacei)
133  {
134  const vector& n = pNormals[localFacei];
135  const labelList& faces = boundaryFaceFaces_[localPyrolysisFacei++];
136  forAll(faces, facei)
137  {
138  const label faceID = faces[facei];
139  nMagSf[faceID] = regionMesh().Sf()[faceID] & n;
140  }
141  }
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
147 
149 {
150  if (regionModel::read())
151  {
152  moveMesh_.readIfPresent("moveMesh", coeffs_);
153 
154  return true;
155  }
156  else
157  {
158  return false;
159  }
160 }
161 
162 
164 {
165  if (regionModel::read(dict))
166  {
167  moveMesh_.readIfPresent("moveMesh", coeffs_);
168 
169  return true;
170  }
171  else
172  {
173  return false;
174  }
175 }
176 
177 
179 (
180  const scalarList& deltaV,
181  const scalar minDelta
182 )
183 {
184  tmp<labelField> tcellMoveMap(new labelField(regionMesh().nCells(), 0));
185  labelField& cellMoveMap = tcellMoveMap.ref();
186 
187  if (!moveMesh_)
188  {
189  return cellMoveMap;
190  }
191 
192  pointField oldPoints = regionMesh().points();
193  pointField newPoints = oldPoints;
194 
195  const polyBoundaryMesh& bm = regionMesh().boundaryMesh();
196 
197  label totalFaceId = 0;
198  forAll(intCoupledPatchIDs_, localPatchi)
199  {
200  label patchi = intCoupledPatchIDs_[localPatchi];
201  const polyPatch pp = bm[patchi];
202  const vectorField& cf = regionMesh().Cf().boundaryField()[patchi];
203 
204  forAll(pp, patchFacei)
205  {
206  const labelList& faces = boundaryFaceFaces_[totalFaceId];
207  const labelList& cells = boundaryFaceCells_[totalFaceId];
208 
209  const vector n = pp.faceNormals()[patchFacei];
210  const vector sf = pp.faceAreas()[patchFacei];
211 
212  List<point> oldCf(faces.size() + 1);
213  oldCf[0] = cf[patchFacei];
214  forAll(faces, i)
215  {
216  oldCf[i + 1] = regionMesh().faceCentres()[faces[i]];
217  }
218 
219  vector newDelta = Zero;
220  point nbrCf = oldCf[0];
221 
222  forAll(faces, i)
223  {
224  const label facei = faces[i];
225  const label celli = cells[i];
226 
227  const face f = regionMesh().faces()[facei];
228 
229  newDelta += (deltaV[celli]/mag(sf))*n;
230 
231  vector localDelta = Zero;
232  forAll(f, pti)
233  {
234  const label pointi = f[pti];
235 
236  if
237  (
238  mag((nbrCf - (oldPoints[pointi] + newDelta)) & n)
239  > minDelta
240  )
241  {
242  newPoints[pointi] = oldPoints[pointi] + newDelta;
243  localDelta = newDelta;
244  cellMoveMap[celli] = 1;
245  }
246  }
247  nbrCf = oldCf[i + 1] + localDelta;
248  }
249  // Modify boundary
250  const label bFacei = boundaryFaceOppositeFace_[totalFaceId];
251  const face f = regionMesh().faces()[bFacei];
252  const label celli = cells[cells.size() - 1];
253  newDelta += (deltaV[celli]/mag(sf))*n;
254  forAll(f, pti)
255  {
256  const label pointi = f[pti];
257  if
258  (
259  mag((nbrCf - (oldPoints[pointi] + newDelta)) & n)
260  > minDelta
261  )
262  {
263  newPoints[pointi] = oldPoints[pointi] + newDelta;
264  cellMoveMap[celli] = 1;
265  }
266  }
267  totalFaceId ++;
268  }
269  }
270  // Move points
271  regionMesh().movePoints(newPoints);
272 
273  return tcellMoveMap;
274 }
275 
276 
277 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
278 
279 Foam::regionModels::regionModel1D::regionModel1D
280 (
281  const fvMesh& mesh,
282  const word& regionType
283 )
284 :
285  regionModel(mesh, regionType),
286  boundaryFaceFaces_(),
287  boundaryFaceCells_(),
288  boundaryFaceOppositeFace_(),
289  nLayers_(0),
290  nMagSfPtr_(nullptr),
291  moveMesh_(false)
292 {}
293 
294 
295 Foam::regionModels::regionModel1D::regionModel1D
296 (
297  const fvMesh& mesh,
298  const word& regionType,
299  const word& modelName,
300  bool readFields
301 )
302 :
303  regionModel(mesh, regionType, modelName, false),
304  boundaryFaceFaces_(regionMesh().nCells()),
305  boundaryFaceCells_(regionMesh().nCells()),
306  boundaryFaceOppositeFace_(regionMesh().nCells()),
307  nLayers_(0),
308  nMagSfPtr_(nullptr),
309  moveMesh_(true)
310 {
311  if (active_)
312  {
313  constructMeshObjects();
314  initialise();
315 
316  if (readFields)
317  {
318  read();
319  }
320  }
321 }
322 
323 
324 Foam::regionModels::regionModel1D::regionModel1D
325 (
326  const fvMesh& mesh,
327  const word& regionType,
328  const word& modelName,
329  const dictionary& dict,
330  bool readFields
331 )
332 :
333  regionModel(mesh, regionType, modelName, dict, readFields),
334  boundaryFaceFaces_(regionMesh().nCells()),
335  boundaryFaceCells_(regionMesh().nCells()),
336  boundaryFaceOppositeFace_(regionMesh().nCells()),
337  nLayers_(0),
338  nMagSfPtr_(nullptr),
339  moveMesh_(false)
340 {
341  if (active_)
342  {
343  constructMeshObjects();
344  initialise();
345 
346  if (readFields)
347  {
348  read(dict);
349  }
350  }
351 }
352 
353 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
354 
356 {}
357 
358 
359 // ************************************************************************* //
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< labelField > moveMesh(const scalarList &deltaV, const scalar minDelta=0.0)
Move mesh points according to change in cell volumes.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
defineTypeNameAndDebug(regionModel, 0)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const cellShapeList & cells
virtual ~regionModel1D()
Destructor.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:141
#define DebugInFunction
Report an information message using Foam::Info.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
const Field< PointType > & faceNormals() const
Return face normals for patch.
Foam::polyBoundaryMesh.
const vectorField::subField faceAreas() const
Return face areas.
Definition: polyPatch.C:290
labelList f(nPoints)
volScalarField sf(fieldObject, mesh)
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read control parameters from dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Base class for region models.
Definition: regionModel.H:55
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.