regionModel1D.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 | Copyright (C) 2011-2016 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 
43  nMagSfPtr_.reset
44  (
46  (
47  IOobject
48  (
49  "nMagSf",
50  time().timeName(),
51  regionMesh(),
54  ),
55  regionMesh(),
56  dimensionedScalar("zero", dimArea, 0.0)
57  )
58  );
59 }
60 
61 
62 void Foam::regionModels::regionModel1D::initialise()
63 {
64  if (debug)
65  {
66  Pout<< "regionModel1D::initialise()" << endl;
67  }
68 
69  // Calculate boundaryFaceFaces and boundaryFaceCells
70 
71  DynamicList<label> faceIDs;
72  DynamicList<label> cellIDs;
73 
74  label localPyrolysisFacei = 0;
75 
76  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
77 
78  forAll(intCoupledPatchIDs_, i)
79  {
80  const label patchi = intCoupledPatchIDs_[i];
81  const polyPatch& ppCoupled = rbm[patchi];
82  forAll(ppCoupled, localFacei)
83  {
84  label facei = ppCoupled.start() + localFacei;
85  label celli = -1;
86  label nFaces = 0;
87  label nCells = 0;
88  do
89  {
90  label ownCelli = regionMesh().faceOwner()[facei];
91  if (ownCelli != celli)
92  {
93  celli = ownCelli;
94  }
95  else
96  {
97  celli = regionMesh().faceNeighbour()[facei];
98  }
99  nCells++;
100  cellIDs.append(celli);
101  const cell& cFaces = regionMesh().cells()[celli];
102  facei = cFaces.opposingFaceLabel(facei, regionMesh().faces());
103  faceIDs.append(facei);
104  nFaces++;
105  } while (regionMesh().isInternalFace(facei));
106 
107  boundaryFaceOppositeFace_[localPyrolysisFacei] = facei;
108  faceIDs.remove(); //remove boundary face.
109  nFaces--;
110 
111  boundaryFaceFaces_[localPyrolysisFacei].transfer(faceIDs);
112  boundaryFaceCells_[localPyrolysisFacei].transfer(cellIDs);
113 
114  localPyrolysisFacei++;
115  nLayers_ = nCells;
116  }
117  }
118 
119  boundaryFaceOppositeFace_.setSize(localPyrolysisFacei);
120  boundaryFaceFaces_.setSize(localPyrolysisFacei);
121  boundaryFaceCells_.setSize(localPyrolysisFacei);
122 
123  surfaceScalarField& nMagSf = nMagSfPtr_();
124 
125  surfaceScalarField::Boundary nMagSfBf =
126  nMagSf.boundaryFieldRef();
127 
128  localPyrolysisFacei = 0;
129 
130  forAll(intCoupledPatchIDs_, i)
131  {
132  const label patchi = intCoupledPatchIDs_[i];
133  const polyPatch& ppCoupled = rbm[patchi];
134  const vectorField& pNormals = ppCoupled.faceNormals();
135  nMagSfBf[patchi] = regionMesh().Sf().boundaryField()[patchi] & pNormals;
136  forAll(pNormals, localFacei)
137  {
138  const vector& n = pNormals[localFacei];
139  const labelList& faces = boundaryFaceFaces_[localPyrolysisFacei++];
140  forAll(faces, facei)
141  {
142  const label faceID = faces[facei];
143  nMagSf[faceID] = regionMesh().Sf()[faceID] & n;
144  }
145  }
146  }
147 }
148 
149 
150 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
151 
153 {
154  if (regionModel::read())
155  {
156  moveMesh_.readIfPresent("moveMesh", coeffs_);
157 
158  return true;
159  }
160  else
161  {
162  return false;
163  }
164 }
165 
166 
168 {
169  if (regionModel::read(dict))
170  {
171  moveMesh_.readIfPresent("moveMesh", coeffs_);
172 
173  return true;
174  }
175  else
176  {
177  return false;
178  }
179 }
180 
181 
183 (
184  const scalarList& deltaV,
185  const scalar minDelta
186 )
187 {
188  tmp<labelField> tcellMoveMap(new labelField(regionMesh().nCells(), 0));
189  labelField& cellMoveMap = tcellMoveMap.ref();
190 
191  if (!moveMesh_)
192  {
193  return cellMoveMap;
194  }
195 
196  pointField oldPoints = regionMesh().points();
197  pointField newPoints = oldPoints;
198 
199  const polyBoundaryMesh& bm = regionMesh().boundaryMesh();
200 
201  label totalFaceId = 0;
202  forAll(intCoupledPatchIDs_, localPatchi)
203  {
204  label patchi = intCoupledPatchIDs_[localPatchi];
205  const polyPatch pp = bm[patchi];
206  const vectorField& cf = regionMesh().Cf().boundaryField()[patchi];
207 
208  forAll(pp, patchFacei)
209  {
210  const labelList& faces = boundaryFaceFaces_[totalFaceId];
211  const labelList& cells = boundaryFaceCells_[totalFaceId];
212 
213  const vector n = pp.faceNormals()[patchFacei];
214  const vector sf = pp.faceAreas()[patchFacei];
215 
216  List<point> oldCf(faces.size() + 1);
217  oldCf[0] = cf[patchFacei];
218  forAll(faces, i)
219  {
220  oldCf[i + 1] = regionMesh().faceCentres()[faces[i]];
221  }
222 
223  vector newDelta = Zero;
224  point nbrCf = oldCf[0];
225 
226  forAll(faces, i)
227  {
228  const label facei = faces[i];
229  const label celli = cells[i];
230 
231  const face f = regionMesh().faces()[facei];
232 
233  newDelta += (deltaV[celli]/mag(sf))*n;
234 
235  vector localDelta = Zero;
236  forAll(f, pti)
237  {
238  const label pointi = f[pti];
239 
240  if
241  (
242  mag((nbrCf - (oldPoints[pointi] + newDelta)) & n)
243  > minDelta
244  )
245  {
246  newPoints[pointi] = oldPoints[pointi] + newDelta;
247  localDelta = newDelta;
248  cellMoveMap[celli] = 1;
249  }
250  }
251  nbrCf = oldCf[i + 1] + localDelta;
252  }
253  // Modify boundary
254  const label bFacei = boundaryFaceOppositeFace_[totalFaceId];
255  const face f = regionMesh().faces()[bFacei];
256  const label celli = cells[cells.size() - 1];
257  newDelta += (deltaV[celli]/mag(sf))*n;
258  forAll(f, pti)
259  {
260  const label pointi = f[pti];
261  if
262  (
263  mag((nbrCf - (oldPoints[pointi] + newDelta)) & n)
264  > minDelta
265  )
266  {
267  newPoints[pointi] = oldPoints[pointi] + newDelta;
268  cellMoveMap[celli] = 1;
269  }
270  }
271  totalFaceId ++;
272  }
273  }
274  // Move points
275  regionMesh().movePoints(newPoints);
276 
277  return tcellMoveMap;
278 }
279 
280 
281 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
282 
283 Foam::regionModels::regionModel1D::regionModel1D
284 (
285  const fvMesh& mesh,
286  const word& regionType
287 )
288 :
289  regionModel(mesh, regionType),
290  boundaryFaceFaces_(),
291  boundaryFaceCells_(),
292  boundaryFaceOppositeFace_(),
293  nLayers_(0),
294  nMagSfPtr_(NULL),
295  moveMesh_(false)
296 {}
297 
298 
299 Foam::regionModels::regionModel1D::regionModel1D
300 (
301  const fvMesh& mesh,
302  const word& regionType,
303  const word& modelName,
304  bool readFields
305 )
306 :
307  regionModel(mesh, regionType, modelName, false),
308  boundaryFaceFaces_(regionMesh().nCells()),
309  boundaryFaceCells_(regionMesh().nCells()),
310  boundaryFaceOppositeFace_(regionMesh().nCells()),
311  nLayers_(0),
312  nMagSfPtr_(NULL),
313  moveMesh_(true)
314 {
315  if (active_)
316  {
317  constructMeshObjects();
318  initialise();
319 
320  if (readFields)
321  {
322  read();
323  }
324  }
325 }
326 
327 
328 Foam::regionModels::regionModel1D::regionModel1D
329 (
330  const fvMesh& mesh,
331  const word& regionType,
332  const word& modelName,
333  const dictionary& dict,
334  bool readFields
335 )
336 :
337  regionModel(mesh, regionType, modelName, dict, readFields),
338  boundaryFaceFaces_(regionMesh().nCells()),
339  boundaryFaceCells_(regionMesh().nCells()),
340  boundaryFaceOppositeFace_(regionMesh().nCells()),
341  nLayers_(0),
342  nMagSfPtr_(NULL),
343  moveMesh_(false)
344 {
345  if (active_)
346  {
347  constructMeshObjects();
348  initialise();
349 
350  if (readFields)
351  {
352  read(dict);
353  }
354  }
355 }
356 
357 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
358 
360 {}
361 
362 
363 // ************************************************************************* //
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:290
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:428
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:137
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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:171
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
Foam::polyBoundaryMesh.
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
labelList f(nPoints)
volScalarField sf(fieldObject, mesh)
label patchi
const Field< PointType > & faceNormals() const
Return face normals for patch.
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: PtrList.H:54
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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.