surfaceDisplacementPointPatchVectorField.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 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 const NamedEnum<surfaceDisplacementPointPatchVectorField::projectMode, 3>
42 surfaceDisplacementPointPatchVectorField::projectModeNames_
43 {
44  "nearest",
45  "pointNormal",
46  "fixedNormal"
47 };
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void surfaceDisplacementPointPatchVectorField::calcProjection
53 (
54  vectorField& displacement
55 ) const
56 {
57  const polyMesh& mesh = patch().boundaryMesh().mesh()();
58  const pointField& localPoints = patch().localPoints();
59  const labelList& meshPoints = patch().meshPoints();
60 
61  // const scalar deltaT = mesh.time().deltaTValue();
62 
63  // Construct large enough vector in direction of projectDir so
64  // we're guaranteed to hit something.
65 
66  //- Per point projection vector:
67  const scalar projectLen = mag(mesh.bounds().max()-mesh.bounds().min());
68 
69  // For case of fixed projection vector:
70  vector projectVec(Zero);
71  if (projectMode_ == FIXEDNORMAL)
72  {
73  vector n = projectDir_/mag(projectDir_);
74  projectVec = projectLen*n;
75  }
76 
77 
78  // Get fixed points (bit of a hack)
79  const pointZone* zonePtr = nullptr;
80 
81  if (frozenPointsZone_.size() > 0)
82  {
83  const pointZoneList& pZones = mesh.pointZones();
84 
85  zonePtr = &pZones[frozenPointsZone_];
86 
87  Pout<< "surfaceDisplacementPointPatchVectorField : Fixing all "
88  << zonePtr->size() << " points in pointZone " << zonePtr->name()
89  << endl;
90  }
91 
92  // Get the motionSolver from the mesh
93  const motionSolver& motion =
94  refCast<const fvMeshMovers::motionSolver>
95  (
96  refCast<const fvMesh>(mesh).mover()
97  ).motion();
98 
99  // Get the starting locations from the motionSolver
100  const pointField& points0 =
101  refCast<const displacementMotionSolver>(motion).points0();
102 
103  pointField start(meshPoints.size());
104  forAll(start, i)
105  {
106  start[i] = points0[meshPoints[i]] + displacement[i];
107  }
108 
109  label nNotProjected = 0;
110 
111  if (projectMode_ == NEAREST)
112  {
113  List<pointIndexHit> nearest;
114  labelList hitSurfaces;
116  (
117  start,
118  scalarField(start.size(), sqr(projectLen)),
119  hitSurfaces,
120  nearest
121  );
122 
123  forAll(nearest, i)
124  {
125  if (zonePtr && (zonePtr->localIndex(meshPoints[i]) >= 0))
126  {
127  // Fixed point. Reset to point0 location.
128  displacement[i] = points0[meshPoints[i]] - localPoints[i];
129  }
130  else if (nearest[i].hit())
131  {
132  displacement[i] =
133  nearest[i].hitPoint()
134  - points0[meshPoints[i]];
135  }
136  else
137  {
138  nNotProjected++;
139 
140  if (debug)
141  {
142  Pout<< " point:" << meshPoints[i]
143  << " coord:" << localPoints[i]
144  << " did not find any surface within " << projectLen
145  << endl;
146  }
147  }
148  }
149  }
150  else
151  {
152  // Do tests on all points. Combine later on.
153 
154  // 1. Check if already on surface
155  List<pointIndexHit> nearest;
156  {
157  labelList nearestSurface;
159  (
160  start,
161  scalarField(start.size(), sqr(small)),
162  nearestSurface,
163  nearest
164  );
165  }
166 
167  // 2. intersection. (combined later on with information from nearest
168  // above)
169  vectorField projectVecs(start.size(), projectVec);
170 
171  if (projectMode_ == POINTNORMAL)
172  {
173  projectVecs = projectLen*patch().pointNormals();
174  }
175 
176  // Knock out any wedge component
177  scalarField offset(start.size(), 0.0);
178  if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
179  {
180  forAll(offset, i)
181  {
182  offset[i] = start[i][wedgePlane_];
183  start[i][wedgePlane_] = 0;
184  projectVecs[i][wedgePlane_] = 0;
185  }
186  }
187 
188  List<pointIndexHit> rightHit;
189  {
190  labelList rightSurf;
192  (
193  start,
194  start+projectVecs,
195  rightSurf,
196  rightHit
197  );
198  }
199 
200  List<pointIndexHit> leftHit;
201  {
202  labelList leftSurf;
204  (
205  start,
206  start-projectVecs,
207  leftSurf,
208  leftHit
209  );
210  }
211 
212  // 3. Choose either -fixed, nearest, right, left.
213  forAll(displacement, i)
214  {
215  if (zonePtr && (zonePtr->localIndex(meshPoints[i]) >= 0))
216  {
217  // Fixed point. Reset to point0 location.
218  displacement[i] = points0[meshPoints[i]] - localPoints[i];
219  }
220  else if (nearest[i].hit())
221  {
222  // Found nearest.
223  displacement[i] =
224  nearest[i].hitPoint()
225  - points0[meshPoints[i]];
226  }
227  else
228  {
229  pointIndexHit interPt;
230 
231  if (rightHit[i].hit())
232  {
233  if (leftHit[i].hit())
234  {
235  if
236  (
237  magSqr(rightHit[i].hitPoint()-start[i])
238  < magSqr(leftHit[i].hitPoint()-start[i])
239  )
240  {
241  interPt = rightHit[i];
242  }
243  else
244  {
245  interPt = leftHit[i];
246  }
247  }
248  else
249  {
250  interPt = rightHit[i];
251  }
252  }
253  else
254  {
255  if (leftHit[i].hit())
256  {
257  interPt = leftHit[i];
258  }
259  }
260 
261 
262  if (interPt.hit())
263  {
264  if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
265  {
266  interPt.rawPoint()[wedgePlane_] += offset[i];
267  }
268  displacement[i] = interPt.rawPoint()-points0[meshPoints[i]];
269  }
270  else
271  {
272  nNotProjected++;
273 
274  if (debug)
275  {
276  Pout<< " point:" << meshPoints[i]
277  << " coord:" << localPoints[i]
278  << " did not find any intersection between"
279  << " ray from " << start[i]-projectVecs[i]
280  << " to " << start[i]+projectVecs[i] << endl;
281  }
282  }
283  }
284  }
285  }
286 
287  reduce(nNotProjected, sumOp<label>());
288 
289  if (nNotProjected > 0)
290  {
291  Info<< "surfaceDisplacement :"
292  << " on patch " << patch().name()
293  << " did not project " << nNotProjected
294  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
295  << " points." << endl;
296  }
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
301 
304 (
305  const pointPatch& p,
307  const dictionary& dict
308 )
309 :
310  fixedValuePointPatchVectorField(p, iF, dict),
311  velocity_(dict.lookup<vector>("velocity", dimVelocity)),
312  surfacesDict_(dict.subDict("geometry")),
313  projectMode_(projectModeNames_.read(dict.lookup("projectMode"))),
314  projectDir_(dict.lookup<vector>("projectDirection", dimless)),
315  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
316  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
317 {
318  if (velocity_.x() < 0 || velocity_.y() < 0 || velocity_.z() < 0)
319  {
321  << "All components of velocity have to be positive : "
322  << velocity_ << nl
323  << "Set velocity components to a great value if no clipping"
324  << " necessary." << exit(FatalError);
325  }
326 }
327 
328 
331 (
333  const pointPatch& p,
335  const fieldMapper& mapper
336 )
337 :
338  fixedValuePointPatchVectorField(ppf, p, iF, mapper),
339  velocity_(ppf.velocity_),
340  surfacesDict_(ppf.surfacesDict_),
341  projectMode_(ppf.projectMode_),
342  projectDir_(ppf.projectDir_),
343  wedgePlane_(ppf.wedgePlane_),
344  frozenPointsZone_(ppf.frozenPointsZone_)
345 {}
346 
347 
350 (
353 )
354 :
355  fixedValuePointPatchVectorField(ppf, iF),
356  velocity_(ppf.velocity_),
357  surfacesDict_(ppf.surfacesDict_),
358  projectMode_(ppf.projectMode_),
359  projectDir_(ppf.projectDir_),
360  wedgePlane_(ppf.wedgePlane_),
361  frozenPointsZone_(ppf.frozenPointsZone_)
362 {}
363 
364 
365 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
366 
369 {
370  if (surfacesPtr_.empty())
371  {
372  surfacesPtr_.reset
373  (
375  (
376  IOobject
377  (
378  "abc", // dummy name
379  db().time().constant(),
380  searchableSurface::geometryDir(db().time()),
381  db().time(),
384  ),
385  surfacesDict_,
386  true // allow single-region shortcut
387  )
388  );
389  }
390  return surfacesPtr_();
391 }
392 
393 
395 {
396  if (this->updated())
397  {
398  return;
399  }
400 
401  const polyMesh& mesh = patch().boundaryMesh().mesh()();
402 
403  vectorField currentDisplacement(this->patchInternalField());
404 
405  // Calculate intersections with surface w.r.t points0.
406  vectorField displacement(currentDisplacement);
407  calcProjection(displacement);
408 
409  // offset wrt current displacement
410  vectorField offset(displacement-currentDisplacement);
411 
412  // Clip offset to maximum displacement possible: velocity*timestep
413 
414  const scalar deltaT = mesh.time().deltaTValue();
415  const vector clipVelocity = velocity_*deltaT;
416 
417  forAll(displacement, i)
418  {
419  vector& d = offset[i];
420 
421  for (direction cmpt = 0; cmpt < vector::nComponents; cmpt++)
422  {
423  if (d[cmpt] < 0)
424  {
425  d[cmpt] = max(d[cmpt], -clipVelocity[cmpt]);
426  }
427  else
428  {
429  d[cmpt] = min(d[cmpt], clipVelocity[cmpt]);
430  }
431  }
432  }
433 
434  this->operator==(currentDisplacement+offset);
435 
436  fixedValuePointPatchVectorField::updateCoeffs();
437 }
438 
439 
441 {
443  writeEntry(os, "velocity", velocity_);
444  writeEntry(os, "geometry", surfacesDict_);
445  writeEntry(os, "projectMode", projectModeNames_[projectMode_]);
446  writeEntry(os, "projectDirection", projectDir_);
447  writeEntry(os, "wedgePlane", wedgePlane_);
448  if (frozenPointsZone_ != word::null)
449  {
450  writeEntry(os, "frozenPointsZone", frozenPointsZone_);
451  }
452 }
453 
454 
455 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
456 
458 (
459  fixedValuePointPatchVectorField,
461 );
462 
463 
464 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
465 
466 } // End namespace Foam
467 
468 // ************************************************************************* //
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
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:105
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
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 Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
const polyMesh & mesh() const
Return the mesh reference.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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 fixed by projection onto triSurface. Use in a displacementMotionSolver as a bc on the po...
const searchableSurfaceList & surfaces() const
Surface to follow. Demand loads surfaceNames.
surfaceDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
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.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
const dimensionSet dimVelocity
error FatalError
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void offset(label &lst, const label o)
static const char nl
Definition: Ostream.H:267
uint8_t direction
Definition: direction.H:45
dictionary dict
volScalarField & p
Spatial transformation functions for primitive fields.