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 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::surfaceSlipDisplacementPointPatchVectorField::projectModeNames_
43({
44 { projectMode::NEAREST, "nearest" },
45 { projectMode::POINTNORMAL, "pointNormal" },
46 { projectMode::FIXEDNORMAL, "fixedNormal" },
47});
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52void 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 = mag(mesh.bounds().max()-mesh.bounds().min());
68
69 // For case of fixed projection vector:
70 vector projectVec(0, 0, 0);
71 if (projectMode_ == FIXEDNORMAL)
72 {
73 vector n = projectDir_/mag(projectDir_);
74 projectVec = projectLen*n;
75 }
76
77
78 // Get fixed points (bit of a hack)
79 const pointZone* zonePtr = nullptr;
80
81 if (frozenPointsZone_.size() > 0)
82 {
84
85 zonePtr = &pZones[frozenPointsZone_];
86
87 Pout<< "surfaceSlipDisplacementPointPatchVectorField : Fixing all "
88 << zonePtr->size() << " points in pointZone " << zonePtr->name()
89 << endl;
90 }
91
92 // Get the starting locations from the motionSolver
93 const pointField& points0 = mesh.lookupObject<displacementMotionSolver>
94 (
95 "dynamicMeshDict"
96 ).points0();
97
98
99 pointField start(meshPoints.size());
100 forAll(start, i)
101 {
102 start[i] = points0[meshPoints[i]] + displacement[i];
103 }
104
105 label nNotProjected = 0;
106
107 if (projectMode_ == NEAREST)
108 {
109 List<pointIndexHit> nearest;
110 labelList hitSurfaces;
112 (
113 start,
114 scalarField(start.size(), sqr(projectLen)),
115 hitSurfaces,
116 nearest
117 );
118
119 forAll(nearest, i)
120 {
121 if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
122 {
123 // Fixed point. Reset to point0 location.
124 displacement[i] = points0[meshPoints[i]] - localPoints[i];
125 }
126 else if (nearest[i].hit())
127 {
128 displacement[i] =
129 nearest[i].hitPoint()
130 - points0[meshPoints[i]];
131 }
132 else
133 {
134 nNotProjected++;
135
136 if (debug)
137 {
138 Pout<< " point:" << meshPoints[i]
139 << " coord:" << localPoints[i]
140 << " did not find any surface within " << projectLen
141 << endl;
142 }
143 }
144 }
145 }
146 else
147 {
148 // Do tests on all points. Combine later on.
149
150 // 1. Check if already on surface
151 List<pointIndexHit> nearest;
152 {
153 labelList nearestSurface;
155 (
156 start,
157 scalarField(start.size(), sqr(SMALL)),
158 nearestSurface,
159 nearest
160 );
161 }
162
163 // 2. intersection. (combined later on with information from nearest
164 // above)
165 vectorField projectVecs(start.size(), projectVec);
166
167 if (projectMode_ == POINTNORMAL)
168 {
169 projectVecs = projectLen*patch().pointNormals();
170 }
171
172 // Knock out any wedge component
173 scalarField offset(start.size(), Zero);
174 if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
175 {
176 forAll(offset, i)
177 {
178 offset[i] = start[i][wedgePlane_];
179 start[i][wedgePlane_] = 0;
180 projectVecs[i][wedgePlane_] = 0;
181 }
182 }
183
184 List<pointIndexHit> rightHit;
185 {
186 labelList rightSurf;
188 (
189 start,
190 start+projectVecs,
191 rightSurf,
192 rightHit
193 );
194 }
195
196 List<pointIndexHit> leftHit;
197 {
198 labelList leftSurf;
200 (
201 start,
202 start-projectVecs,
203 leftSurf,
204 leftHit
205 );
206 }
207
208 // 3. Choose either -fixed, nearest, right, left.
209 forAll(displacement, i)
210 {
211 if (zonePtr && (zonePtr->whichPoint(meshPoints[i]) >= 0))
212 {
213 // Fixed point. Reset to point0 location.
214 displacement[i] = points0[meshPoints[i]] - localPoints[i];
215 }
216 else if (nearest[i].hit())
217 {
218 // Found nearest.
219 displacement[i] =
220 nearest[i].hitPoint()
221 - points0[meshPoints[i]];
222 }
223 else
224 {
225 pointIndexHit interPt;
226
227 if (rightHit[i].hit())
228 {
229 if (leftHit[i].hit())
230 {
231 if
232 (
233 magSqr(rightHit[i].hitPoint()-start[i])
234 < magSqr(leftHit[i].hitPoint()-start[i])
235 )
236 {
237 interPt = rightHit[i];
238 }
239 else
240 {
241 interPt = leftHit[i];
242 }
243 }
244 else
245 {
246 interPt = rightHit[i];
247 }
248 }
249 else
250 {
251 if (leftHit[i].hit())
252 {
253 interPt = leftHit[i];
254 }
255 }
256
257
258 if (interPt.hit())
259 {
260 if (wedgePlane_ >= 0 && wedgePlane_ < vector::nComponents)
261 {
262 interPt.rawPoint()[wedgePlane_] += offset[i];
263 }
264 displacement[i] = interPt.rawPoint()-points0[meshPoints[i]];
265 }
266 else
267 {
268 nNotProjected++;
269
270 if (debug)
271 {
272 Pout<< " point:" << meshPoints[i]
273 << " coord:" << localPoints[i]
274 << " did not find any intersection between"
275 << " ray from " << start[i]-projectVecs[i]
276 << " to " << start[i]+projectVecs[i] << endl;
277 }
278 }
279 }
280 }
281 }
282
283 reduce(nNotProjected, sumOp<label>());
284
285 if (nNotProjected > 0)
286 {
287 Info<< "surfaceSlipDisplacement :"
288 << " on patch " << patch().name()
289 << " did not project " << nNotProjected
290 << " out of " << returnReduce(localPoints.size(), sumOp<label>())
291 << " points." << endl;
292 }
293}
294
295
296// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
297
298Foam::surfaceSlipDisplacementPointPatchVectorField::
299surfaceSlipDisplacementPointPatchVectorField
300(
301 const pointPatch& p,
303)
304:
306 projectMode_(NEAREST),
307 projectDir_(Zero),
308 wedgePlane_(-1)
309{}
310
311
312Foam::surfaceSlipDisplacementPointPatchVectorField::
313surfaceSlipDisplacementPointPatchVectorField
314(
315 const pointPatch& p,
317 const dictionary& dict
318)
319:
321 surfacesDict_(dict.subDict("geometry")),
322 projectMode_(projectModeNames_.get("projectMode", dict)),
323 projectDir_(dict.get<vector>("projectDirection")),
324 wedgePlane_(dict.getOrDefault("wedgePlane", -1)),
325 frozenPointsZone_(dict.getOrDefault("frozenPointsZone", word::null))
326{}
327
328
329Foam::surfaceSlipDisplacementPointPatchVectorField::
330surfaceSlipDisplacementPointPatchVectorField
331(
333 const pointPatch& p,
336)
337:
339 surfacesDict_(ppf.surfacesDict_),
340 projectMode_(ppf.projectMode_),
341 projectDir_(ppf.projectDir_),
342 wedgePlane_(ppf.wedgePlane_),
343 frozenPointsZone_(ppf.frozenPointsZone_)
344{}
345
346
347Foam::surfaceSlipDisplacementPointPatchVectorField::
348surfaceSlipDisplacementPointPatchVectorField
349(
351)
352:
354 surfacesDict_(ppf.surfacesDict_),
355 projectMode_(ppf.projectMode_),
356 projectDir_(ppf.projectDir_),
357 wedgePlane_(ppf.wedgePlane_),
358 frozenPointsZone_(ppf.frozenPointsZone_)
359{}
360
361
362Foam::surfaceSlipDisplacementPointPatchVectorField::
363surfaceSlipDisplacementPointPatchVectorField
364(
367)
368:
369 pointPatchVectorField(ppf, iF),
370 surfacesDict_(ppf.surfacesDict_),
371 projectMode_(ppf.projectMode_),
372 projectDir_(ppf.projectDir_),
373 wedgePlane_(ppf.wedgePlane_),
374 frozenPointsZone_(ppf.frozenPointsZone_)
375{}
376
377
378// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
379
382{
383 if (!surfacesPtr_)
384 {
385 surfacesPtr_.reset
386 (
388 (
390 (
391 "abc", // dummy name
392 db().time().constant(), // directory
393 "triSurface", // instance
394 db().time(), // registry
397 ),
398 surfacesDict_,
399 true // use single region naming shortcut
400 )
401 );
402 }
403
404 return *surfacesPtr_;
405}
406
407
409(
410 const Pstream::commsTypes commsType
411)
412{
413 vectorField displacement(this->patchInternalField());
414
415 // Calculate displacement to project points onto surface
416 calcProjection(displacement);
417
418 // Get internal field to insert values into
419 Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
420
421 //setInInternalField(iF, motionU);
422 setInInternalField(iF, displacement);
423
425}
426
427
429(
430 Ostream& os
431) const
432{
434 os.writeEntry("geometry", surfacesDict_);
435 os.writeEntry("projectMode", projectModeNames_[projectMode_]);
436 os.writeEntry("projectDirection", projectDir_);
437 os.writeEntry("wedgePlane", wedgePlane_);
438
440 (
441 "frozenPointsZone",
443 frozenPointsZone_
444 );
445}
446
447
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450namespace Foam
451{
452
454(
457);
458
459} // End namespace Foam
460
461// ************************************************************************* //
label n
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:61
void evaluate()
Evaluate boundary conditions.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
commsTypes
Types of communications.
Definition: UPstream.H:67
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
const Type & lookupObject(const word &name, const bool recursive=false) const
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
constant condensation/saturation model.
const pointMesh & mesh() const noexcept
Return the mesh reference.
Foam::pointPatchFieldMapper.
const pointPatch & patch() const
Return patch.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
virtual const labelList & meshPoints() const =0
Return mesh points.
virtual const word & name() const =0
Return name.
virtual const vectorField & localPoints() const =0
Return pointField of points in patch.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: pointPatch.H:141
virtual const vectorField & pointNormals() const =0
Return point unit normals.
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:492
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 follows a triSurface. Use in a displacementMotionSolver as a bc on the pointDisplacement...
const searchableSurfaces & surfaces() const
Surface to follow. Demand loads surfaceNames.
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:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
dynamicFvMesh & mesh
IOporosityModelList pZones(mesh)
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Spatial transformation functions for primitive fields.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))