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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "regionModel1D.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace regionModels
36{
38}
39}
40
41// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42
43void Foam::regionModels::regionModel1D::constructMeshObjects()
44{
45
46 nMagSfPtr_.reset
47 (
49 (
50 IOobject
51 (
52 "nMagSf",
53 time().timeName(),
54 regionMesh(),
57 ),
58 regionMesh(),
60 )
61 );
62}
63
64
65void Foam::regionModels::regionModel1D::initialise()
66{
67 if (debug)
68 {
69 Pout<< "regionModel1D::initialise()" << endl;
70 }
71
72 // Calculate boundaryFaceFaces and boundaryFaceCells
73
74 DynamicList<label> faceIDs;
75 DynamicList<label> cellIDs;
76
77 label localPyrolysisFacei = 0;
78
79 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
80
81 forAll(intCoupledPatchIDs_, i)
82 {
83 const label patchi = intCoupledPatchIDs_[i];
84 const polyPatch& ppCoupled = rbm[patchi];
85 localPyrolysisFacei += ppCoupled.size();
86 }
87
88 boundaryFaceOppositeFace_.setSize(localPyrolysisFacei);
89 boundaryFaceFaces_.setSize(localPyrolysisFacei);
90 boundaryFaceCells_.setSize(localPyrolysisFacei);
91
92 localPyrolysisFacei = 0;
93
94 forAll(intCoupledPatchIDs_, i)
95 {
96 const label patchi = intCoupledPatchIDs_[i];
97 const polyPatch& ppCoupled = rbm[patchi];
98 forAll(ppCoupled, localFacei)
99 {
100 label facei = ppCoupled.start() + localFacei;
101 label celli = -1;
102 label nCells = 0;
103 do
104 {
105 label ownCelli = regionMesh().faceOwner()[facei];
106 if (ownCelli != celli)
107 {
108 celli = ownCelli;
109 }
110 else
111 {
112 celli = regionMesh().faceNeighbour()[facei];
113 }
114 nCells++;
115 cellIDs.append(celli);
116 const cell& cFaces = regionMesh().cells()[celli];
117 faceIDs.append(facei);
118 label face0 =
119 cFaces.opposingFaceLabel(facei, regionMesh().faces());
120 facei = face0;
121 } while (regionMesh().isInternalFace(facei));
122
123 boundaryFaceOppositeFace_[localPyrolysisFacei] = facei;
124 //faceIDs.remove(); //remove boundary face.
125
126 boundaryFaceFaces_[localPyrolysisFacei].transfer(faceIDs);
127 boundaryFaceCells_[localPyrolysisFacei].transfer(cellIDs);
128
129 localPyrolysisFacei++;
130 nLayers_ = nCells;
131 }
132 }
133 faceIDs.clear();
134 cellIDs.clear();
135
136 surfaceScalarField& nMagSf = nMagSfPtr_();
137
138 surfaceScalarField::Boundary& nMagSfBf = nMagSf.boundaryFieldRef();
139
140 localPyrolysisFacei = 0;
141
142 forAll(intCoupledPatchIDs_, i)
143 {
144 const label patchi = intCoupledPatchIDs_[i];
145 const polyPatch& ppCoupled = rbm[patchi];
146 const vectorField& pNormals = ppCoupled.faceNormals();
147
148 nMagSfBf[patchi] = regionMesh().Sf().boundaryField()[patchi] & pNormals;
149
150 forAll(pNormals, localFacei)
151 {
152 const vector n = pNormals[localFacei];
153 const labelList& faces = boundaryFaceFaces_[localPyrolysisFacei++];
154 forAll(faces, facei)
155 {
156 // facei = 0 is on boundary
157 if (facei > 0)
158 {
159 const label faceID = faces[facei];
160 nMagSf[faceID] = regionMesh().Sf()[faceID] & n;
161 }
162 }
163 }
164 }
165}
166
167
168// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
169
171{
172 if (regionModel::read())
173 {
174 return true;
175 }
176
177 return false;
178}
179
180
182{
184 {
185 moveMesh_.readIfPresent("moveMesh", coeffs_);
186
187 return true;
188 }
189
190 return false;
191}
192
193
195(
196 const scalarList& deltaV,
197 const scalar minDelta
198)
199{
200 tmp<labelField> tcellMoveMap(new labelField(regionMesh().nCells(), Zero));
201 labelField& cellMoveMap = tcellMoveMap.ref();
202
203 if (!moveMesh_)
204 {
205 return cellMoveMap;
206 }
207
208 pointField oldPoints = regionMesh().points();
209 pointField newPoints = oldPoints;
210
211 const polyBoundaryMesh& bm = regionMesh().boundaryMesh();
212
213 label totalFacei = 0;
214 forAll(intCoupledPatchIDs_, localPatchi)
215 {
216 label patchi = intCoupledPatchIDs_[localPatchi];
217 const polyPatch& pp = bm[patchi];
218
219 forAll(pp, patchFacei)
220 {
221 const labelList& faces = boundaryFaceFaces_[totalFacei];
222 const labelList& cells = boundaryFaceCells_[totalFacei];
223 const label oFace = boundaryFaceOppositeFace_[totalFacei];
224
225 const vector n = pp.faceNormals()[patchFacei];
226 const vector sf = pp.faceAreas()[patchFacei];
227
228 List<point> oldCf(faces.size() + 1, Zero);
229 List<bool> frozen(faces.size(), false);
230
231 forAll(faces, i)
232 {
233 oldCf[i] = regionMesh().faceCentres()[faces[i]];
234 }
235
236 oldCf[faces.size()] = regionMesh().faceCentres()[oFace];
237
238 forAll(faces, i)
239 {
240 const label celli = cells[i];
241
242 if (mag(oldCf[i + 1] - oldCf[i]) < minDelta)
243 {
244 frozen[i] = true;
245 cellMoveMap[celli] = 1;
246 }
247 }
248
249 vectorField newDelta(cells.size() + 1, Zero);
250
251 label j = 0;
252 forAll(cells, i)
253 {
254 const label celli = cells[i];
255 newDelta[j+1] = (deltaV[celli]/mag(sf))*n + newDelta[j];
256 j++;
257 }
258
259 // Move the back face first
260 const face of = regionMesh().faces()[oFace];
261 {
262 scalar omagV = mag(newDelta[newDelta.size()-1]);
263
264 if (!frozen[cells.size()-1] && (omagV > ROOTVSMALL))
265 {
266 forAll(of, pti)
267 {
268 const label pointi = of[pti];
269 newPoints[pointi] =
270 oldPoints[pointi] - newDelta[newDelta.size()-1];
271 }
272 }
273 }
274 // Do internal faces
275 for (label i=0; i < faces.size(); i++)
276 {
277 const label facei = faces[i];
278 const face f = regionMesh().faces()[facei];
279
280 scalar magV = mag(newDelta[i]);
281 if (!frozen[i] && magV > 0)
282 {
283 forAll(f, pti)
284 {
285 const label pointi = f[pti];
286 newPoints[pointi] = oldPoints[pointi] - newDelta[i];
287 }
288 }
289 }
290
291 totalFacei++;
292 }
293 }
294
295 // Move points
296 regionMesh().movePoints(newPoints);
297
298 return tcellMoveMap;
299}
300
301
302// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
303
305(
306 const fvMesh& mesh,
307 const word& regionType
308)
309:
310 regionModel(mesh, regionType),
311 boundaryFaceFaces_(),
312 boundaryFaceCells_(),
313 boundaryFaceOppositeFace_(),
314 nLayers_(0),
315 nMagSfPtr_(nullptr),
316 moveMesh_(false)
317{}
318
319
321(
322 const fvMesh& mesh,
323 const word& regionType,
324 const word& modelName,
325 bool readFields
326)
327:
328 regionModel(mesh, regionType, modelName, false),
329 boundaryFaceFaces_(regionMesh().nCells()),
330 boundaryFaceCells_(regionMesh().nCells()),
331 boundaryFaceOppositeFace_(regionMesh().nCells()),
332 nLayers_(0),
333 nMagSfPtr_(nullptr),
334 moveMesh_(true)
335{
336 if (active_)
337 {
338 constructMeshObjects();
339 initialise();
340 if (readFields)
341 {
342 moveMesh_.readIfPresent("moveMesh", coeffs_);
343 }
344 }
345}
346
347
349(
350 const fvMesh& mesh,
351 const word& regionType,
352 const word& modelName,
353 const dictionary& dict,
354 bool readFields
355)
356:
357 regionModel(mesh, regionType, modelName, dict, readFields),
358 boundaryFaceFaces_(regionMesh().nCells()),
359 boundaryFaceCells_(regionMesh().nCells()),
360 boundaryFaceOppositeFace_(regionMesh().nCells()),
361 nLayers_(0),
362 nMagSfPtr_(nullptr),
363 moveMesh_(false)
364{
365 if (active_)
366 {
367 constructMeshObjects();
368 initialise();
369 if (readFields)
370 {
371 moveMesh_.readIfPresent("moveMesh", coeffs_);
372 }
373 }
374}
375
376// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
377
379{}
380
381
382// ************************************************************************* //
label n
labelList cellIDs
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
bool readIfPresent(const word &key, const dictionary &dict)
Update the value of the Switch if it is found in the dictionary.
Definition: Switch.C:335
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual void moveMesh()
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
Base class for 1-D region models.
Definition: regionModel1D.H:57
autoPtr< surfaceScalarField > nMagSfPtr_
Face area magnitude normal to patch.
Definition: regionModel1D.H:95
virtual ~regionModel1D()
Destructor.
Switch moveMesh_
Flag to allow mesh movement.
Definition: regionModel1D.H:98
virtual bool read()
Read control parameters from dictionary.
Base class for region models.
Definition: regionModel.H:63
Switch active_
Active flag.
Definition: regionModel.H:93
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
dictionary coeffs_
Model coefficients dictionary.
Definition: regionModel.H:102
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:149
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< label > labelList
A List of labels.
Definition: List.H:66
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333