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-2025 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 const NamedEnum<surfaceSlipDisplacementPointPatchVectorField::projectMode, 3>
41 surfaceSlipDisplacementPointPatchVectorField::projectModeNames_
42 {
43  "nearest",
44  "pointNormal",
45  "fixedNormal"
46 };
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void surfaceSlipDisplacementPointPatchVectorField::calcProjection
52 (
53  vectorField& displacement
54 ) const
55 {
56  const polyMesh& mesh = patch().boundaryMesh().mesh()();
57  const pointField& localPoints = patch().localPoints();
58  const labelList& meshPoints = patch().meshPoints();
59 
60  // const scalar deltaT = mesh.time().deltaTValue();
61 
62  // Construct large enough vector in direction of projectDir so
63  // we're guaranteed to hit something.
64 
65  //- Per point projection vector:
66  const scalar projectLen = mag(mesh.bounds().max()-mesh.bounds().min());
67 
68  // For case of fixed projection vector:
69  vector projectVec(0, 0, 0);
70  if (projectMode_ == FIXEDNORMAL)
71  {
72  vector n = projectDir_/mag(projectDir_);
73  projectVec = projectLen*n;
74  }
75 
76 
77  // Get fixed points (bit of a hack)
78  const pointZone* zonePtr = nullptr;
79 
80  if (frozenPointsZone_.size() > 0)
81  {
82  const pointZoneList& pZones = mesh.pointZones();
83 
84  zonePtr = &pZones[frozenPointsZone_];
85 
86  Pout<< "surfaceSlipDisplacementPointPatchVectorField : Fixing all "
87  << zonePtr->size() << " points in pointZone " << zonePtr->name()
88  << endl;
89  }
90 
91  // Get the motionSolver from the mesh
92  const motionSolver& motion =
93  refCast<const fvMeshMovers::motionSolver>
94  (
95  refCast<const fvMesh>(mesh).mover()
96  ).motion();
97 
98  // Get the starting locations from the motionSolver
99  const pointField& points0 =
100  refCast<const displacementMotionSolver>(motion).points0();
101 
102  pointField start(meshPoints.size());
103  forAll(start, i)
104  {
105  start[i] = points0[meshPoints[i]] + displacement[i];
106  }
107 
108  label nNotProjected = 0;
109 
110  if (projectMode_ == NEAREST)
111  {
112  List<pointIndexHit> nearest;
113  labelList hitSurfaces;
115  (
116  start,
117  scalarField(start.size(), sqr(projectLen)),
118  hitSurfaces,
119  nearest
120  );
121 
122  forAll(nearest, i)
123  {
124  if (zonePtr && (zonePtr->localIndex(meshPoints[i]) >= 0))
125  {
126  // Fixed point. Reset to point0 location.
127  displacement[i] = points0[meshPoints[i]] - localPoints[i];
128  }
129  else if (nearest[i].hit())
130  {
131  displacement[i] =
132  nearest[i].hitPoint()
133  - points0[meshPoints[i]];
134  }
135  else
136  {
137  nNotProjected++;
138 
139  if (debug)
140  {
141  Pout<< " point:" << meshPoints[i]
142  << " coord:" << localPoints[i]
143  << " did not find any surface within " << projectLen
144  << endl;
145  }
146  }
147  }
148  }
149  else
150  {
151  // Do tests on all points. Combine later on.
152 
153  // 1. Check if already on surface
154  List<pointIndexHit> nearest;
155  {
156  labelList nearestSurface;
158  (
159  start,
160  scalarField(start.size(), sqr(small)),
161  nearestSurface,
162  nearest
163  );
164  }
165 
166  // 2. intersection. (combined later on with information from nearest
167  // above)
168  vectorField projectVecs(start.size(), projectVec);
169 
170  if (projectMode_ == POINTNORMAL)
171  {
172  projectVecs = projectLen*patch().pointNormals();
173  }
174 
175  // Knock out any wedge component
176  scalarField offset(start.size(), 0.0);
177  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
178  {
179  forAll(offset, i)
180  {
181  offset[i] = start[i][wedgePlane_];
182  start[i][wedgePlane_] = 0;
183  projectVecs[i][wedgePlane_] = 0;
184  }
185  }
186 
187  List<pointIndexHit> rightHit;
188  {
189  labelList rightSurf;
191  (
192  start,
193  start+projectVecs,
194  rightSurf,
195  rightHit
196  );
197  }
198 
199  List<pointIndexHit> leftHit;
200  {
201  labelList leftSurf;
203  (
204  start,
205  start-projectVecs,
206  leftSurf,
207  leftHit
208  );
209  }
210 
211  // 3. Choose either -fixed, nearest, right, left.
212  forAll(displacement, i)
213  {
214  if (zonePtr && (zonePtr->localIndex(meshPoints[i]) >= 0))
215  {
216  // Fixed point. Reset to point0 location.
217  displacement[i] = points0[meshPoints[i]] - localPoints[i];
218  }
219  else if (nearest[i].hit())
220  {
221  // Found nearest.
222  displacement[i] =
223  nearest[i].hitPoint()
224  - points0[meshPoints[i]];
225  }
226  else
227  {
228  pointIndexHit interPt;
229 
230  if (rightHit[i].hit())
231  {
232  if (leftHit[i].hit())
233  {
234  if
235  (
236  magSqr(rightHit[i].hitPoint()-start[i])
237  < magSqr(leftHit[i].hitPoint()-start[i])
238  )
239  {
240  interPt = rightHit[i];
241  }
242  else
243  {
244  interPt = leftHit[i];
245  }
246  }
247  else
248  {
249  interPt = rightHit[i];
250  }
251  }
252  else
253  {
254  if (leftHit[i].hit())
255  {
256  interPt = leftHit[i];
257  }
258  }
259 
260 
261  if (interPt.hit())
262  {
263  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
264  {
265  interPt.rawPoint()[wedgePlane_] += offset[i];
266  }
267  displacement[i] = interPt.rawPoint()-points0[meshPoints[i]];
268  }
269  else
270  {
271  nNotProjected++;
272 
273  if (debug)
274  {
275  Pout<< " point:" << meshPoints[i]
276  << " coord:" << localPoints[i]
277  << " did not find any intersection between"
278  << " ray from " << start[i]-projectVecs[i]
279  << " to " << start[i]+projectVecs[i] << endl;
280  }
281  }
282  }
283  }
284  }
285 
286  reduce(nNotProjected, sumOp<label>());
287 
288  if (nNotProjected > 0)
289  {
290  Info<< "surfaceSlipDisplacement :"
291  << " on patch " << patch().name()
292  << " did not project " << nNotProjected
293  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
294  << " points." << endl;
295  }
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
300 
303 (
304  const pointPatch& p,
306  const dictionary& dict
307 )
308 :
310  surfacesDict_(dict.subDict("geometry")),
311  projectMode_(projectModeNames_.read(dict.lookup("projectMode"))),
312  projectDir_(dict.lookup<vector>("projectDirection", dimless)),
313  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
314  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
315 {}
316 
317 
320 (
322  const pointPatch& p,
324  const fieldMapper&
325 )
326 :
328  surfacesDict_(ppf.surfacesDict_),
329  projectMode_(ppf.projectMode_),
330  projectDir_(ppf.projectDir_),
331  wedgePlane_(ppf.wedgePlane_),
332  frozenPointsZone_(ppf.frozenPointsZone_)
333 {}
334 
335 
338 (
341 )
342 :
343  pointPatchVectorField(ppf, iF),
344  surfacesDict_(ppf.surfacesDict_),
345  projectMode_(ppf.projectMode_),
346  projectDir_(ppf.projectDir_),
347  wedgePlane_(ppf.wedgePlane_),
348  frozenPointsZone_(ppf.frozenPointsZone_)
349 {}
350 
351 
352 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
353 
356 {
357  if (surfacesPtr_.empty())
358  {
359  surfacesPtr_.reset
360  (
362  (
363  IOobject
364  (
365  "abc", // dummy name
366  db().time().constant(),
368  db().time(),
371  ),
372  surfacesDict_,
373  true // use single region naming shortcut
374  )
375  );
376  }
377  return surfacesPtr_();
378 }
379 
380 
382 (
383  const Pstream::commsTypes commsType
384 )
385 {
386  vectorField displacement(this->patchInternalField());
387 
388  // Calculate displacement to project points onto surface
389  calcProjection(displacement);
390 
391  // Get internal field to insert values into
392  Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
393 
394  // setInternalField(iF, motionU);
395  setInternalField(iF, displacement);
396 
398 }
399 
400 
402 {
404  writeEntry(os, "geometry", surfacesDict_);
405  writeEntry(os, "projectMode", projectModeNames_[projectMode_]);
406  writeEntry(os, "projectDirection", projectDir_);
407  writeEntry(os, "wedgePlane", wedgePlane_);
408  if (frozenPointsZone_ != word::null)
409  {
410  writeEntry(os, "frozenPointsZone", frozenPointsZone_);
411  }
412 }
413 
414 
415 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
416 
418 (
421 );
422 
423 
424 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
425 
426 } // End namespace Foam
427 
428 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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
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:105
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
const pointMesh & mesh() const
Return the mesh reference.
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.
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:407
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.
static const word & geometryDir()
Return the geometry directory name.
Displacement follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
const searchableSurfaceList & surfaces() const
Surface to follow. Demand loads surfaceNames.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Update the patch field.
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)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Namespace for OpenFOAM.
makePointPatchTypeField(pointPatchVectorField, angularOscillatingDisplacementPointPatchVectorField)
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
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
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)
Field< vector > vectorField
Specialisation of Field<T> for vector.
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void offset(label &lst, const label o)
dictionary dict
volScalarField & p
Spatial transformation functions for primitive fields.