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-2022 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*
42 NamedEnum<surfaceSlipDisplacementPointPatchVectorField::projectMode, 3>::
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 )
312 :
313  pointPatchVectorField(p, iF),
314  projectMode_(NEAREST),
315  projectDir_(Zero),
316  wedgePlane_(-1)
317 {}
318 
319 
322 (
323  const pointPatch& p,
325  const dictionary& dict
326 )
327 :
328  pointPatchVectorField(p, iF, dict),
329  surfacesDict_(dict.subDict("geometry")),
330  projectMode_(projectModeNames_.read(dict.lookup("projectMode"))),
331  projectDir_(dict.lookup("projectDirection")),
332  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
333  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
334 {}
335 
336 
339 (
341  const pointPatch& p,
343  const pointPatchFieldMapper&
344 )
345 :
346  pointPatchVectorField(p, iF),
347  surfacesDict_(ppf.surfacesDict_),
348  projectMode_(ppf.projectMode_),
349  projectDir_(ppf.projectDir_),
350  wedgePlane_(ppf.wedgePlane_),
351  frozenPointsZone_(ppf.frozenPointsZone_)
352 {}
353 
354 
357 (
360 )
361 :
362  pointPatchVectorField(ppf, iF),
363  surfacesDict_(ppf.surfacesDict_),
364  projectMode_(ppf.projectMode_),
365  projectDir_(ppf.projectDir_),
366  wedgePlane_(ppf.wedgePlane_),
367  frozenPointsZone_(ppf.frozenPointsZone_)
368 {}
369 
370 
371 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
372 
373 const searchableSurfaces&
375 {
376  if (surfacesPtr_.empty())
377  {
378  surfacesPtr_.reset
379  (
381  (
382  IOobject
383  (
384  "abc", // dummy name
385  db().time().constant(),
387  db().time(),
390  ),
391  surfacesDict_,
392  true // use single region naming shortcut
393  )
394  );
395  }
396  return surfacesPtr_();
397 }
398 
399 
401 (
402  const Pstream::commsTypes commsType
403 )
404 {
405  vectorField displacement(this->patchInternalField());
406 
407  // Calculate displacement to project points onto surface
408  calcProjection(displacement);
409 
410  // Get internal field to insert values into
411  Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
412 
413  // setInternalField(iF, motionU);
414  setInternalField(iF, displacement);
415 
417 }
418 
419 
421 {
423  writeEntry(os, "geometry", surfacesDict_);
424  writeEntry(os, "projectMode", projectModeNames_[projectMode_]);
425  writeEntry(os, "projectDirection", projectDir_);
426  writeEntry(os, "wedgePlane", wedgePlane_);
427  if (frozenPointsZone_ != word::null)
428  {
429  writeEntry(os, "frozenPointsZone", frozenPointsZone_);
430  }
431 }
432 
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
437 (
440 );
441 
442 
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444 
445 } // End namespace Foam
446 
447 // ************************************************************************* //
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Displacement follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Update the patch field.
commsTypes
Types of communications.
Definition: UPstream.H:64
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Foam::pointPatchFieldMapper.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Virtual base class for mesh motion solver.
Definition: motionSolver.H:56
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
fvMesh & mesh
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
const Field< vector > & primitiveField() const
Return internal field reference.
Macros for easy insertion into run-time selection tables.
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:489
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
Spatial transformation functions for primitive fields.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
const objectRegistry & db() const
Return local objectRegistry.
static const word & geometryDir()
Return the geometry directory name.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const pointPatch & patch() const
Return patch.
bool hit() const
Is there a hit.
static const word null
An empty word.
Definition: word.H:77
Container for searchableSurfaces.
static const zero Zero
Definition: zero.H:97
virtual const word & name() const =0
Return name.
const word & name() const
Return name.
Definition: zone.H:147
label whichPoint(const label globalPointID) const
Helper function to re-direct to zone::localID(...)
Definition: pointZone.C:126
tmp< Field< vector > > patchInternalField() const
Return field created from appropriate internal field values.
const Point & rawPoint() const
Return point with no checking.
const pointMesh & mesh() const
Return the mesh reference.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
pointPatchField< vector > pointPatchVectorField
virtual const vectorField & pointNormals() const =0
Return point normals.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: pointPatch.H:133
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual const vectorField & localPoints() const =0
Return mesh points.
virtual void write(Ostream &) const
Write.
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:60
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:459
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list...
Definition: pointZone.H:62
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void setInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelList &meshPoints) const
Given the internal field and a patch field,.
messageStream Info
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
surfaceSlipDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
virtual const labelList & meshPoints() const =0
Return mesh points.
IOporosityModelList pZones(mesh)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864