powerLawLopesdaCosta.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) 2018-2024 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 
27 #include "powerLawLopesdaCosta.H"
28 #include "geometricOneField.H"
29 #include "fvMatrices.H"
30 #include "triSurfaceMesh.H"
31 #include "triSurfaceTools.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace porosityModels
38  {
41  }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const word& name,
50  const word& modelType,
51  const fvMesh& mesh,
52  const dictionary& dict
53 )
54 :
55  zoneName_(name + ":porousZone")
56 {
57  dictionary coeffs(dict.optionalSubDict(modelType + "Coeffs"));
58 
59  // Vertical direction within porous region
60  vector zDir(coeffs.lookup("zDir"));
61 
62  // Span of the search for the cells in the porous region
63  scalar searchSpan(coeffs.lookup<scalar>("searchSpan"));
64 
65  // Top surface file name defining the extent of the porous region
66  word topSurfaceFileName(coeffs.lookup("topSurface"));
67 
68  // List of the ground patches defining the lower surface
69  // of the porous region
70  List<wordRe> groundPatches(coeffs.lookup("groundPatches"));
71 
72  // Functional form of the porosity surface area per unit volume as a
73  // function of the normalised vertical position
75  (
77  (
78  dict.lookupEntryBackwardsCompatible
79  (
80  {"Av", "Sigma"},
81  false,
82  true
83  ).keyword(),
84  dimLength,
86  dict
87  )
88  );
89 
90  // Searchable triSurface for the top of the porous region
91  triSurfaceMesh searchSurf
92  (
93  IOobject
94  (
95  topSurfaceFileName,
96  mesh.time().constant(),
98  mesh.time()
99  )
100  );
101 
102  // Limit the porous cell search to those within the lateral and vertical
103  // extent of the top surface
104 
105  boundBox surfBounds(searchSurf.points());
106  boundBox searchBounds
107  (
108  surfBounds.min() - searchSpan*zDir, surfBounds.max()
109  );
110 
111  const pointField& C = mesh.cellCentres();
112 
113  // List of cells within the porous region
114  labelList porousCells(C.size());
115  label porousCelli = 0;
116 
117  forAll(C, celli)
118  {
119  if (searchBounds.contains(C[celli]))
120  {
121  porousCells[porousCelli++] = celli;
122  }
123  }
124 
125  porousCells.setSize(porousCelli);
126 
127  pointField start(porousCelli);
128  pointField end(porousCelli);
129 
130  forAll(porousCells, porousCelli)
131  {
132  start[porousCelli] = C[porousCells[porousCelli]];
133  end[porousCelli] = start[porousCelli] + searchSpan*zDir;
134  }
135 
136  // Field of distance between the cell centre and the corresponding point
137  // on the porous region top surface
138  scalarField zTop(porousCelli);
139 
140  // Search the porous cells for those with a corresponding point on the
141  // porous region top surface
142  List<pointIndexHit> hitInfo;
143  searchSurf.findLine(start, end, hitInfo);
144 
145  porousCelli = 0;
146 
147  forAll(porousCells, celli)
148  {
149  const pointIndexHit& hit = hitInfo[celli];
150 
151  if (hit.hit())
152  {
153  porousCells[porousCelli] = porousCells[celli];
154 
155  zTop[porousCelli] =
156  (hit.hitPoint() - C[porousCells[porousCelli]]) & zDir;
157 
158  porousCelli++;
159  }
160  }
161 
162  // Resize the porous cells list and height field
163  porousCells.setSize(porousCelli);
164  zTop.setSize(porousCelli);
165 
166  // Create a triSurface for the ground patch(s)
167  triSurface groundSurface
168  (
170  (
171  mesh.boundaryMesh(),
172  mesh.boundaryMesh().patchSet(groundPatches),
173  searchBounds
174  )
175  );
176 
177  // Combine the ground triSurfaces across all processors
178  if (Pstream::parRun())
179  {
180  List<List<labelledTri>> groundSurfaceProcTris(Pstream::nProcs());
181  List<pointField> groundSurfaceProcPoints(Pstream::nProcs());
182 
183  groundSurfaceProcTris[Pstream::myProcNo()] = groundSurface;
184  groundSurfaceProcPoints[Pstream::myProcNo()] = groundSurface.points();
185 
186  Pstream::gatherList(groundSurfaceProcTris);
187  Pstream::scatterList(groundSurfaceProcTris);
188  Pstream::gatherList(groundSurfaceProcPoints);
189  Pstream::scatterList(groundSurfaceProcPoints);
190 
191  label nTris = 0;
192  forAll(groundSurfaceProcTris, i)
193  {
194  nTris += groundSurfaceProcTris[i].size();
195  }
196 
197  List<labelledTri> groundSurfaceTris(nTris);
198  label trii = 0;
199  label offset = 0;
200  forAll(groundSurfaceProcTris, i)
201  {
202  forAll(groundSurfaceProcTris[i], j)
203  {
204  groundSurfaceTris[trii] = groundSurfaceProcTris[i][j];
205  groundSurfaceTris[trii][0] += offset;
206  groundSurfaceTris[trii][1] += offset;
207  groundSurfaceTris[trii][2] += offset;
208  trii++;
209  }
210  offset += groundSurfaceProcPoints[i].size();
211  }
212 
213  label nPoints = 0;
214  forAll(groundSurfaceProcPoints, i)
215  {
216  nPoints += groundSurfaceProcPoints[i].size();
217  }
218 
219  pointField groundSurfacePoints(nPoints);
220 
221  label pointi = 0;
222  forAll(groundSurfaceProcPoints, i)
223  {
224  forAll(groundSurfaceProcPoints[i], j)
225  {
226  groundSurfacePoints[pointi++] = groundSurfaceProcPoints[i][j];
227  }
228  }
229 
230  groundSurface = triSurface(groundSurfaceTris, groundSurfacePoints);
231  }
232 
233  // Create a searchable triSurface for the ground
234  triSurfaceSearch groundSearch(groundSurface);
235 
236  // Search the porous cells for the corresponding point on the ground
237 
238  start.setSize(porousCelli);
239  end.setSize(porousCelli);
240 
241  forAll(porousCells, porousCelli)
242  {
243  start[porousCelli] = C[porousCells[porousCelli]];
244  end[porousCelli] = start[porousCelli] - searchSpan*zDir;
245  }
246 
247  groundSearch.findLine(start, end, hitInfo);
248 
249  scalarField zBottom(porousCelli);
250 
251  forAll(porousCells, porousCelli)
252  {
253  const pointIndexHit& hit = hitInfo[porousCelli];
254 
255  if (hit.hit())
256  {
257  zBottom[porousCelli] =
258  (C[porousCells[porousCelli]] - hit.hitPoint()) & zDir;
259  }
260  }
261 
262  // Create the normalised height field
263  const scalarField zNorm(zBottom/(zBottom + zTop));
264 
265  // Create the porosity surface area per unit volume zone field
266  Av_ = AvFunc->value(zNorm);
267 
268  // Create the porous region cellZone and add to the mesh cellZones
269  if (!mesh.cellZones().found(zoneName_))
270  {
271  mesh.cellZones().append
272  (
273  new cellZone
274  (
275  zoneName_,
276  porousCells,
277  mesh.cellZones()
278  )
279  );
280  }
281  else
282  {
284  << "Unable to create porous cellZone " << zoneName_
285  << ": zone already exists"
286  << abort(FatalError);
287  }
288 }
289 
290 
292 (
293  const word& name,
294  const word& modelType,
295  const fvMesh& mesh,
296  const dictionary& dict,
297  const word& dummyCellZoneName
298 )
299 :
300  powerLawLopesdaCostaZone(name, modelType, mesh, dict),
302  (
303  name,
304  modelType,
305  mesh,
306  dict,
307  powerLawLopesdaCostaZone::zoneName_
308  ),
309  Cd_(coeffs_.lookup<scalar>("Cd")),
310  C1_(coeffs_.lookup<scalar>("C1")),
311  rhoName_(coeffs_.lookupOrDefault<word>("rho", "rho"))
312 {}
313 
314 
315 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
316 
318 {}
319 
320 
321 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
322 
323 const Foam::scalarField&
325 {
326  return Av_;
327 }
328 
329 
331 {}
332 
333 
335 (
336  const volVectorField& U,
337  const volScalarField& rho,
338  const volScalarField& mu,
339  vectorField& force
340 ) const
341 {
342  scalarField Udiag(U.size(), 0.0);
343  const scalarField& V = mesh_.V();
344 
345  apply(Udiag, V, rho, U);
346 
347  force = Udiag*U;
348 }
349 
350 
352 (
354 ) const
355 {
356  const vectorField& U = UEqn.psi();
357  const scalarField& V = mesh_.V();
358  scalarField& Udiag = UEqn.diag();
359 
360  if (UEqn.dimensions() == dimForce)
361  {
362  const volScalarField& rho =
363  mesh_.lookupObject<volScalarField>(rhoName_);
364 
365  apply(Udiag, V, rho, U);
366  }
367  else
368  {
369  apply(Udiag, V, geometricOneField(), U);
370  }
371 }
372 
373 
375 (
377  const volScalarField& rho,
378  const volScalarField& mu
379 ) const
380 {
381  const vectorField& U = UEqn.psi();
382  const scalarField& V = mesh_.V();
383  scalarField& Udiag = UEqn.diag();
384 
385  apply(Udiag, V, rho, U);
386 }
387 
388 
390 (
391  const fvVectorMatrix& UEqn,
392  volTensorField& AU
393 ) const
394 {
395  const vectorField& U = UEqn.psi();
396 
397  if (UEqn.dimensions() == dimForce)
398  {
399  const volScalarField& rho =
400  mesh_.lookupObject<volScalarField>(rhoName_);
401 
402  apply(AU, rho, U);
403  }
404  else
405  {
406  apply(AU, geometricOneField(), U);
407  }
408 }
409 
410 
412 {
413  os << indent << name_ << endl;
414  dict_.write(os);
415 
416  return true;
417 }
418 
419 
420 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:51
Run-time selectable general function of one variable.
Definition: Function1.H:125
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
const Point & hitPoint() const
Return hit point.
bool hit() const
Is there a hit.
const Field< PointType > & points() const
Return reference to global points.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
void append(ZoneType *) const
Append or update a zone.
Definition: ZoneList.C:176
bool found(const label objectIndex) const
Return true if objectIndex is in any zone.
Definition: ZoneList.C:122
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
bool contains(const point &) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:176
A subset of mesh cells.
Definition: cellZone.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch indices corresponding to the given names.
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
Top level model for porosity models.
Definition: porosityModel.H:57
const word zoneName_
Automatically generated zone name for this porous zone.
const scalarField & Av() const
Return the porosity surface area per unit volume zone field.
scalarField Av_
Porosity surface area per unit volume zone field.
powerLawLopesdaCostaZone(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Constructor.
Variant of the power law porosity model with spatially varying drag coefficient.
powerLawLopesdaCosta(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName)
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
const vectorField & cellCentres() const
static const word & geometryDir()
Return the geometry directory name.
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals,...
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual tmp< pointField > points() const
Get the points that define the surface.
Helper class to search on triSurface.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Triangulated surface description with patch information.
Definition: triSurface.H:69
A class for handling words, derived from string.
Definition: word.H:62
fvVectorMatrix & UEqn
Definition: UEqn.H:11
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
A special matrix type and solver, designed for finite volume solutions of scalar equations.
label nPoints
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
Namespace for OpenFOAM.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimLength
const dimensionSet dimForce
const dimensionSet dimVolume
error FatalError
const dimensionSet dimArea
void offset(label &lst, const label o)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
dictionary dict