wallDistAddressing.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) 2023 OpenCFD Ltd.
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 
28 #include "FaceCellWave.H"
29 #include "wallDistAddressing.H"
30 #include "wallPointAddressing.H"
32 #include "surfaceFields.H"
33 #include "wallPolyPatch.H"
34 #include "patchDistMethod.H"
35 #include "OBJstream.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
48 (
50  const List<wallPointAddressing>& allCellInfo,
51  const List<wallPointAddressing>& allFaceInfo,
53 ) const
54 {
55  // Extract volField. See patchWave::getValues
56 
57  label nIllegal = 0;
58  forAll(allCellInfo, celli)
59  {
60  const scalar dist = allCellInfo[celli].distSqr();
61 
62  if (allCellInfo[celli].valid(wave.data()))
63  {
64  y[celli] = Foam::sqrt(dist);
65  }
66  else
67  {
68  y[celli] = dist;
69  ++nIllegal;
70  }
71  }
72 
73  // Copy boundary values
74  auto& bfld = y.boundaryFieldRef();
75 
76  for (auto& pfld : bfld)
77  {
78  scalarField patchField(pfld.size(), Zero);
79 
80  forAll(pfld, patchFacei)
81  {
82  const label meshFacei = pfld.patch().start() + patchFacei;
83  const scalar dist = allFaceInfo[meshFacei].distSqr();
84 
85  if (allFaceInfo[meshFacei].valid(wave.data()))
86  {
87  // Adding SMALL to avoid problems with /0 in the turbulence
88  // models
89  patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
90  }
91  else
92  {
93  patchField[patchFacei] = dist;
94  ++nIllegal;
95  }
96  }
97  pfld = patchField;
98  }
99 
100  return nIllegal;
101 }
102 
103 
105 (
106  const label item,
107  const labelPair& data,
108  label& untransformi,
109  label& transformi,
110  labelPairList& transformedWallInfo
111 )
112 {
113  const auto& gt = mesh_.globalData().globalTransforms();
114 
115  if (gt.transformIndex(data) == gt.nullTransformIndex())
116  {
117  // Store item and global index of untransformed data
118  const label proci = gt.processor(data);
119  const label index = gt.index(data);
120  untransformedItems_[untransformi] = item;
121  untransformedSlots_[untransformi++] =
122  globalWallsPtr_().toGlobal(proci, index);
123  }
124  else
125  {
126  // Store item and transformed data
127  transformedItems_[transformi] = item;
128  transformedWallInfo[transformi++] = data;
129  }
130 }
131 
132 
134 {
135  y = dimensionedScalar("yWall", dimLength, GREAT);
136 
137  const auto& C = mesh_.C();
138  const auto& gt = mesh_.globalData().globalTransforms();
139 
140  label nWalls = 0;
141  for (const label patchi : patchIDs_)
142  {
143  nWalls += C.boundaryField()[patchi].size();
144  }
145 
146  globalWallsPtr_.reset(new globalIndex(nWalls));
147  const globalIndex& globalWalls = globalWallsPtr_();
148 
149  labelList seedFaces(nWalls);
150  List<wallPointAddressing> seedInfo(nWalls);
151 
152 
153  nWalls = 0;
154  for (const label patchi : patchIDs_)
155  {
156  const auto& fc = C.boundaryField()[patchi];
157 
158  forAll(fc, patchFacei)
159  {
160  seedFaces[nWalls] = fc.patch().start()+patchFacei;
161  seedInfo[nWalls] = wallPointAddressing
162  (
163  fc[patchFacei],
164  gt.encode
165  (
167  nWalls,
168  gt.nullTransformIndex()
169  ),
170  scalar(0.0)
171  );
172  nWalls++;
173  }
174  }
175 
176  List<wallPointAddressing> allCellInfo(mesh_.nCells());
177  // Initialise the passive data (since no mesh reference in the constructor
178  // of wallPointAddressing)
179  forAll(allCellInfo, celli)
180  {
181  allCellInfo[celli].data() = gt.encode
182  (
184  celli,
185  gt.nullTransformIndex()
186  );
187  }
188  List<wallPointAddressing> allFaceInfo(mesh_.nFaces());
189  forAll(allFaceInfo, facei)
190  {
191  allFaceInfo[facei].data() = gt.encode
192  (
194  facei,
195  gt.nullTransformIndex()
196  );
197  }
198 
199  // Propagate information inwards
200  FaceCellWave<wallPointAddressing> wave
201  (
202  mesh_,
203  seedFaces,
204  seedInfo,
205  allFaceInfo,
206  allCellInfo,
207  mesh_.globalData().nTotalCells()+1
208  );
209 
210 
211  // Extract wall distance
212  // ~~~~~~~~~~~~~~~~~~~~~
213 
214  getValues(wave, allCellInfo, allFaceInfo, y);
215 
216 
217  const labelHashSet patchSet(patchIDs_);
218 
219  // Correct wall cells for true distance
220  // - store for corrected cells the nearest wall face
221  // - update distance field (y) for these cells
222  Map<label> cellToWallFace;
223  if (correctWalls_)
224  {
225  cellToWallFace.reserve(nWalls);
226  correctBoundaryFaceCells
227  (
228  patchSet,
229  y,
230  cellToWallFace
231  );
232 
233  correctBoundaryPointCells
234  (
235  patchSet,
236  y,
237  cellToWallFace
238  );
239  }
240 
241  // Make sure boundary values are up to date
242  y.correctBoundaryConditions();
243 
244 
245  // Extract all addressing
246  // ~~~~~~~~~~~~~~~~~~~~~~
247  // - untransformed : -local cell/boundary, -nearest (global) wall
248  // - transformed : -local cell/boundary, -nearest (global) wall
249 
250  label untransformi = 0;
251  untransformedItems_.resize(mesh_.nCells()+mesh_.nBoundaryFaces());
252  untransformedSlots_.resize(untransformedItems_.size());
253 
254  label transformi = 0;
255  transformedItems_.resize(mesh_.nCells()+mesh_.nBoundaryFaces());
256  labelPairList transformedWallInfo(transformedItems_.size());
257 
258  for (label celli = 0; celli < mesh_.nCells(); celli++)
259  {
260  const auto cellFnd = cellToWallFace.find(celli);
261  if (cellFnd)
262  {
263  const label wallFacei = cellFnd();
264  const auto& data = allFaceInfo[wallFacei].data();
265  addItem(celli, data, untransformi, transformi, transformedWallInfo);
266  }
267  else if (allCellInfo[celli].valid(wave.data()))
268  {
269  const auto& data = allCellInfo[celli].data();
270  addItem(celli, data, untransformi, transformi, transformedWallInfo);
271  }
272  }
273  untransformedPatchStarts_.resize(mesh_.boundary().size()+1);
274  transformedPatchStarts_.resize(mesh_.boundary().size()+1);
275  for (const auto& pp : mesh_.boundary())
276  {
277  untransformedPatchStarts_[pp.index()] = untransformi;
278  transformedPatchStarts_[pp.index()] = transformi;
279 
280  forAll(pp, patchFacei)
281  {
282  const label facei = pp.start()+patchFacei;
283  if (allFaceInfo[facei].valid(wave.data()))
284  {
285  const auto& data = allFaceInfo[facei].data();
286  addItem
287  (
288  facei,
289  data,
290  untransformi,
291  transformi,
292  transformedWallInfo
293  );
294  }
295  }
296  }
297 
298  untransformedItems_.resize(untransformi);
299  untransformedSlots_.resize(untransformi);
300  untransformedPatchStarts_.back() = untransformi;
301  transformedItems_.resize(transformi);
302  transformedWallInfo.resize(transformi);
303  transformedPatchStarts_.back() = transformi;
304 
305  if (debug)
306  {
307  Pout<< typeName
308  << " : untransformed:" << untransformi
309  << " transformed:" << transformi
310  << " out of nWalls:" << nWalls
311  << " untransformStart:" << flatOutput(untransformedPatchStarts_)
312  << " transformStart:" << flatOutput(transformedPatchStarts_)
313  << endl;
314  }
315 
316  List<Map<label>> compactMap;
317  mapPtr_.reset
318  (
319  new mapDistribute
320  (
321  globalWalls, // globalIndex
322  untransformedSlots_, // untransformedElements
323 
324  gt, // globalIndexAndTransform
325  transformedWallInfo, // transformedElements
326  transformedSlots_, // transformedIndices
327  compactMap
328  )
329  );
330 
331  if (debug & 2)
332  {
333  // Use the obtained map to write nearest
334  {
335  OBJstream os(mesh_.polyMesh::path()/"nearest.obj");
336  Pout<< typeName
337  << " : writing line from cell/face to nearest wall face to "
338  << os.name() << endl;
339 
340  // Seed all wall faces with centroid and override all other
341  // cells/faces with that of nearest wall
342  volVectorField wallCentre(mesh_.C());
343  this->map(wallCentre, mapDistribute::transformPosition());
344 
345  // Seed all wall faces with normal and override all other
346  // cells/faces with that of nearest wall
348  (
349  IOobject
350  (
351  "n" & patchTypeName_,
352  mesh_.time().timeName(),
353  mesh_,
355  ),
356  mesh_,
358  patchDistMethod::patchTypes<vector>(mesh_, patchSet)
359  );
360  // Calculate normal
361  for (const label patchi : patchIDs_)
362  {
363  auto& pnf = wallNormal.boundaryFieldRef()[patchi];
364  pnf == pnf.patch().nf();
365  }
366  this->map(wallNormal, mapDistribute::transform());
367 
368  forAll(wallCentre, celli)
369  {
370  const point& cc = mesh_.C()[celli];
371  const point& wallC = wallCentre[celli];
372  const vector& wallN = wallNormal[celli];
373  os.write(linePointRef(cc, wallC), wallN, wallN);
374  }
375  }
376  }
377 }
378 
379 
380 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
381 
383 (
384  const fvMesh& mesh,
385  const bool correctWalls,
386  const label updateInterval
387 )
388 :
389  // Register as "wallDistAddressing"
390  MeshObject<fvMesh, Foam::UpdateableMeshObject, wallDistAddressing>(mesh),
391  cellDistFuncs(mesh),
392  patchIDs_(mesh.boundaryMesh().findPatchIDs<wallPolyPatch>().sortedToc()),
393  patchTypeName_("wall"),
394  updateInterval_(updateInterval),
395  correctWalls_(correctWalls),
396  requireUpdate_(true),
397  y_
398  (
399  IOobject
400  (
401  "y" & patchTypeName_,
402  mesh.time().timeName(),
403  mesh.thisDb(),
404  IOobjectOption::NO_REGISTER
405  ),
406  mesh,
407  dimensionedScalar("y", dimLength, GREAT) // Same as wallDistData
408  )
409 {
410  movePoints();
411 }
412 
413 
415 (
416  const word& patchTypeName,
417  const fvMesh& mesh,
418  const labelList& patchIDs,
419  const bool correctWalls,
420  const label updateInterval
421 )
422 :
424  (
425  patchTypeName,
426  mesh
427  ),
429  patchIDs_(patchIDs),
430  patchTypeName_(patchTypeName),
431  updateInterval_(updateInterval),
432  correctWalls_(correctWalls),
433  requireUpdate_(true),
434  y_
435  (
436  IOobject
437  (
438  "y" & patchTypeName_,
439  mesh.time().timeName(),
440  mesh.thisDb(),
441  IOobjectOption::NO_REGISTER
442  ),
443  mesh,
444  dimensionedScalar("y", dimLength, GREAT) // Same as wallDistData
445  )
446 {
447  movePoints();
448 }
449 
450 
451 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
454 {} // mapDistribute was forward declared
455 
456 
457 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
458 
460 {
461  if
462  (
463  (updateInterval_ > 0)
464  && ((mesh_.time().timeIndex() % updateInterval_) == 0)
465  )
466  {
467  requireUpdate_ = true;
468  }
469 
470  if (requireUpdate_)
471  {
472  DebugInfo<< "Updating wall distance" << endl;
473 
474  requireUpdate_ = false;
475 
476  correct(y_);
477 
478  return true;
479  }
480 
481  return false;
482 }
483 
484 
485 void Foam::wallDistAddressing::updateMesh(const mapPolyMesh& mpm)
486 {
487  // Force update if performing topology change
488  // Note: needed?
489  // - field would have been mapped, so if using updateInterval option (!= 1)
490  // live with error associated of not updating and use mapped values?
491  requireUpdate_ = true;
492  movePoints();
493 }
494 
495 
496 // ************************************************************************* //
Foam::surfaceFields.
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field. Unvisited.
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:128
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
wallDistAddressing(const wallDistAddressing &)=delete
No copy construct.
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:272
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:206
label getValues(const FaceCellWave< wallPointAddressing > &wave, const List< wallPointAddressing > &allCellInfo, const List< wallPointAddressing > &allFaceInfo, volScalarField &y) const
Extract FaceCellWave data.
void correct(volScalarField &y)
Extract nearest-patch distance data.
const TrackingData & data() const noexcept
Additional data to be passed into container.
Definition: FaceCellWave.H:506
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:59
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
word timeName
Definition: getTimeIndex.H:3
scalar y
dynamicFvMesh & mesh
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual ~wallDistAddressing()
Destructor.
Vector< scalar > vector
Definition: vector.H:57
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
#define DebugInfo
Report an information message using Foam::Info.
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
vector point
Point is a vector.
Definition: point.H:37
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
vector wallNormal(Zero)
List< label > labelList
A List of labels.
Definition: List.H:62
virtual bool movePoints()
Update the y-field when the mesh moves.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
void addItem(const label item, const labelPair &data, label &untransformi, label &transformi, labelPairList &transformedWallInfo)
Store nearest-data to cell or boundary face.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
Variant of wallDist that uses meshWave and stores the addressing.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127