uniformDistance.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) 2012-2015 OpenFOAM Foundation
9  Copyright (C) 2018-2022 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 
29 #include "uniformDistance.H"
31 #include "volumeType.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(uniformDistance, 0);
38  addToRunTimeSelectionTable(cellSizeFunction, uniformDistance, dictionary);
39 }
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& initialPointsDict,
46  const searchableSurface& surface,
47  const scalar& defaultCellSize,
48  const labelList regionIndices
49 )
50 :
51  cellSizeFunction
52  (
53  typeName,
54  initialPointsDict,
55  surface,
56  defaultCellSize,
57  regionIndices
58  ),
59  distance_
60  (
61  coeffsDict().get<scalar>("distanceCoeff") * defaultCellSize
62  ),
63  distanceSqr_(sqr(distance_))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 
71 (
72  const pointIndexHit& hitPt,
73  const vector& n,
74  pointField& shapePts,
75  scalarField& shapeSizes
76 ) const
77 {
78  const Foam::point& pt = hitPt.hitPoint();
79 
80  const scalar distanceCellSize =
81  surfaceCellSizeFunction_().interpolate(pt, hitPt.index());
82 
83  if (sideMode_ == rmBothsides)
84  {
85  shapePts.resize(2);
86  shapeSizes.resize(2);
87 
88  shapePts[0] = pt - n*distance_;
89  shapePts[1] = pt + n*distance_;
90 
91  shapeSizes[0] = distanceCellSize;
92  shapeSizes[1] = distanceCellSize;
93  }
94  else if (sideMode_ == smInside)
95  {
96  shapePts.resize(1);
97  shapeSizes.resize(1);
98 
99  shapePts[0] = pt - n*distance_;
100  shapeSizes[0] = distanceCellSize;
101  }
102  else if (sideMode_ == smOutside)
103  {
104  shapePts.resize(1);
105  shapeSizes.resize(1);
106 
107  shapePts[0] = pt - n*distance_;
108  shapeSizes[0] = distanceCellSize;
109  }
110 
111  return false;
112 }
113 
114 
116 (
117  const point& pt,
118  scalar& size
119 ) const
120 {
121  size = 0;
122 
123  List<pointIndexHit> hits;
124 
125  surface_.findNearest
126  (
127  pointField(1, pt),
128  scalarField(1, distanceSqr_),
129  regionIndices_,
130  hits
131  );
132 
133  const pointIndexHit& hitInfo = hits[0];
134 
135  if (hitInfo.hit())
136  {
137  const point& hitPt = hitInfo.point();
138  const label index = hitInfo.index();
139 
140  if (sideMode_ == rmBothsides)
141  {
142  size = surfaceCellSizeFunction_().interpolate(hitPt, index);
143 
144  return true;
145  }
146 
147  // If the nearest point is essentially on the surface, do not do a
148  // getVolumeType calculation, as it will be prone to error.
149  if (hitInfo.point().dist(pt) < snapToSurfaceTol_)
150  {
151  size = surfaceCellSizeFunction_().interpolate(hitPt, index);
152 
153  return true;
154  }
155 
156  pointField ptF(1, pt);
157  List<volumeType> vTL;
158 
159  surface_.getVolumeType(ptF, vTL);
160 
161  bool functionApplied = false;
162 
163  if
164  (
165  sideMode_ == smInside
166  && vTL[0] == volumeType::INSIDE
167  )
168  {
169  size = surfaceCellSizeFunction_().interpolate(hitPt, index);
170 
171  functionApplied = true;
172  }
173  else if
174  (
175  sideMode_ == smOutside
176  && vTL[0] == volumeType::OUTSIDE
177  )
178  {
179  size = surfaceCellSizeFunction_().interpolate(hitPt, index);
180 
181  functionApplied = true;
182  }
183 
184  return functionApplied;
185  }
186 
187  return false;
188 }
189 
190 
192 (
193  const pointField& pts
194 )
195 {
196 // labelHashSet surfaceAlreadyHit(surface_.size());
197 //
198 // for (const Foam::point& pt : pts)
199 // {
200 // List<pointIndexHit> hits;
201 //
202 // surface_.findNearest
203 // (
204 // pointField(1, pt),
205 // scalarField(1, distanceSqr_),
206 // regionIndices_,
207 // hits
208 // );
209 //
210 // if (hits[0].hit() && !surfaceAlreadyHit.found(hits[0].index()))
211 // {
212 // surfaceCellSizeFunction_().refineSurfaceSize(hits[0].index());
213 //
214 // surfaceAlreadyHit.insert(hits[0].index());
215 // }
216 // }
217 //
218 // return true;
219 
220  return false;
221 }
222 
223 
224 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
uniformDistance(const dictionary &initialPointsDict, const searchableSurface &surface, const scalar &defaultCellSize, const labelList regionIndices)
Construct from components.
virtual bool setCellSize(const pointField &pts)
Adapt local cell size. Return true if anything changed.
virtual bool cellSize(const point &pt, scalar &size) const
Modify scalar argument to the cell size specified by function.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
virtual bool sizeLocations(const pointIndexHit &hitPt, const vector &n, pointField &shapePts, scalarField &shapeSizes) const
Macros for easy insertion into run-time selection tables.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
Definition: vector.H:57
A location inside the volume.
Definition: volumeType.H:65
A location outside the volume.
Definition: volumeType.H:66
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineTypeNameAndDebug(combustionModel, 0)
vector point
Point is a vector.
Definition: point.H:37
label n
List< label > labelList
A List of labels.
Definition: List.H:62
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const pointField & pts