surfaceDisplacementPointPatchVectorField.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020-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 
31 #include "Time.H"
32 #include "transformField.H"
33 #include "fvMesh.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 const Foam::Enum
39 <
41 >
42 Foam::surfaceDisplacementPointPatchVectorField::projectModeNames_
43 ({
44  { projectMode::NEAREST, "nearest" },
45  { projectMode::POINTNORMAL, "pointNormal" },
46  { projectMode::FIXEDNORMAL, "fixedNormal" },
47 });
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::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 = mesh.bounds().mag();
68 
69  // For case of fixed projection vector:
70  vector projectVec(Zero);
71  if (projectMode_ == FIXEDNORMAL)
72  {
73  projectVec = projectLen * normalised(projectDir_);
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 pointZoneMesh& pZones = mesh.pointZones();
83 
84  zonePtr = &pZones[frozenPointsZone_];
85 
86  Pout<< "surfaceDisplacementPointPatchVectorField : Fixing all "
87  << zonePtr->size() << " points in pointZone " << zonePtr->name()
88  << endl;
89  }
90 
91  // Get the starting locations from the motionSolver
92  const pointField& points0 = mesh.lookupObject<displacementMotionSolver>
93  (
94  "dynamicMeshDict"
95  ).points0();
96 
97 
98  pointField start(meshPoints.size());
99  forAll(start, i)
100  {
101  start[i] = points0[meshPoints[i]] + displacement[i];
102  }
103 
104  label nNotProjected = 0;
105 
106  if (projectMode_ == NEAREST)
107  {
108  List<pointIndexHit> nearest;
109  labelList hitSurfaces;
111  (
112  start,
113  scalarField(start.size(), sqr(projectLen)),
114  hitSurfaces,
115  nearest
116  );
117 
118  forAll(nearest, i)
119  {
120  if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
121  {
122  // Fixed point. Reset to point0 location.
123  displacement[i] = points0[meshPoints[i]] - localPoints[i];
124  }
125  else if (nearest[i].hit())
126  {
127  displacement[i] =
128  nearest[i].point()
129  - points0[meshPoints[i]];
130  }
131  else
132  {
133  nNotProjected++;
134 
135  if (debug)
136  {
137  Pout<< " point:" << meshPoints[i]
138  << " coord:" << localPoints[i]
139  << " did not find any surface within " << projectLen
140  << endl;
141  }
142  }
143  }
144  }
145  else
146  {
147  // Do tests on all points. Combine later on.
148 
149  // 1. Check if already on surface
150  List<pointIndexHit> nearest;
151  {
152  labelList nearestSurface;
154  (
155  start,
156  scalarField(start.size(), sqr(SMALL)),
157  nearestSurface,
158  nearest
159  );
160  }
161 
162  // 2. intersection. (combined later on with information from nearest
163  // above)
164  vectorField projectVecs(start.size(), projectVec);
165 
166  if (projectMode_ == POINTNORMAL)
167  {
168  projectVecs = projectLen*patch().pointNormals();
169  }
170 
171  // Knock out any wedge component
172  scalarField offset(start.size(), Zero);
173  if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
174  {
175  forAll(offset, i)
176  {
177  offset[i] = start[i][wedgePlane_];
178  start[i][wedgePlane_] = 0;
179  projectVecs[i][wedgePlane_] = 0;
180  }
181  }
182 
183  List<pointIndexHit> rightHit;
184  {
185  labelList rightSurf;
187  (
188  start,
189  start+projectVecs,
190  rightSurf,
191  rightHit
192  );
193  }
194 
195  List<pointIndexHit> leftHit;
196  {
197  labelList leftSurf;
199  (
200  start,
201  start-projectVecs,
202  leftSurf,
203  leftHit
204  );
205  }
206 
207  // 3. Choose either -fixed, nearest, right, left.
208  forAll(displacement, i)
209  {
210  if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
211  {
212  // Fixed point. Reset to point0 location.
213  displacement[i] = points0[meshPoints[i]] - localPoints[i];
214  }
215  else if (nearest[i].hit())
216  {
217  // Found nearest.
218  displacement[i] =
219  nearest[i].point()
220  - points0[meshPoints[i]];
221  }
222  else
223  {
224  pointIndexHit interPt;
225 
226  if (rightHit[i].hit())
227  {
228  if
229  (
230  !leftHit[i].hit()
231  ||
232  (
233  start[i].distSqr(rightHit[i].point())
234  < start[i].distSqr(leftHit[i].point())
235  )
236  )
237  {
238  interPt = rightHit[i];
239  }
240  else
241  {
242  interPt = leftHit[i];
243  }
244  }
245  else
246  {
247  if (leftHit[i].hit())
248  {
249  interPt = leftHit[i];
250  }
251  }
252 
253 
254  if (interPt.hit())
255  {
256  if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
257  {
258  interPt.point()[wedgePlane_] += offset[i];
259  }
260  displacement[i] = interPt.point() - points0[meshPoints[i]];
261  }
262  else
263  {
264  nNotProjected++;
265 
266  if (debug)
267  {
268  Pout<< " point:" << meshPoints[i]
269  << " coord:" << localPoints[i]
270  << " did not find any intersection between"
271  << " ray from " << start[i]-projectVecs[i]
272  << " to " << start[i]+projectVecs[i] << endl;
273  }
274  }
275  }
276  }
277  }
278 
279  reduce(nNotProjected, sumOp<label>());
280 
281  if (nNotProjected > 0)
282  {
283  Info<< "surfaceDisplacement :"
284  << " on patch " << patch().name()
285  << " did not project " << nNotProjected
286  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
287  << " points." << endl;
288  }
289 }
290 
292 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
293 
296 (
297  const pointPatch& p,
299 )
300 :
301  fixedValuePointPatchVectorField(p, iF),
302  velocity_(Zero),
303  projectMode_(NEAREST),
304  projectDir_(Zero),
305  wedgePlane_(-1)
306 {}
307 
308 
311 (
312  const pointPatch& p,
314  const dictionary& dict
315 )
316 :
317  fixedValuePointPatchVectorField(p, iF, dict),
318  velocity_(dict.get<vector>("velocity")),
319  surfacesDict_(dict.subDict("geometry")),
320  projectMode_(projectModeNames_.get("projectMode", dict)),
321  projectDir_(dict.get<vector>("projectDirection")),
322  wedgePlane_(dict.getOrDefault("wedgePlane", -1)),
323  frozenPointsZone_(dict.getOrDefault("frozenPointsZone", word::null))
324 {
325  if (velocity_.x() < 0 || velocity_.y() < 0 || velocity_.z() < 0)
326  {
328  << "All components of velocity have to be positive : "
329  << velocity_ << nl
330  << "Set velocity components to a great value if no clipping"
331  << " necessary." << exit(FatalError);
332  }
333 }
334 
335 
338 (
340  const pointPatch& p,
342  const pointPatchFieldMapper& mapper
343 )
344 :
345  fixedValuePointPatchVectorField(ppf, p, iF, mapper),
346  velocity_(ppf.velocity_),
347  surfacesDict_(ppf.surfacesDict_),
348  projectMode_(ppf.projectMode_),
349  projectDir_(ppf.projectDir_),
350  wedgePlane_(ppf.wedgePlane_),
351  frozenPointsZone_(ppf.frozenPointsZone_)
352 {}
353 
354 
357 (
359 )
360 :
361  fixedValuePointPatchVectorField(ppf),
362  velocity_(ppf.velocity_),
363  surfacesDict_(ppf.surfacesDict_),
364  projectMode_(ppf.projectMode_),
365  projectDir_(ppf.projectDir_),
366  wedgePlane_(ppf.wedgePlane_),
367  frozenPointsZone_(ppf.frozenPointsZone_)
368 {}
369 
370 
373 (
376 )
377 :
378  fixedValuePointPatchVectorField(ppf, iF),
379  velocity_(ppf.velocity_),
380  surfacesDict_(ppf.surfacesDict_),
381  projectMode_(ppf.projectMode_),
382  projectDir_(ppf.projectDir_),
383  wedgePlane_(ppf.wedgePlane_),
384  frozenPointsZone_(ppf.frozenPointsZone_)
385 {}
387 
388 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
389 
392 {
393  if (!surfacesPtr_)
394  {
395  surfacesPtr_.reset
396  (
398  (
399  IOobject
400  (
401  "abc", // dummy name
402  db().time().constant(), // directory
403  "triSurface", // instance
404  db().time(), // registry
407  ),
408  surfacesDict_,
409  true // allow single-region shortcut
410  )
411  );
412  }
414  return *surfacesPtr_;
415 }
416 
417 
419 {
420  if (this->updated())
421  {
422  return;
423  }
424 
425  const polyMesh& mesh = patch().boundaryMesh().mesh()();
426 
427  vectorField currentDisplacement(this->patchInternalField());
428 
429  // Calculate intersections with surface w.r.t points0.
430  vectorField displacement(currentDisplacement);
431  calcProjection(displacement);
432 
433  // offset wrt current displacement
434  vectorField offset(displacement-currentDisplacement);
435 
436  // Clip offset to maximum displacement possible: velocity*timestep
437 
438  const scalar deltaT = mesh.time().deltaTValue();
439  const vector clipVelocity = velocity_*deltaT;
440 
441  forAll(displacement, i)
442  {
443  vector& d = offset[i];
444 
445  for (direction cmpt = 0; cmpt < vector::nComponents; cmpt++)
446  {
447  if (d[cmpt] < 0)
448  {
449  d[cmpt] = max(d[cmpt], -clipVelocity[cmpt]);
450  }
451  else
452  {
453  d[cmpt] = min(d[cmpt], clipVelocity[cmpt]);
454  }
455  }
456  }
457 
458  this->operator==(currentDisplacement+offset);
460  fixedValuePointPatchVectorField::updateCoeffs();
461 }
462 
463 
465 {
467  os.writeEntry("velocity", velocity_);
468  os.writeEntry("geometry", surfacesDict_);
469  os.writeEntry("projectMode", projectModeNames_[projectMode_]);
470  os.writeEntry("projectDirection", projectDir_);
471  os.writeEntry("wedgePlane", wedgePlane_);
472 
474  (
475  "frozenPointsZone",
476  word::null,
477  frozenPointsZone_
478  );
479 }
480 
481 
482 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
483 
484 namespace Foam
485 {
486 
488 (
489  fixedValuePointPatchVectorField,
490  surfaceDisplacementPointPatchVectorField
491 );
492 
493 } // End namespace Foam
494 
495 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
uint8_t direction
Definition: direction.H:46
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Foam::pointPatchFieldMapper.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
surfaceDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
virtual void write(Ostream &) const
Write.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Spatial transformation functions for primitive fields.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const word null
An empty word.
Definition: word.H:84
Container for searchableSurfaces. The collection is specified as a dictionary. For example...
Vector< scalar > vector
Definition: vector.H:57
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
Displacement fixed by projection onto triSurface. Use in a displacementMotionSolver as a bc on the po...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
vector point
Point is a vector.
Definition: point.H:37
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
IOporosityModelList pZones(mesh)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127