surfaceSlipDisplacementPointPatchVectorField.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::surfaceSlipDisplacementPointPatchVectorField::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::surfaceSlipDisplacementPointPatchVectorField::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<< "surfaceSlipDisplacementPointPatchVectorField : 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] = nearest[i].point() - points0[meshPoints[i]];
128  }
129  else
130  {
131  nNotProjected++;
132 
133  if (debug)
134  {
135  Pout<< " point:" << meshPoints[i]
136  << " coord:" << localPoints[i]
137  << " did not find any surface within " << projectLen
138  << endl;
139  }
140  }
141  }
142  }
143  else
144  {
145  // Do tests on all points. Combine later on.
146 
147  // 1. Check if already on surface
148  List<pointIndexHit> nearest;
149  {
150  labelList nearestSurface;
152  (
153  start,
154  scalarField(start.size(), sqr(SMALL)),
155  nearestSurface,
156  nearest
157  );
158  }
159 
160  // 2. intersection. (combined later on with information from nearest
161  // above)
162  vectorField projectVecs(start.size(), projectVec);
163 
164  if (projectMode_ == POINTNORMAL)
165  {
166  projectVecs = projectLen*patch().pointNormals();
167  }
168 
169  // Knock out any wedge component
170  scalarField offset(start.size(), Zero);
171  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
172  {
173  forAll(offset, i)
174  {
175  offset[i] = start[i][wedgePlane_];
176  start[i][wedgePlane_] = 0;
177  projectVecs[i][wedgePlane_] = 0;
178  }
179  }
180 
181  List<pointIndexHit> rightHit;
182  {
183  labelList rightSurf;
185  (
186  start,
187  start+projectVecs,
188  rightSurf,
189  rightHit
190  );
191  }
192 
193  List<pointIndexHit> leftHit;
194  {
195  labelList leftSurf;
197  (
198  start,
199  start-projectVecs,
200  leftSurf,
201  leftHit
202  );
203  }
204 
205  // 3. Choose either -fixed, nearest, right, left.
206  forAll(displacement, i)
207  {
208  if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
209  {
210  // Fixed point. Reset to point0 location.
211  displacement[i] = points0[meshPoints[i]] - localPoints[i];
212  }
213  else if (nearest[i].hit())
214  {
215  // Found nearest.
216  displacement[i] = nearest[i].point() - points0[meshPoints[i]];
217  }
218  else
219  {
220  pointIndexHit interPt;
221 
222  if (rightHit[i].hit())
223  {
224  if
225  (
226  !leftHit[i].hit()
227  ||
228  (
229  start[i].distSqr(rightHit[i].point())
230  < start[i].distSqr(leftHit[i].point())
231  )
232  )
233  {
234  interPt = rightHit[i];
235  }
236  else
237  {
238  interPt = leftHit[i];
239  }
240  }
241  else
242  {
243  if (leftHit[i].hit())
244  {
245  interPt = leftHit[i];
246  }
247  }
248 
249 
250  if (interPt.hit())
251  {
252  if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
253  {
254  interPt.point()[wedgePlane_] += offset[i];
255  }
256  displacement[i] = interPt.point() - points0[meshPoints[i]];
257  }
258  else
259  {
260  nNotProjected++;
261 
262  if (debug)
263  {
264  Pout<< " point:" << meshPoints[i]
265  << " coord:" << localPoints[i]
266  << " did not find any intersection between"
267  << " ray from " << start[i]-projectVecs[i]
268  << " to " << start[i]+projectVecs[i] << endl;
269  }
270  }
271  }
272  }
273  }
274 
275  reduce(nNotProjected, sumOp<label>());
276 
277  if (nNotProjected > 0)
278  {
279  Info<< "surfaceSlipDisplacement :"
280  << " on patch " << patch().name()
281  << " did not project " << nNotProjected
282  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
283  << " points." << endl;
284  }
285 }
286 
288 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
289 
292 (
293  const pointPatch& p,
295 )
296 :
298  projectMode_(NEAREST),
299  projectDir_(Zero),
300  wedgePlane_(-1)
301 {}
302 
303 
306 (
307  const pointPatch& p,
309  const dictionary& dict
310 )
311 :
313  surfacesDict_(dict.subDict("geometry")),
314  projectMode_(projectModeNames_.get("projectMode", dict)),
315  projectDir_(dict.get<vector>("projectDirection")),
316  wedgePlane_(dict.getOrDefault("wedgePlane", -1)),
317  frozenPointsZone_(dict.getOrDefault("frozenPointsZone", word::null))
318 {}
319 
320 
323 (
325  const pointPatch& p,
327  const pointPatchFieldMapper&
328 )
329 :
331  surfacesDict_(ppf.surfacesDict_),
332  projectMode_(ppf.projectMode_),
333  projectDir_(ppf.projectDir_),
334  wedgePlane_(ppf.wedgePlane_),
335  frozenPointsZone_(ppf.frozenPointsZone_)
336 {}
337 
338 
341 (
343 )
344 :
346  surfacesDict_(ppf.surfacesDict_),
347  projectMode_(ppf.projectMode_),
348  projectDir_(ppf.projectDir_),
349  wedgePlane_(ppf.wedgePlane_),
350  frozenPointsZone_(ppf.frozenPointsZone_)
351 {}
352 
353 
356 (
359 )
360 :
361  pointPatchVectorField(ppf, iF),
362  surfacesDict_(ppf.surfacesDict_),
363  projectMode_(ppf.projectMode_),
364  projectDir_(ppf.projectDir_),
365  wedgePlane_(ppf.wedgePlane_),
366  frozenPointsZone_(ppf.frozenPointsZone_)
367 {}
369 
370 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
371 
374 {
375  if (!surfacesPtr_)
376  {
377  surfacesPtr_.reset
378  (
380  (
381  IOobject
382  (
383  "abc", // dummy name
384  db().time().constant(), // directory
385  "triSurface", // instance
386  db().time(), // registry
389  ),
390  surfacesDict_,
391  true // use single region naming shortcut
392  )
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  //setInInternalField(iF, motionU);
414  setInInternalField(iF, displacement);
415 
417 }
418 
419 
421 (
422  Ostream& os
423 ) const
424 {
426  os.writeEntry("geometry", surfacesDict_);
427  os.writeEntry("projectMode", projectModeNames_[projectMode_]);
428  os.writeEntry("projectDirection", projectDir_);
429  os.writeEntry("wedgePlane", wedgePlane_);
430 
432  (
433  "frozenPointsZone",
434  word::null,
435  frozenPointsZone_
436  );
437 }
438 
439 
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441 
442 namespace Foam
443 {
444 
446 (
448  surfaceSlipDisplacementPointPatchVectorField
449 );
450 
451 } // End namespace Foam
452 
453 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
dictionary dict
Displacement follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
const pointPatch & patch() const noexcept
Return the patch.
commsTypes
Communications types.
Definition: UPstream.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Update the patch field.
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)
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
virtual void write(Ostream &os) const
Write.
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
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.
const pointMesh & mesh() const noexcept
Return the mesh reference.
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...
virtual const word & name() const =0
Return name.
Vector< scalar > vector
Definition: vector.H:57
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
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
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)
pointPatchField< vector > pointPatchVectorField
virtual const vectorField & pointNormals() const =0
Return point unit normals.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: pointPatch.H:171
virtual const vectorField & localPoints() const =0
Return pointField of points in patch.
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
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
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
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
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.
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)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127