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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "powerLawLopesdaCosta.H"
30 #include "geometricOneField.H"
31 #include "fvMatrices.H"
32 #include "triSurfaceMesh.H"
33 #include "triSurfaceTools.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace porosityModels
40  {
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const word& modelType,
53  const fvMesh& mesh,
54  const dictionary& dict
55 )
56 :
57  zoneName_(name + ":porousZone")
58 {
59  dictionary coeffs(dict.optionalSubDict(modelType + "Coeffs"));
60 
61  // Vertical direction within porous region
62  vector zDir(coeffs.lookup("zDir"));
63 
64  // Span of the search for the cells in the porous region
65  scalar searchSpan(coeffs.get<scalar>("searchSpan"));
66 
67  // Top surface file name defining the extent of the porous region
68  const word topSurfaceFileName(coeffs.get<word>("topSurface"));
69 
70  // List of the ground patches defining the lower surface
71  // of the porous region
72  List<wordRe> groundPatches(coeffs.lookup("groundPatches"));
73 
74  // Functional form of the porosity surface area per unit volume as a
75  // function of the normalized vertical position
76  autoPtr<Function1<scalar>> SigmaFunc
77  (
79  );
80 
81  // Searchable triSurface for the top of the porous region
82  triSurfaceMesh searchSurf
83  (
84  IOobject
85  (
86  topSurfaceFileName,
87  mesh.time().constant(),
88  "triSurface",
89  mesh.time()
90  )
91  );
92 
93  // Limit the porous cell search to those within the lateral and vertical
94  // extent of the top surface
95 
96  boundBox surfBounds(searchSurf.points());
97  boundBox searchBounds
98  (
99  surfBounds.min() - searchSpan*zDir, surfBounds.max()
100  );
101 
102  const pointField& C = mesh.cellCentres();
103 
104  // List of cells within the porous region
105  labelList porousCells(C.size());
106  label porousCelli = 0;
107 
108  forAll(C, celli)
109  {
110  if (searchBounds.contains(C[celli]))
111  {
112  porousCells[porousCelli++] = celli;
113  }
114  }
115 
116  porousCells.setSize(porousCelli);
117 
118  pointField start(porousCelli);
119  pointField end(porousCelli);
120 
121  forAll(porousCells, porousCelli)
122  {
123  start[porousCelli] = C[porousCells[porousCelli]];
124  end[porousCelli] = start[porousCelli] + searchSpan*zDir;
125  }
126 
127  // Field of distance between the cell centre and the corresponding point
128  // on the porous region top surface
129  scalarField zTop(porousCelli);
130 
131  // Search the porous cells for those with a corresponding point on the
132  // porous region top surface
133  List<pointIndexHit> hitInfo;
134  searchSurf.findLine(start, end, hitInfo);
135 
136  porousCelli = 0;
137 
138  forAll(porousCells, celli)
139  {
140  const pointIndexHit& hit = hitInfo[celli];
141 
142  if (hit.hit())
143  {
144  porousCells[porousCelli] = porousCells[celli];
145 
146  zTop[porousCelli] =
147  (hit.hitPoint() - C[porousCells[porousCelli]]) & zDir;
148 
149  porousCelli++;
150  }
151  }
152 
153  // Resize the porous cells list and height field
154  porousCells.setSize(porousCelli);
155  zTop.setSize(porousCelli);
156 
157  // Create a triSurface for the ground patch(s)
158  triSurface groundSurface
159  (
160  triSurfaceTools::triangulate
161  (
162  mesh.boundaryMesh(),
163  mesh.boundaryMesh().patchSet(groundPatches),
164  searchBounds
165  )
166  );
167 
168  // Combine the ground triSurfaces across all processors
169  if (Pstream::parRun())
170  {
171  List<List<labelledTri>> groundSurfaceProcTris(Pstream::nProcs());
172  List<pointField> groundSurfaceProcPoints(Pstream::nProcs());
173 
174  groundSurfaceProcTris[Pstream::myProcNo()] = groundSurface;
175  groundSurfaceProcPoints[Pstream::myProcNo()] = groundSurface.points();
176 
177  Pstream::gatherList(groundSurfaceProcTris);
178  Pstream::scatterList(groundSurfaceProcTris);
179  Pstream::gatherList(groundSurfaceProcPoints);
180  Pstream::scatterList(groundSurfaceProcPoints);
181 
182  label nTris = 0;
183  forAll(groundSurfaceProcTris, i)
184  {
185  nTris += groundSurfaceProcTris[i].size();
186  }
187 
188  List<labelledTri> groundSurfaceTris(nTris);
189  label trii = 0;
190  label offset = 0;
191  forAll(groundSurfaceProcTris, i)
192  {
193  forAll(groundSurfaceProcTris[i], j)
194  {
195  groundSurfaceTris[trii] = groundSurfaceProcTris[i][j];
196  groundSurfaceTris[trii][0] += offset;
197  groundSurfaceTris[trii][1] += offset;
198  groundSurfaceTris[trii][2] += offset;
199  trii++;
200  }
201  offset += groundSurfaceProcPoints[i].size();
202  }
203 
204  label nPoints = 0;
205  forAll(groundSurfaceProcPoints, i)
206  {
207  nPoints += groundSurfaceProcPoints[i].size();
208  }
209 
210  pointField groundSurfacePoints(nPoints);
211 
212  label pointi = 0;
213  forAll(groundSurfaceProcPoints, i)
214  {
215  forAll(groundSurfaceProcPoints[i], j)
216  {
217  groundSurfacePoints[pointi++] = groundSurfaceProcPoints[i][j];
218  }
219  }
220 
221  groundSurface = triSurface(groundSurfaceTris, groundSurfacePoints);
222  }
223 
224  // Create a searchable triSurface for the ground
225  triSurfaceSearch groundSearch(groundSurface);
226 
227  // Search the porous cells for the corresponding point on the ground
228 
229  start.setSize(porousCelli);
230  end.setSize(porousCelli);
231 
232  forAll(porousCells, porousCelli)
233  {
234  start[porousCelli] = C[porousCells[porousCelli]];
235  end[porousCelli] = start[porousCelli] - searchSpan*zDir;
236  }
237 
238  groundSearch.findLine(start, end, hitInfo);
239 
240  scalarField zBottom(porousCelli);
241 
242  forAll(porousCells, porousCelli)
243  {
244  const pointIndexHit& hit = hitInfo[porousCelli];
245 
246  if (hit.hit())
247  {
248  zBottom[porousCelli] =
249  (C[porousCells[porousCelli]] - hit.hitPoint()) & zDir;
250  }
251  }
252 
253  // Create the normalized height field
254  scalarField zNorm(zBottom/(zBottom + zTop));
255 
256  // Create the porosity surface area per unit volume zone field
257  Sigma_ = SigmaFunc->value(zNorm);
258 
259  // Create the porous region cellZone and add to the mesh cellZones
260 
261  cellZoneMesh& cellZones = const_cast<cellZoneMesh&>(mesh.cellZones());
262 
263  label zoneID = cellZones.findZoneID(zoneName_);
264 
265  if (zoneID == -1)
266  {
267  zoneID = cellZones.size();
268  cellZones.setSize(zoneID + 1);
269 
270  cellZones.set
271  (
272  zoneID,
273  new cellZone
274  (
275  zoneName_,
276  porousCells,
277  zoneID,
278  cellZones
279  )
280  );
281  }
282  else
283  {
285  << "Unable to create porous cellZone " << zoneName_
286  << ": zone already exists"
287  << abort(FatalError);
288  }
289 }
290 
291 
292 Foam::porosityModels::powerLawLopesdaCosta::powerLawLopesdaCosta
293 (
294  const word& name,
295  const word& modelType,
296  const fvMesh& mesh,
297  const dictionary& dict,
298  const word& dummyCellZoneName
299 )
300 :
301  powerLawLopesdaCostaZone(name, modelType, mesh, dict),
303  (
304  name,
305  modelType,
306  mesh,
307  dict,
308  powerLawLopesdaCostaZone::zoneName_
309  ),
310  Cd_(coeffs_.get<scalar>("Cd")),
311  C1_(coeffs_.get<scalar>("C1")),
312  rhoName_(coeffs_.lookupOrDefault<word>("rho", "rho"))
313 {}
314 
315 
316 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
317 
319 {}
320 
321 
322 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
323 
324 const Foam::scalarField&
326 {
327  return Sigma_;
328 }
329 
330 
332 {}
333 
334 
336 (
337  const volVectorField& U,
338  const volScalarField& rho,
339  const volScalarField& mu,
340  vectorField& force
341 ) const
342 {
343  scalarField Udiag(U.size(), Zero);
344  const scalarField& V = mesh_.V();
345 
346  apply(Udiag, V, rho, U);
347 
348  force = Udiag*U;
349 }
350 
351 
353 (
355 ) const
356 {
357  const vectorField& U = UEqn.psi();
358  const scalarField& V = mesh_.V();
359  scalarField& Udiag = UEqn.diag();
360 
361  if (UEqn.dimensions() == dimForce)
362  {
363  const volScalarField& rho =
364  mesh_.lookupObject<volScalarField>(rhoName_);
365 
366  apply(Udiag, V, rho, U);
367  }
368  else
369  {
370  apply(Udiag, V, geometricOneField(), U);
371  }
372 }
373 
374 
376 (
378  const volScalarField& rho,
379  const volScalarField& mu
380 ) const
381 {
382  const vectorField& U = UEqn.psi();
383  const scalarField& V = mesh_.V();
384  scalarField& Udiag = UEqn.diag();
385 
386  apply(Udiag, V, rho, U);
387 }
388 
389 
391 (
392  const fvVectorMatrix& UEqn,
393  volTensorField& AU
394 ) const
395 {
396  const vectorField& U = UEqn.psi();
397 
398  if (UEqn.dimensions() == dimForce)
399  {
400  const volScalarField& rho =
401  mesh_.lookupObject<volScalarField>(rhoName_);
402 
403  apply(AU, rho, U);
404  }
405  else
406  {
407  apply(AU, geometricOneField(), U);
408  }
409 }
410 
411 
413 {
414  dict_.writeEntry(name_, os);
415 
416  return true;
417 }
418 
419 
420 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:336
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:107
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:54
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:100
Foam::dimForce
const dimensionSet dimForce
rho
rho
Definition: readInitialConditions.H:96
Foam::porosityModels::powerLawLopesdaCosta::writeData
bool writeData(Ostream &os) const
Write.
Definition: powerLawLopesdaCosta.C:412
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:50
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:353
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:55
Foam::label
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:62
Foam::Field< vector >
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:70
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:119
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C: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:933
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:318
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:121
Foam::porosityModels::powerLawLopesdaCostaZone::Sigma
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
Definition: powerLawLopesdaCosta.C:325
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:115
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::triSurfaceMesh::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: triSurfaceMesh.C:677
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name, return -1 if not found.
Definition: ZoneMesh.C:484
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:355
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< wordRe >
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::porosityModels::powerLawLopesdaCosta::calcTransformModelData
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
Definition: powerLawLopesdaCosta.C:331
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
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