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 -------------------------------------------------------------------------------
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::surfaceSlipDisplacementPointPatchVectorField::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::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 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] =
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<< "surfaceSlipDisplacement :"
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 :
305  projectMode_(NEAREST),
306  projectDir_(Zero),
307  wedgePlane_(-1)
308 {}
309 
310 
313 (
314  const pointPatch& p,
316  const dictionary& dict
317 )
318 :
320  surfacesDict_(dict.subDict("geometry")),
321  projectMode_(projectModeNames_.get("projectMode", dict)),
322  projectDir_(dict.get<vector>("projectDirection")),
323  wedgePlane_(dict.lookupOrDefault("wedgePlane", -1)),
324  frozenPointsZone_(dict.lookupOrDefault("frozenPointsZone", word::null))
325 {}
326 
327 
330 (
332  const pointPatch& p,
334  const pointPatchFieldMapper&
335 )
336 :
338  surfacesDict_(ppf.surfacesDict_),
339  projectMode_(ppf.projectMode_),
340  projectDir_(ppf.projectDir_),
341  wedgePlane_(ppf.wedgePlane_),
342  frozenPointsZone_(ppf.frozenPointsZone_)
343 {}
344 
345 
348 (
350 )
351 :
353  surfacesDict_(ppf.surfacesDict_),
354  projectMode_(ppf.projectMode_),
355  projectDir_(ppf.projectDir_),
356  wedgePlane_(ppf.wedgePlane_),
357  frozenPointsZone_(ppf.frozenPointsZone_)
358 {}
359 
360 
363 (
366 )
367 :
368  pointPatchVectorField(ppf, iF),
369  surfacesDict_(ppf.surfacesDict_),
370  projectMode_(ppf.projectMode_),
371  projectDir_(ppf.projectDir_),
372  wedgePlane_(ppf.wedgePlane_),
373  frozenPointsZone_(ppf.frozenPointsZone_)
374 {}
375 
376 
377 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
378 
381 {
382  if (surfacesPtr_.empty())
383  {
384  surfacesPtr_.reset
385  (
387  (
388  IOobject
389  (
390  "abc", // dummy name
391  db().time().constant(), // directory
392  "triSurface", // instance
393  db().time(), // registry
396  ),
397  surfacesDict_,
398  true // use single region naming shortcut
399  )
400  );
401  }
402 
403  return *surfacesPtr_;
404 }
405 
406 
408 (
409  const Pstream::commsTypes commsType
410 )
411 {
412  vectorField displacement(this->patchInternalField());
413 
414  // Calculate displacement to project points onto surface
415  calcProjection(displacement);
416 
417  // Get internal field to insert values into
418  Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
419 
420  //setInInternalField(iF, motionU);
421  setInInternalField(iF, displacement);
422 
424 }
425 
426 
428 (
429  Ostream& os
430 ) const
431 {
433  os.writeEntry("geometry", surfacesDict_);
434  os.writeEntry("projectMode", projectModeNames_[projectMode_]);
435  os.writeEntry("projectDirection", projectDir_);
436  os.writeEntry("wedgePlane", wedgePlane_);
437 
439  (
440  "frozenPointsZone",
441  word::null,
442  frozenPointsZone_
443  );
444 }
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 namespace Foam
450 {
451 
453 (
456 );
457 
458 } // End namespace Foam
459 
460 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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
Foam::surfaceSlipDisplacementPointPatchVectorField::POINTNORMAL
Definition: surfaceSlipDisplacementPointPatchVectorField.H:81
Foam::surfaceSlipDisplacementPointPatchVectorField::projectMode
projectMode
Definition: surfaceSlipDisplacementPointPatchVectorField.H:78
p
volScalarField & p
Definition: createFieldRefs.H:8
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::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::pointPatch::pointNormals
virtual const vectorField & pointNormals() const =0
Return point normals.
Foam::pointPatch::localPoints
virtual const vectorField & localPoints() const =0
Return mesh points.
Foam::pointPatchField< vector >::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
Definition: pointPatchField.C:294
Foam::surfaceSlipDisplacementPointPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: surfaceSlipDisplacementPointPatchVectorField.C:428
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::pointPatchField< vector >
transformField.H
Spatial transformation functions for primitive fields.
Foam::pointPatchField< vector >::patch
const pointPatch & patch() const
Return patch.
Definition: pointPatchField.H:268
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)
Foam::surfaceSlipDisplacementPointPatchVectorField::surfaces
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
Definition: surfaceSlipDisplacementPointPatchVectorField.C:380
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::surfaceSlipDisplacementPointPatchVectorField::surfaceSlipDisplacementPointPatchVectorField
surfaceSlipDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Definition: surfaceSlipDisplacementPointPatchVectorField.C:299
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
Foam::pointPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: pointPatchField.C:116
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::pointBoundaryMesh::mesh
const pointMesh & mesh() const
Return the mesh reference.
Definition: pointBoundaryMesh.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
surfaceSlipDisplacementPointPatchVectorField.H
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::pointPatch::name
virtual const word & name() const =0
Return name.
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
displacementMotionSolver.H
Time.H
Foam::surfaceSlipDisplacementPointPatchVectorField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Update the patch field.
Definition: surfaceSlipDisplacementPointPatchVectorField.C:408
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::Vector< scalar >
Foam::pointPatch::meshPoints
virtual const labelList & meshPoints() const =0
Return mesh points.
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::pointPatch::boundaryMesh
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: pointPatch.H:135
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::pointPatchVectorField
pointPatchField< vector > pointPatchVectorField
Definition: pointPatchFieldsFwd.H:43
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
constant
constant condensation/saturation model.
Foam::surfaceSlipDisplacementPointPatchVectorField::FIXEDNORMAL
Definition: surfaceSlipDisplacementPointPatchVectorField.H:82
Foam::surfaceSlipDisplacementPointPatchVectorField::NEAREST
Definition: surfaceSlipDisplacementPointPatchVectorField.H:80
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::surfaceSlipDisplacementPointPatchVectorField
Displacement follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
Definition: surfaceSlipDisplacementPointPatchVectorField.H:69
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