pointDist.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) 2013-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 
26 #include "pointDist.H"
27 #include "pointEdgeDist.H"
28 #include "pointMesh.H"
29 #include "PointEdgeWave.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const pointMesh& pMesh,
36  const labelHashSet& patchIDs,
37  const pointField& points,
38  const scalar maxDist
39 )
40 :
42  (
43  IOobject
44  (
45  "pointDistance",
46  pMesh.db().time().name(),
47  pMesh.db()
48  ),
49  pMesh,
51  ),
52  points_(points),
53  patchIndices_(patchIDs),
54  maxDist_(maxDist)
55 {
56  correct();
57 }
58 
59 
61 (
62  const pointMesh& pMesh,
63  const labelHashSet& patchIDs,
64  const labelHashSet& zoneIDs,
65  const pointField& points,
66  const scalar maxDist
67 )
68 :
70  (
71  IOobject
72  (
73  "pointDistance",
74  pMesh.db().time().name(),
75  pMesh.db()
76  ),
77  pMesh,
79  ),
80  points_(points),
81  patchIndices_(patchIDs),
82  zoneIndices_(zoneIDs),
83  maxDist_(maxDist)
84 {
85  correct();
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  const pointBoundaryMesh& pbm = mesh().boundary();
100 
101  label nPatchPoints = 0;
102 
103  forAllConstIter(labelHashSet, patchIndices_, iter)
104  {
105  const label patchi = iter.key();
106  nPatchPoints += pbm[patchi].meshPoints().size();
107  }
108 
109  forAllConstIter(labelHashSet, zoneIndices_, iter)
110  {
111  const label zonei = iter.key();
112  nPatchPoints += mesh()().pointZones()[zonei].size();
113  }
114 
115  pointEdgeDist::data pointEdgeData(points_, maxDist_);
116 
117  // Set initial changed points to all the patch points(if patch present)
118  List<pointEdgeDist> patchPointsInfo(nPatchPoints);
119  labelList patchPoints(nPatchPoints);
120  nPatchPoints = 0;
121 
122  // Add the patch points to the patchPointsInfo
123  forAllConstIter(labelHashSet, patchIndices_, iter)
124  {
125  const label patchi = iter.key();
126  const labelList& mp = pbm[patchi].meshPoints();
127 
128  forAll(mp, ppi)
129  {
130  const label meshPointi = mp[ppi];
131  patchPoints[nPatchPoints] = meshPointi;
132  patchPointsInfo[nPatchPoints] = pointEdgeDist
133  (
134  pointEdgeData.points[meshPointi],
135  0
136  );
137  nPatchPoints++;
138  }
139  }
140 
141  // Add the zone points to the patchPointsInfo
142  forAllConstIter(labelHashSet, zoneIndices_, iter)
143  {
144  const label zonei = iter.key();
145  const labelList& zonePoints = mesh()().pointZones()[zonei];
146 
147  forAll(zonePoints, j)
148  {
149  const label meshPointi = zonePoints[j];
150  patchPoints[nPatchPoints] = meshPointi;
151  patchPointsInfo[nPatchPoints] = pointEdgeDist
152  (
153  pointEdgeData.points[meshPointi],
154  0
155  );
156  nPatchPoints++;
157  }
158  }
159 
160  // Current info on points
161  List<pointEdgeDist> allPointInfo(mesh()().nPoints());
162 
163  // Current info on edges
164  List<pointEdgeDist> allEdgeInfo(mesh()().nEdges());
165 
167  <
170  > patchCalc
171  (
172  mesh()(),
173  patchPoints,
174  patchPointsInfo,
175 
176  allPointInfo,
177  allEdgeInfo,
178  mesh().globalData().nTotalPoints(), // max iterations
179  pointEdgeData
180  );
181 
182  pointScalarField& psf = *this;
183 
184  forAll(allPointInfo, pointi)
185  {
186  if (allPointInfo[pointi].valid(pointEdgeData))
187  {
188  psf[pointi] = sqrt(allPointInfo[pointi].distSqr());
189  }
190  else
191  {
192  psf[pointi] = maxDist_;
193  }
194  }
195 
196  forAllConstIter(labelHashSet, zoneIndices_, iter)
197  {
198  const label zonei = iter.key();
199  const labelList& zonePoints = mesh()().pointZones()[zonei];
200 
201  forAll(zonePoints, j)
202  {
203  const label meshPointi = zonePoints[j];
204  psf[meshPointi] = 0;
205  }
206  }
207 }
208 
209 
210 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:89
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Foam::pointBoundaryMesh.
pointDist(const pointMesh &pMesh, const labelHashSet &patchIDs, const pointField &points, const scalar maxDist=rootVGreat)
Construct from mesh and set of patches and points.
Definition: pointDist.C:34
void correct()
Correct for mesh geom/topo changes.
Definition: pointDist.C:97
virtual ~pointDist()
Destructor.
Definition: pointDist.C:91
Class used to pass data into container.
Definition: pointEdgeDist.H:80
const pointField & points
Definition: pointEdgeDist.H:82
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgeDist.H:66
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
label patchi
const pointField & points
label nPoints
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar mp
Proton mass.
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
dimensionedScalar sqrt(const dimensionedScalar &ds)