Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
38const Foam::Enum
39<
41>
42Foam::surfaceDisplacementPointPatchVectorField::projectModeNames_
43({
44 { projectMode::NEAREST, "nearest" },
45 { projectMode::POINTNORMAL, "pointNormal" },
46 { projectMode::FIXEDNORMAL, "fixedNormal" },
47});
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52void Foam::surfaceDisplacementPointPatchVectorField::calcProjection
53(
54 vectorField& displacement
55) const
56{
57 const polyMesh& mesh = patch().boundaryMesh().mesh()();
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 = mesh.bounds().mag();
67
68 // For case of fixed projection vector:
69 vector projectVec(Zero);
70 if (projectMode_ == FIXEDNORMAL)
71 {
72 projectVec = projectLen * normalised(projectDir_);
73 }
74
75
76 // Get fixed points (bit of a hack)
77 const pointZone* zonePtr = nullptr;
78
79 if (frozenPointsZone_.size() > 0)
80 {
81 const pointZoneMesh& pZones = mesh.pointZones();
82
83 zonePtr = &pZones[frozenPointsZone_];
84
85 Pout<< "surfaceDisplacementPointPatchVectorField : Fixing all "
86 << zonePtr->size() << " points in pointZone " << zonePtr->name()
87 << endl;
88 }
89
90 // Get the starting locations from the motionSolver
91 const pointField& points0 = mesh.lookupObject<displacementMotionSolver>
92 (
93 "dynamicMeshDict"
94 ).points0();
95
96
97 pointField start(meshPoints.size());
98 forAll(start, i)
99 {
100 start[i] = points0[meshPoints[i]] + displacement[i];
101 }
102
103 label nNotProjected = 0;
104
105 if (projectMode_ == NEAREST)
106 {
107 List<pointIndexHit> nearest;
108 labelList hitSurfaces;
110 (
111 start,
112 scalarField(start.size(), sqr(projectLen)),
113 hitSurfaces,
114 nearest
115 );
116
117 forAll(nearest, i)
118 {
119 const label meshPointi = meshPoints[i];
120 const point& pt = mesh.points()[meshPointi];
121
122 if (zonePtr && (zonePtr->whichPoint(meshPointi) >= 0))
123 {
124 // Fixed point. Reset to point0 location.
125 displacement[i] = points0[meshPointi] - pt;
126 }
127 else if (nearest[i].hit())
128 {
129 displacement[i] =
130 nearest[i].point()
131 - points0[meshPointi];
132 }
133 else
134 {
135 nNotProjected++;
136
137 if (debug)
138 {
139 Pout<< " point:" << meshPointi
140 << " coord:" << pt
141 << " did not find any surface within " << projectLen
142 << endl;
143 }
144 }
145 }
146 }
147 else
148 {
149 // Do tests on all points. Combine later on.
150
151 // 1. Check if already on surface
152 List<pointIndexHit> nearest;
153 {
154 labelList nearestSurface;
156 (
157 start,
158 scalarField(start.size(), sqr(SMALL)),
159 nearestSurface,
160 nearest
161 );
162 }
163
164 // 2. intersection. (combined later on with information from nearest
165 // above)
166 vectorField projectVecs(start.size(), projectVec);
167
168 if (projectMode_ == POINTNORMAL)
169 {
170 projectVecs = projectLen*patch().pointNormals();
171 }
172
173 // Knock out any wedge component
174 scalarField offset(start.size(), Zero);
175 if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
176 {
177 forAll(offset, i)
178 {
179 offset[i] = start[i][wedgePlane_];
180 start[i][wedgePlane_] = 0;
181 projectVecs[i][wedgePlane_] = 0;
182 }
183 }
184
185 List<pointIndexHit> rightHit;
186 {
187 labelList rightSurf;
189 (
190 start,
191 start+projectVecs,
192 rightSurf,
193 rightHit
194 );
195 }
196
197 List<pointIndexHit> leftHit;
198 {
199 labelList leftSurf;
201 (
202 start,
203 start-projectVecs,
204 leftSurf,
205 leftHit
206 );
207 }
208
209 // 3. Choose either -fixed, nearest, right, left.
210 forAll(displacement, i)
211 {
212 const label meshPointi = meshPoints[i];
213 const point& pt = mesh.points()[meshPointi];
214
215 if (zonePtr && (zonePtr->whichPoint(meshPointi) >= 0))
216 {
217 // Fixed point. Reset to point0 location.
218 displacement[i] = points0[meshPointi] - pt;
219 }
220 else if (nearest[i].hit())
221 {
222 // Found nearest.
223 displacement[i] =
224 nearest[i].point()
225 - points0[meshPointi];
226 }
227 else
228 {
229 pointIndexHit interPt;
230
231 if (rightHit[i].hit())
232 {
233 if
234 (
235 !leftHit[i].hit()
236 ||
237 (
238 start[i].distSqr(rightHit[i].point())
239 < start[i].distSqr(leftHit[i].point())
240 )
241 )
242 {
243 interPt = rightHit[i];
244 }
245 else
246 {
247 interPt = leftHit[i];
248 }
249 }
250 else
251 {
252 if (leftHit[i].hit())
253 {
254 interPt = leftHit[i];
255 }
256 }
257
258
259 if (interPt.hit())
260 {
261 if (wedgePlane_ >= 0 && wedgePlane_ <= vector::nComponents)
262 {
263 interPt.point()[wedgePlane_] += offset[i];
264 }
265 displacement[i] = interPt.point() - points0[meshPointi];
266 }
267 else
268 {
269 nNotProjected++;
270
271 if (debug)
272 {
273 Pout<< " point:" << meshPointi
274 << " coord:" << pt
275 << " did not find any intersection between"
276 << " ray from " << start[i]-projectVecs[i]
277 << " to " << start[i]+projectVecs[i] << endl;
278 }
279 }
280 }
281 }
282 }
283
284 reduce(nNotProjected, sumOp<label>());
285
286 if (nNotProjected > 0)
287 {
288 Info<< "surfaceDisplacement :"
289 << " on patch " << patch().name()
290 << " did not project " << nNotProjected
291 << " out of " << returnReduce(meshPoints.size(), sumOp<label>())
292 << " points." << endl;
293 }
295
296
297// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
298
301(
302 const pointPatch& p,
304)
305:
306 fixedValuePointPatchVectorField(p, iF),
307 velocity_(Zero),
308 projectMode_(NEAREST),
309 projectDir_(Zero),
310 wedgePlane_(-1)
311{}
312
313
316(
317 const pointPatch& p,
319 const dictionary& dict
320)
321:
322 fixedValuePointPatchVectorField(p, iF, dict),
323 velocity_(dict.get<vector>("velocity")),
324 surfacesDict_(dict.subDict("geometry")),
325 projectMode_(projectModeNames_.get("projectMode", dict)),
326 projectDir_
327 (
328 (projectMode_ == FIXEDNORMAL)
329 ? dict.get<vector>("projectDirection")
330 : Zero
331 ),
332 wedgePlane_(dict.getOrDefault("wedgePlane", -1)),
333 frozenPointsZone_(dict.getOrDefault("frozenPointsZone", word::null))
334{
335 if (velocity_.x() < 0 || velocity_.y() < 0 || velocity_.z() < 0)
336 {
337 FatalErrorInFunction
338 << "All components of velocity have to be positive : "
339 << velocity_ << nl
340 << "Set velocity components to a great value if no clipping"
341 << " necessary." << exit(FatalError);
342 }
343}
344
345
348(
350 const pointPatch& p,
352 const pointPatchFieldMapper& mapper
353)
354:
355 fixedValuePointPatchVectorField(ppf, p, iF, mapper),
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
367(
369)
370:
371 fixedValuePointPatchVectorField(ppf),
372 velocity_(ppf.velocity_),
373 surfacesDict_(ppf.surfacesDict_),
374 projectMode_(ppf.projectMode_),
375 projectDir_(ppf.projectDir_),
376 wedgePlane_(ppf.wedgePlane_),
377 frozenPointsZone_(ppf.frozenPointsZone_)
378{}
379
380
383(
386)
387:
388 fixedValuePointPatchVectorField(ppf, iF),
389 velocity_(ppf.velocity_),
390 surfacesDict_(ppf.surfacesDict_),
391 projectMode_(ppf.projectMode_),
392 projectDir_(ppf.projectDir_),
393 wedgePlane_(ppf.wedgePlane_),
394 frozenPointsZone_(ppf.frozenPointsZone_)
395{}
397
398// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
399
402{
403 if (!surfacesPtr_)
404 {
405 surfacesPtr_.reset
406 (
408 (
410 (
411 "abc", // dummy name
412 db().time().constant(), // directory
413 "triSurface", // instance
414 db().time(), // registry
417 ),
418 surfacesDict_,
419 true // allow single-region shortcut
420 )
421 );
422 }
424 return *surfacesPtr_;
425}
426
427
429{
430 if (this->updated())
431 {
432 return;
433 }
434
435 const polyMesh& mesh = patch().boundaryMesh().mesh()();
436
437 const vectorField currentDisplacement(this->patchInternalField());
438
439 // Calculate intersections with surface w.r.t points0.
440 vectorField displacement(currentDisplacement);
441 calcProjection(displacement);
442
443
444 // offset wrt current displacement
445 vectorField offset(displacement-currentDisplacement);
446
447 // Clip offset to maximum displacement possible: velocity*timestep
448
449 const scalar deltaT = mesh.time().deltaTValue();
450 const vector clipVelocity = velocity_*deltaT;
451
452 forAll(displacement, i)
453 {
454 vector& d = offset[i];
455
456 const scalar magD(mag(d));
457 if (magD > ROOTVSMALL)
458 {
459 d /= magD;
460 d *= min(magD, mag(clipVelocity));
461 }
462 }
463
464 this->operator==(currentDisplacement+offset);
465 fixedValuePointPatchVectorField::updateCoeffs();
466}
467
468
470{
472 os.writeEntry("velocity", velocity_);
473 os.writeEntry("geometry", surfacesDict_);
474 os.writeEntry("projectMode", projectModeNames_[projectMode_]);
475 os.writeEntryIfDifferent<vector>
476 (
477 "projectDirection",
478 Zero,
479 projectDir_
480 );
481 os.writeEntryIfDifferent<label>("wedgePlane", -1, wedgePlane_);
482 os.writeEntryIfDifferent<word>
483 (
484 "frozenPointsZone",
486 frozenPointsZone_
487 );
488}
489
490
491// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
492
493namespace Foam
494{
495
497(
498 fixedValuePointPatchVectorField,
500);
501
502} // End namespace Foam
503
504// ************************************************************************* //
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...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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:134
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static constexpr direction nComponents
Number of components in this vector space.
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Foam::pointPatchFieldMapper.
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
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.
Displacement fixed by projection onto triSurface. Use in a displacementMotionSolver as a bc on the po...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
surfaceDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
virtual void write(Ostream &) const
Write.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
volScalarField & p
dynamicFvMesh & mesh
IOporosityModelList pZones(mesh)
OBJstream os(runTime.globalPath()/outputName)
Different types of constants.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))