surfaceSlipDisplacementPointPatchVectorField.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) 2011-2023 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 
28 #include "Time.H"
29 #include "transformField.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 template<>
41 const char*
43 names[] =
44 {
45  "nearest",
46  "pointNormal",
47  "fixedNormal"
48 };
49 
51  surfaceSlipDisplacementPointPatchVectorField::projectModeNames_;
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void surfaceSlipDisplacementPointPatchVectorField::calcProjection
57 (
58  vectorField& displacement
59 ) const
60 {
61  const polyMesh& mesh = patch().boundaryMesh().mesh()();
62  const pointField& localPoints = patch().localPoints();
63  const labelList& meshPoints = patch().meshPoints();
64 
65  // const scalar deltaT = mesh.time().deltaTValue();
66 
67  // Construct large enough vector in direction of projectDir so
68  // we're guaranteed to hit something.
69 
70  //- Per point projection vector:
71  const scalar projectLen = mag(mesh.bounds().max()-mesh.bounds().min());
72 
73  // For case of fixed projection vector:
74  vector projectVec(0, 0, 0);
75  if (projectMode_ == FIXEDNORMAL)
76  {
77  vector n = projectDir_/mag(projectDir_);
78  projectVec = projectLen*n;
79  }
80 
81 
82  // Get fixed points (bit of a hack)
83  const pointZone* zonePtr = nullptr;
84 
85  if (frozenPointsZone_.size() > 0)
86  {
87  const meshPointZones& pZones = mesh.pointZones();
88 
89  zonePtr = &pZones[frozenPointsZone_];
90 
91  Pout<< "surfaceSlipDisplacementPointPatchVectorField : Fixing all "
92  << zonePtr->size() << " points in pointZone " << zonePtr->name()
93  << endl;
94  }
95 
96  // Get the motionSolver from the mesh
97  const motionSolver& motion =
98  refCast<const fvMeshMovers::motionSolver>
99  (
100  refCast<const fvMesh>(mesh).mover()
101  ).motion();
102 
103  // Get the starting locations from the motionSolver
104  const pointField& points0 =
105  refCast<const displacementMotionSolver>(motion).points0();
106 
107  pointField start(meshPoints.size());
108  forAll(start, i)
109  {
110  start[i] = points0[meshPoints[i]] + displacement[i];
111  }
112 
113  label nNotProjected = 0;
114 
115  if (projectMode_ == NEAREST)
116  {
117  List<pointIndexHit> nearest;
118  labelList hitSurfaces;
120  (
121  start,
122  scalarField(start.size(), sqr(projectLen)),
123  hitSurfaces,
124  nearest
125  );
126 
127  forAll(nearest, i)
128  {
129  if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
130  {
131  // Fixed point. Reset to point0 location.
132  displacement[i] = points0[meshPoints[i]] - localPoints[i];
133  }
134  else if (nearest[i].hit())
135  {
136  displacement[i] =
137  nearest[i].hitPoint()
138  - points0[meshPoints[i]];
139  }
140  else
141  {
142  nNotProjected++;
143 
144  if (debug)
145  {
146  Pout<< " point:" << meshPoints[i]
147  << " coord:" << localPoints[i]
148  << " did not find any surface within " << projectLen
149  << endl;
150  }
151  }
152  }
153  }
154  else
155  {
156  // Do tests on all points. Combine later on.
157 
158  // 1. Check if already on surface
159  List<pointIndexHit> nearest;
160  {
161  labelList nearestSurface;
163  (
164  start,
165  scalarField(start.size(), sqr(small)),
166  nearestSurface,
167  nearest
168  );
169  }
170 
171  // 2. intersection. (combined later on with information from nearest
172  // above)
173  vectorField projectVecs(start.size(), projectVec);
174 
175  if (projectMode_ == POINTNORMAL)
176  {
177  projectVecs = projectLen*patch().pointNormals();
178  }
179 
180  // Knock out any wedge component
181  scalarField offset(start.size(), 0.0);
182  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
183  {
184  forAll(offset, i)
185  {
186  offset[i] = start[i][wedgePlane_];
187  start[i][wedgePlane_] = 0;
188  projectVecs[i][wedgePlane_] = 0;
189  }
190  }
191 
192  List<pointIndexHit> rightHit;
193  {
194  labelList rightSurf;
196  (
197  start,
198  start+projectVecs,
199  rightSurf,
200  rightHit
201  );
202  }
203 
204  List<pointIndexHit> leftHit;
205  {
206  labelList leftSurf;
208  (
209  start,
210  start-projectVecs,
211  leftSurf,
212  leftHit
213  );
214  }
215 
216  // 3. Choose either -fixed, nearest, right, left.
217  forAll(displacement, i)
218  {
219  if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
220  {
221  // Fixed point. Reset to point0 location.
222  displacement[i] = points0[meshPoints[i]] - localPoints[i];
223  }
224  else if (nearest[i].hit())
225  {
226  // Found nearest.
227  displacement[i] =
228  nearest[i].hitPoint()
229  - points0[meshPoints[i]];
230  }
231  else
232  {
233  pointIndexHit interPt;
234 
235  if (rightHit[i].hit())
236  {
237  if (leftHit[i].hit())
238  {
239  if
240  (
241  magSqr(rightHit[i].hitPoint()-start[i])
242  < magSqr(leftHit[i].hitPoint()-start[i])
243  )
244  {
245  interPt = rightHit[i];
246  }
247  else
248  {
249  interPt = leftHit[i];
250  }
251  }
252  else
253  {
254  interPt = rightHit[i];
255  }
256  }
257  else
258  {
259  if (leftHit[i].hit())
260  {
261  interPt = leftHit[i];
262  }
263  }
264 
265 
266  if (interPt.hit())
267  {
268  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
269  {
270  interPt.rawPoint()[wedgePlane_] += offset[i];
271  }
272  displacement[i] = interPt.rawPoint()-points0[meshPoints[i]];
273  }
274  else
275  {
276  nNotProjected++;
277 
278  if (debug)
279  {
280  Pout<< " point:" << meshPoints[i]
281  << " coord:" << localPoints[i]
282  << " did not find any intersection between"
283  << " ray from " << start[i]-projectVecs[i]
284  << " to " << start[i]+projectVecs[i] << endl;
285  }
286  }
287  }
288  }
289  }
290 
291  reduce(nNotProjected, sumOp<label>());
292 
293  if (nNotProjected > 0)
294  {
295  Info<< "surfaceSlipDisplacement :"
296  << " on patch " << patch().name()
297  << " did not project " << nNotProjected
298  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
299  << " points." << endl;
300  }
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
305 
308 (
309  const pointPatch& p,
311  const dictionary& dict
312 )
313 :
315  surfacesDict_(dict.subDict("geometry")),
316  projectMode_(projectModeNames_.read(dict.lookup("projectMode"))),
317  projectDir_(dict.lookup("projectDirection")),
318  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
319  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
320 {}
321 
322 
325 (
327  const pointPatch& p,
329  const pointPatchFieldMapper&
330 )
331 :
333  surfacesDict_(ppf.surfacesDict_),
334  projectMode_(ppf.projectMode_),
335  projectDir_(ppf.projectDir_),
336  wedgePlane_(ppf.wedgePlane_),
337  frozenPointsZone_(ppf.frozenPointsZone_)
338 {}
339 
340 
343 (
346 )
347 :
348  pointPatchVectorField(ppf, iF),
349  surfacesDict_(ppf.surfacesDict_),
350  projectMode_(ppf.projectMode_),
351  projectDir_(ppf.projectDir_),
352  wedgePlane_(ppf.wedgePlane_),
353  frozenPointsZone_(ppf.frozenPointsZone_)
354 {}
355 
356 
357 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
358 
359 const searchableSurfaces&
361 {
362  if (surfacesPtr_.empty())
363  {
364  surfacesPtr_.reset
365  (
367  (
368  IOobject
369  (
370  "abc", // dummy name
371  db().time().constant(),
373  db().time(),
376  ),
377  surfacesDict_,
378  true // use single region naming shortcut
379  )
380  );
381  }
382  return surfacesPtr_();
383 }
384 
385 
387 (
388  const Pstream::commsTypes commsType
389 )
390 {
391  vectorField displacement(this->patchInternalField());
392 
393  // Calculate displacement to project points onto surface
394  calcProjection(displacement);
395 
396  // Get internal field to insert values into
397  Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
398 
399  // setInternalField(iF, motionU);
400  setInternalField(iF, displacement);
401 
403 }
404 
405 
407 {
409  writeEntry(os, "geometry", surfacesDict_);
410  writeEntry(os, "projectMode", projectModeNames_[projectMode_]);
411  writeEntry(os, "projectDirection", projectDir_);
412  writeEntry(os, "wedgePlane", wedgePlane_);
413  if (frozenPointsZone_ != word::null)
414  {
415  writeEntry(os, "frozenPointsZone", frozenPointsZone_);
416  }
417 }
418 
419 
420 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421 
423 (
426 );
427 
428 
429 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
430 
431 } // End namespace Foam
432 
433 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
commsTypes
Types of communications.
Definition: UPstream.H:65
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:54
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const pointMesh & mesh() const
Return the mesh reference.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
virtual void write(Ostream &) const
Write.
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
const Field< Type > & primitiveField() const
Return internal field reference.
const objectRegistry & db() const
Return local objectRegistry.
const pointPatch & patch() const
Return patch.
void setInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
virtual const word & name() const =0
Return name.
virtual const vectorField & pointNormals() const =0
Return point normals.
virtual const labelList & meshPoints() const =0
Return mesh points.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: pointPatch.H:133
virtual const vectorField & localPoints() const =0
Return mesh points.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:439
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:409
static const word & geometryDir()
Return the geometry directory name.
Container for searchableSurfaces.
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
Displacement follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Update the patch field.
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
surfaceSlipDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
IOporosityModelList pZones(mesh)
Namespace for OpenFOAM.
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
MeshZones< pointZone, polyMesh > meshPointZones
A MeshZones with the type pointZone.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
void offset(label &lst, const label o)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Spatial transformation functions for primitive fields.