powerLawLopesdaCosta.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) 2018 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 
30 #include "powerLawLopesdaCosta.H"
31 #include "geometricOneField.H"
32 #include "fvMatrices.H"
33 #include "triSurfaceMesh.H"
34 #include "triSurfaceTools.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  namespace porosityModels
41  {
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const word& modelType,
54  const fvMesh& mesh,
55  const dictionary& dict
56 )
57 :
58  zoneName_(name + ":porousZone")
59 {
60  dictionary coeffs(dict.optionalSubDict(modelType + "Coeffs"));
61 
62  // Vertical direction within porous region
63  vector zDir(coeffs.get<vector>("zDir"));
64 
65  // Span of the search for the cells in the porous region
66  scalar searchSpan(coeffs.get<scalar>("searchSpan"));
67 
68  // Top surface file name defining the extent of the porous region
69  const word topSurfaceFileName(coeffs.get<word>("topSurface"));
70 
71  // List of the ground patches defining the lower surface
72  // of the porous region
73  wordRes groundPatches(coeffs.get<wordRes>("groundPatches"));
74 
75  // Functional form of the porosity surface area per unit volume as a
76  // function of the normalized vertical position
77  autoPtr<Function1<scalar>> SigmaFunc
78  (
79  Function1<scalar>::New("Sigma", dict, &mesh)
80  );
81 
82  // Searchable triSurface for the top of the porous region
83  triSurfaceMesh searchSurf
84  (
85  IOobject
86  (
87  topSurfaceFileName,
88  mesh.time().constant(),
89  "triSurface",
90  mesh.time()
91  )
92  );
93 
94  // Limit the porous cell search to those within the lateral and vertical
95  // extent of the top surface
96 
97  boundBox surfBounds(searchSurf.points());
98  boundBox searchBounds
99  (
100  surfBounds.min() - searchSpan*zDir, surfBounds.max()
101  );
102 
103  const pointField& C = mesh.cellCentres();
104 
105  // List of cells within the porous region
106  labelList porousCells(C.size());
107  label porousCelli = 0;
108 
109  forAll(C, celli)
110  {
111  if (searchBounds.contains(C[celli]))
112  {
113  porousCells[porousCelli++] = celli;
114  }
115  }
116 
117  porousCells.setSize(porousCelli);
118 
119  pointField start(porousCelli);
120  pointField end(porousCelli);
121 
122  forAll(porousCells, porousCelli)
123  {
124  start[porousCelli] = C[porousCells[porousCelli]];
125  end[porousCelli] = start[porousCelli] + searchSpan*zDir;
126  }
127 
128  // Field of distance between the cell centre and the corresponding point
129  // on the porous region top surface
130  scalarField zTop(porousCelli);
131 
132  // Search the porous cells for those with a corresponding point on the
133  // porous region top surface
134  List<pointIndexHit> hitInfo;
135  searchSurf.findLine(start, end, hitInfo);
136 
137  porousCelli = 0;
138 
139  forAll(porousCells, celli)
140  {
141  const pointIndexHit& hit = hitInfo[celli];
142 
143  if (hit.hit())
144  {
145  porousCells[porousCelli] = porousCells[celli];
146 
147  zTop[porousCelli] =
148  (hit.hitPoint() - C[porousCells[porousCelli]]) & zDir;
149 
150  porousCelli++;
151  }
152  }
153 
154  // Resize the porous cells list and height field
155  porousCells.setSize(porousCelli);
156  zTop.setSize(porousCelli);
157 
158  // Create a triSurface for the ground patch(s)
159  triSurface groundSurface
160  (
161  triSurfaceTools::triangulate
162  (
163  mesh.boundaryMesh(),
164  mesh.boundaryMesh().patchSet(groundPatches),
165  searchBounds
166  )
167  );
168 
169  // Combine the ground triSurfaces across all processors
170  if (Pstream::parRun())
171  {
172  List<List<labelledTri>> groundSurfaceProcTris(Pstream::nProcs());
173  List<pointField> groundSurfaceProcPoints(Pstream::nProcs());
174 
175  groundSurfaceProcTris[Pstream::myProcNo()] = groundSurface;
176  groundSurfaceProcPoints[Pstream::myProcNo()] = groundSurface.points();
177 
178  Pstream::gatherList(groundSurfaceProcTris);
179  Pstream::scatterList(groundSurfaceProcTris);
180  Pstream::gatherList(groundSurfaceProcPoints);
181  Pstream::scatterList(groundSurfaceProcPoints);
182 
183  label nTris = 0;
184  forAll(groundSurfaceProcTris, i)
185  {
186  nTris += groundSurfaceProcTris[i].size();
187  }
188 
189  List<labelledTri> groundSurfaceTris(nTris);
190  label trii = 0;
191  label offset = 0;
192  forAll(groundSurfaceProcTris, i)
193  {
194  forAll(groundSurfaceProcTris[i], j)
195  {
196  groundSurfaceTris[trii] = groundSurfaceProcTris[i][j];
197  groundSurfaceTris[trii][0] += offset;
198  groundSurfaceTris[trii][1] += offset;
199  groundSurfaceTris[trii][2] += offset;
200  trii++;
201  }
202  offset += groundSurfaceProcPoints[i].size();
203  }
204 
205  label nPoints = 0;
206  forAll(groundSurfaceProcPoints, i)
207  {
208  nPoints += groundSurfaceProcPoints[i].size();
209  }
210 
211  pointField groundSurfacePoints(nPoints);
212 
213  label pointi = 0;
214  forAll(groundSurfaceProcPoints, i)
215  {
216  forAll(groundSurfaceProcPoints[i], j)
217  {
218  groundSurfacePoints[pointi++] = groundSurfaceProcPoints[i][j];
219  }
220  }
221 
222  groundSurface = triSurface(groundSurfaceTris, groundSurfacePoints);
223  }
224 
225  // Create a searchable triSurface for the ground
226  triSurfaceSearch groundSearch(groundSurface);
227 
228  // Search the porous cells for the corresponding point on the ground
229 
230  start.setSize(porousCelli);
231  end.setSize(porousCelli);
232 
233  forAll(porousCells, porousCelli)
234  {
235  start[porousCelli] = C[porousCells[porousCelli]];
236  end[porousCelli] = start[porousCelli] - searchSpan*zDir;
237  }
238 
239  groundSearch.findLine(start, end, hitInfo);
240 
241  scalarField zBottom(porousCelli);
242 
243  forAll(porousCells, porousCelli)
244  {
245  const pointIndexHit& hit = hitInfo[porousCelli];
246 
247  if (hit.hit())
248  {
249  zBottom[porousCelli] =
250  (C[porousCells[porousCelli]] - hit.hitPoint()) & zDir;
251  }
252  }
253 
254  // Create the normalized height field
255  scalarField zNorm(zBottom/(zBottom + zTop));
256 
257  // Create the porosity surface area per unit volume zone field
258  Sigma_ = SigmaFunc->value(zNorm);
259 
260  // Create the porous region cellZone and add to the mesh cellZones
261 
262  cellZoneMesh& cellZones = const_cast<cellZoneMesh&>(mesh.cellZones());
263 
264  label zoneID = cellZones.findZoneID(zoneName_);
265 
266  if (zoneID == -1)
267  {
268  zoneID = cellZones.size();
269  cellZones.setSize(zoneID + 1);
270 
271  cellZones.set
272  (
273  zoneID,
274  new cellZone
275  (
276  zoneName_,
277  porousCells,
278  zoneID,
279  cellZones
280  )
281  );
282  }
283  else
284  {
286  << "Unable to create porous cellZone " << zoneName_
287  << ": zone already exists"
288  << abort(FatalError);
289  }
290 }
291 
292 
293 Foam::porosityModels::powerLawLopesdaCosta::powerLawLopesdaCosta
294 (
295  const word& name,
296  const word& modelType,
297  const fvMesh& mesh,
298  const dictionary& dict,
299  const word& dummyCellZoneName
300 )
301 :
302  powerLawLopesdaCostaZone(name, modelType, mesh, dict),
304  (
305  name,
306  modelType,
307  mesh,
308  dict,
309  powerLawLopesdaCostaZone::zoneName_
310  ),
311  Cd_(coeffs_.get<scalar>("Cd")),
312  C1_(coeffs_.get<scalar>("C1")),
313  rhoName_(coeffs_.getOrDefault<word>("rho", "rho"))
314 {}
315 
316 
317 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
318 
320 {}
321 
322 
323 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
324 
325 const Foam::scalarField&
327 {
328  return Sigma_;
329 }
330 
331 
333 {}
334 
335 
337 (
338  const volVectorField& U,
339  const volScalarField& rho,
340  const volScalarField& mu,
341  vectorField& force
342 ) const
343 {
344  scalarField Udiag(U.size(), Zero);
345  const scalarField& V = mesh_.V();
346 
347  apply(Udiag, V, rho, U);
348 
349  force = Udiag*U;
350 }
351 
352 
354 (
356 ) const
357 {
358  const vectorField& U = UEqn.psi();
359  const scalarField& V = mesh_.V();
360  scalarField& Udiag = UEqn.diag();
361 
362  if (UEqn.dimensions() == dimForce)
363  {
364  const volScalarField& rho =
365  mesh_.lookupObject<volScalarField>(rhoName_);
366 
367  apply(Udiag, V, rho, U);
368  }
369  else
370  {
371  apply(Udiag, V, geometricOneField(), U);
372  }
373 }
374 
375 
377 (
379  const volScalarField& rho,
380  const volScalarField& mu
381 ) const
382 {
383  const vectorField& U = UEqn.psi();
384  const scalarField& V = mesh_.V();
385  scalarField& Udiag = UEqn.diag();
386 
387  apply(Udiag, V, rho, U);
388 }
389 
390 
392 (
393  const fvVectorMatrix& UEqn,
394  volTensorField& AU
395 ) const
396 {
397  const vectorField& U = UEqn.psi();
398 
399  if (UEqn.dimensions() == dimForce)
400  {
401  const volScalarField& rho =
402  mesh_.lookupObject<volScalarField>(rhoName_);
403 
404  apply(AU, rho, U);
405  }
406  else
407  {
408  apply(AU, geometricOneField(), U);
409  }
410 }
411 
412 
414 {
415  dict_.writeEntry(name_, os);
416 
417  return true;
418 }
419 
420 
421 // ************************************************************************* //
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
powerLawLopesdaCosta.H
Foam::porosityModels::powerLawLopesdaCosta::calcForce
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
Definition: powerLawLopesdaCosta.C:337
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:106
Foam::PointIndexHit::hitPoint
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
Definition: PointIndexHit.H:154
Foam::dimForce
const dimensionSet dimForce
rho
rho
Definition: readInitialConditions.H:88
Foam::porosityModels::powerLawLopesdaCosta::writeData
bool writeData(Ostream &os) const
Write.
Definition: powerLawLopesdaCosta.C:413
Foam::triSurfaceSearch
Helper class to search on triSurface.
Definition: triSurfaceSearch.H:58
Foam::porosityModels::powerLawLopesdaCostaZone::powerLawLopesdaCostaZone
powerLawLopesdaCostaZone(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Constructor.
Definition: powerLawLopesdaCosta.C:51
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::porosityModels::powerLawLopesdaCosta::correct
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
Definition: powerLawLopesdaCosta.C:354
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::PointIndexHit::hit
bool hit() const noexcept
Is there a hit?
Definition: PointIndexHit.H:130
Foam::Field< vector >
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:76
Foam::porosityModels::defineTypeNameAndDebug
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
Foam::triSurfaceMesh::findLine
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
Definition: triSurfaceMesh.C:852
Foam::ZoneMesh< cellZone, polyMesh >
Foam::porosityModels::powerLawLopesdaCostaZone
Definition: powerLawLopesdaCosta.H:89
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::porosityModels::powerLawLopesdaCosta::~powerLawLopesdaCosta
virtual ~powerLawLopesdaCosta()
Destructor.
Definition: powerLawLopesdaCosta.C:319
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::porosityModels::powerLawLopesdaCostaZone::Sigma
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
Definition: powerLawLopesdaCosta.C:326
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::triSurfaceMesh::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: triSurfaceMesh.C:596
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:519
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::porosityModel
Top level model for porosity models.
Definition: porosityModel.H:57
Foam::triSurfaceSearch::findLine
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:333
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
geometricOneField.H
Foam::Vector< scalar >
Foam::apply
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Definition: parcelSelectionDetail.C:82
Foam::List< label >
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::porosityModels::powerLawLopesdaCosta::calcTransformModelData
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
Definition: powerLawLopesdaCosta.C:332
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::porosityModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::GeometricField< vector, fvPatchField, volMesh >
triSurfaceMesh.H
triSurfaceTools.H
Foam::porosityModels::powerLawLopesdaCosta
Variant of the power law porosity model with spatially varying drag coefficient.
Definition: powerLawLopesdaCosta.H:124