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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "Time.H"
31 #include "transformField.H"
32 #include "fvMesh.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 const Foam::Enum
38 <
40 >
41 Foam::surfaceDisplacementPointPatchVectorField::projectModeNames_
42 ({
43  { projectMode::NEAREST, "nearest" },
44  { projectMode::POINTNORMAL, "pointNormal" },
45  { projectMode::FIXEDNORMAL, "fixedNormal" },
46 });
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::surfaceDisplacementPointPatchVectorField::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(Zero);
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 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].hitPoint()
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].hitPoint()
220  - points0[meshPoints[i]];
221  }
222  else
223  {
224  pointIndexHit interPt;
225 
226  if (rightHit[i].hit())
227  {
228  if (leftHit[i].hit())
229  {
230  if
231  (
232  magSqr(rightHit[i].hitPoint()-start[i])
233  < magSqr(leftHit[i].hitPoint()-start[i])
234  )
235  {
236  interPt = rightHit[i];
237  }
238  else
239  {
240  interPt = leftHit[i];
241  }
242  }
243  else
244  {
245  interPt = rightHit[i];
246  }
247  }
248  else
249  {
250  if (leftHit[i].hit())
251  {
252  interPt = leftHit[i];
253  }
254  }
255 
256 
257  if (interPt.hit())
258  {
259  if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
260  {
261  interPt.rawPoint()[wedgePlane_] += offset[i];
262  }
263  displacement[i] = interPt.rawPoint()-points0[meshPoints[i]];
264  }
265  else
266  {
267  nNotProjected++;
268 
269  if (debug)
270  {
271  Pout<< " point:" << meshPoints[i]
272  << " coord:" << localPoints[i]
273  << " did not find any intersection between"
274  << " ray from " << start[i]-projectVecs[i]
275  << " to " << start[i]+projectVecs[i] << endl;
276  }
277  }
278  }
279  }
280  }
281 
282  reduce(nNotProjected, sumOp<label>());
283 
284  if (nNotProjected > 0)
285  {
286  Info<< "surfaceDisplacement :"
287  << " on patch " << patch().name()
288  << " did not project " << nNotProjected
289  << " out of " << returnReduce(localPoints.size(), sumOp<label>())
290  << " points." << endl;
291  }
292 }
293 
294 
295 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
296 
299 (
300  const pointPatch& p,
302 )
303 :
304  fixedValuePointPatchVectorField(p, iF),
305  velocity_(Zero),
306  projectMode_(NEAREST),
307  projectDir_(Zero),
308  wedgePlane_(-1)
309 {}
310 
311 
314 (
315  const pointPatch& p,
317  const dictionary& dict
318 )
319 :
320  fixedValuePointPatchVectorField(p, iF, dict),
321  velocity_(dict.get<vector>("velocity")),
322  surfacesDict_(dict.subDict("geometry")),
323  projectMode_(projectModeNames_.get("projectMode", dict)),
324  projectDir_(dict.get<vector>("projectDirection")),
325  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
326  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
327 {
328  if (velocity_.x() < 0 || velocity_.y() < 0 || velocity_.z() < 0)
329  {
331  << "All components of velocity have to be positive : "
332  << velocity_ << nl
333  << "Set velocity components to a great value if no clipping"
334  << " necessary." << exit(FatalError);
335  }
336 }
337 
338 
341 (
343  const pointPatch& p,
345  const pointPatchFieldMapper& mapper
346 )
347 :
348  fixedValuePointPatchVectorField(ppf, p, iF, mapper),
349  velocity_(ppf.velocity_),
350  surfacesDict_(ppf.surfacesDict_),
351  projectMode_(ppf.projectMode_),
352  projectDir_(ppf.projectDir_),
353  wedgePlane_(ppf.wedgePlane_),
354  frozenPointsZone_(ppf.frozenPointsZone_)
355 {}
356 
357 
360 (
362 )
363 :
364  fixedValuePointPatchVectorField(ppf),
365  velocity_(ppf.velocity_),
366  surfacesDict_(ppf.surfacesDict_),
367  projectMode_(ppf.projectMode_),
368  projectDir_(ppf.projectDir_),
369  wedgePlane_(ppf.wedgePlane_),
370  frozenPointsZone_(ppf.frozenPointsZone_)
371 {}
372 
373 
376 (
379 )
380 :
381  fixedValuePointPatchVectorField(ppf, iF),
382  velocity_(ppf.velocity_),
383  surfacesDict_(ppf.surfacesDict_),
384  projectMode_(ppf.projectMode_),
385  projectDir_(ppf.projectDir_),
386  wedgePlane_(ppf.wedgePlane_),
387  frozenPointsZone_(ppf.frozenPointsZone_)
388 {}
389 
390 
391 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
392 
395 {
396  if (surfacesPtr_.empty())
397  {
398  surfacesPtr_.reset
399  (
401  (
402  IOobject
403  (
404  "abc", // dummy name
405  db().time().constant(), // directory
406  "triSurface", // instance
407  db().time(), // registry
410  ),
411  surfacesDict_,
412  true // allow single-region shortcut
413  )
414  );
415  }
416 
417  return *surfacesPtr_;
418 }
419 
420 
422 {
423  if (this->updated())
424  {
425  return;
426  }
427 
428  const polyMesh& mesh = patch().boundaryMesh().mesh()();
429 
430  vectorField currentDisplacement(this->patchInternalField());
431 
432  // Calculate intersections with surface w.r.t points0.
433  vectorField displacement(currentDisplacement);
434  calcProjection(displacement);
435 
436  // offset wrt current displacement
437  vectorField offset(displacement-currentDisplacement);
438 
439  // Clip offset to maximum displacement possible: velocity*timestep
440 
441  const scalar deltaT = mesh.time().deltaTValue();
442  const vector clipVelocity = velocity_*deltaT;
443 
444  forAll(displacement, i)
445  {
446  vector& d = offset[i];
447 
448  for (direction cmpt = 0; cmpt < vector::nComponents; cmpt++)
449  {
450  if (d[cmpt] < 0)
451  {
452  d[cmpt] = max(d[cmpt], -clipVelocity[cmpt]);
453  }
454  else
455  {
456  d[cmpt] = min(d[cmpt], clipVelocity[cmpt]);
457  }
458  }
459  }
460 
461  this->operator==(currentDisplacement+offset);
462 
463  fixedValuePointPatchVectorField::updateCoeffs();
464 }
465 
466 
468 {
470  os.writeEntry("velocity", velocity_);
471  os.writeEntry("geometry", surfacesDict_);
472  os.writeEntry("projectMode", projectModeNames_[projectMode_]);
473  os.writeEntry("projectDirection", projectDir_);
474  os.writeEntry("wedgePlane", wedgePlane_);
475 
477  (
478  "frozenPointsZone",
479  word::null,
480  frozenPointsZone_
481  );
482 }
483 
484 
485 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486 
487 namespace Foam
488 {
489 
491 (
492  fixedValuePointPatchVectorField,
494 );
495 
496 } // End namespace Foam
497 
498 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::surfaceDisplacementPointPatchVectorField::surfaceDisplacementPointPatchVectorField
surfaceDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Definition: surfaceDisplacementPointPatchVectorField.C:299
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::Ostream::writeEntryIfDifferent
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:231
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::surfaceDisplacementPointPatchVectorField::projectMode
projectMode
Definition: surfaceDisplacementPointPatchVectorField.H:85
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::surfaceDisplacementPointPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: surfaceDisplacementPointPatchVectorField.C:421
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:58
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::pointZoneMesh
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
Definition: pointZoneMeshFwd.H:44
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
transformField.H
Spatial transformation functions for primitive fields.
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:48
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
pZones
IOporosityModelList pZones(mesh)
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::searchableSurfaces::findAnyIntersection
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
Definition: searchableSurfaces.C:331
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::makePointPatchTypeField
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
Foam::label
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:62
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::searchableSurfaces::findNearest
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
Definition: searchableSurfaces.C:396
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:91
surfaceDisplacementPointPatchVectorField.H
Foam::surfaceDisplacementPointPatchVectorField::NEAREST
Definition: surfaceDisplacementPointPatchVectorField.H:87
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:523
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
displacementMotionSolver.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::surfaceDisplacementPointPatchVectorField::FIXEDNORMAL
Definition: surfaceDisplacementPointPatchVectorField.H:89
Time.H
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::TimeState::deltaTValue
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:42
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:92
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::surfaceDisplacementPointPatchVectorField
Displacement fixed by projection onto triSurface. Use in a displacementMotionSolver as a bc on the po...
Definition: surfaceDisplacementPointPatchVectorField.H:76
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::surfaceDisplacementPointPatchVectorField::POINTNORMAL
Definition: surfaceDisplacementPointPatchVectorField.H:88
constant
constant condensation/saturation model.
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::surfaceDisplacementPointPatchVectorField::surfaces
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
Definition: surfaceDisplacementPointPatchVectorField.C:394
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::surfaceDisplacementPointPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: surfaceDisplacementPointPatchVectorField.C:467