Loading...
Searching...
No Matches
pointAttractionDisplacementPointPatchVectorField.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) 2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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"
33#include "featureEdgeMesh.H"
35
36// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37
38void Foam::pointAttractionDisplacementPointPatchVectorField::calcProjection
39(
40 vectorField& displacement
41) const
42{
43 const polyMesh& mesh = patch().boundaryMesh().mesh()();
44 const labelList& meshPoints = patch().meshPoints();
45
46 //const scalar deltaT = mesh.time().deltaTValue();
47
48 // Construct large enough vector in direction of projectDir so
49 // we're guaranteed to hit something.
50
51 //- Per point projection vector:
52 const scalar projectLen = mesh.bounds().mag();
53
54
55
56 // Get fixed points (bit of a hack)
57 const pointZone* zonePtr = nullptr;
58
59 if (frozenPointsZone_.size() > 0)
60 {
61 const pointZoneMesh& pZones = mesh.pointZones();
62
63 zonePtr = &pZones[frozenPointsZone_];
64
65 Pout<< "pointAttractionDisplacementPointPatchVectorField : Fixing all "
66 << zonePtr->size() << " points in pointZone " << zonePtr->name()
67 << endl;
68 }
69
70 // Get the starting locations from the motionSolver
71 const pointField& points0 = mesh.lookupObject<displacementMotionSolver>
72 (
73 "dynamicMeshDict"
74 ).points0();
75
76
77 pointField start(meshPoints.size());
78 forAll(start, i)
79 {
80 start[i] = points0[meshPoints[i]] + displacement[i];
81 }
82
83 const auto& tree = pointTree();
84
85 label nNotProjected = 0;
86 forAll(meshPoints, i)
87 {
88 const label meshPointi = meshPoints[i];
89 const point& pt = mesh.points()[meshPointi];
90
91 if (zonePtr && (zonePtr->whichPoint(meshPointi) >= 0))
92 {
93 // Fixed point. Reset to point0 location.
94 displacement[i] = points0[meshPointi] - pt;
95 }
96 else
97 {
98 pointIndexHit nearest = tree.findNearest(start[i], sqr(projectLen));
99 if (nearest.hit())
100 {
101 displacement[i] = nearest.point() - points0[meshPointi];
102 }
103 else
104 {
105 nNotProjected++;
106
107 if (debug)
108 {
109 Pout<< " point:" << meshPointi
110 << " coord:" << pt
111 << " did not find any surface within " << projectLen
112 << endl;
113 }
114 }
115 }
116 }
117
118 reduce(nNotProjected, sumOp<label>());
119
120 if (nNotProjected > 0)
121 {
122 Info<< "pointAttractionDisplacement :"
123 << " on patch " << patch().name()
124 << " did not project " << nNotProjected
125 << " out of " << returnReduce(meshPoints.size(), sumOp<label>())
126 << " points." << endl;
127 }
129
130
131// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132
135(
136 const pointPatch& p,
138)
139:
141 velocity_(Zero)
142{}
143
144
147(
148 const pointPatch& p,
150 const dictionary& dict
151)
152:
154 velocity_(dict.get<vector>("velocity")),
155 featFileName_(dict.get<fileName>("file", keyType::LITERAL)),
156 frozenPointsZone_(dict.getOrDefault("frozenPointsZone", word::null))
157{
158 // Read&store edge mesh on registry
160 (
161 this->patch().boundaryMesh().mesh().time(),
177 velocity_(ppf.velocity_),
178 featFileName_(ppf.featFileName_),
179 frozenPointsZone_(ppf.frozenPointsZone_)
180{}
181
182
185(
187)
188:
190 velocity_(ppf.velocity_),
191 featFileName_(ppf.featFileName_),
192 frozenPointsZone_(ppf.frozenPointsZone_)
193{}
194
195
198(
201)
202:
203 pointPatchVectorField(ppf, iF),
204 velocity_(ppf.velocity_),
205 featFileName_(ppf.featFileName_),
206 frozenPointsZone_(ppf.frozenPointsZone_)
207{}
209
210// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211
214{
215 if (!pointTreePtr_)
216 {
217 const Time& tm = this->patch().boundaryMesh().mesh().time();
218 const auto& eMesh = tm.lookupObject<edgeMesh>(featFileName_);
219
220 const pointField& points = eMesh.points();
221
222 // Calculate bb of all points
224
225 // Random number generator. Bit dodgy since not exactly random ;-)
226 Random rndGen(65431);
227
228 // Slightly extended bb. Slightly off-centred just so on symmetric
229 // geometry there are less face/edge aligned items.
230 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
231
232 pointTreePtr_.reset
233 (
235 (
236 treeDataPoint(points), // All edges
237
238 bb, // overall search domain
239 8, // maxLevel
240 10, // leafsize
241 3.0 // duplicity
242 )
243 );
244 }
246 return pointTreePtr_();
247}
248
249
251{
252 if (this->updated())
253 {
254 return;
255 }
256
257 const vectorField currentDisplacement(this->patchInternalField());
258
259 // Calculate displacement to project points onto surface
260 vectorField displacement(currentDisplacement);
261 calcProjection(displacement);
262
263
264 // offset wrt current displacement
265 vectorField offset(displacement-currentDisplacement);
266
267 // Clip offset to maximum displacement possible: velocity*timestep
268
269 const Time& tm = this->patch().boundaryMesh().mesh().time();
270 const scalar deltaT = tm.deltaTValue();
271 const vector clipVelocity = velocity_*deltaT;
272
273 forAll(displacement, i)
274 {
275 vector& d = offset[i];
276
277 const scalar magD(mag(d));
278 if (magD > ROOTVSMALL)
279 {
280 d /= magD;
281 d *= min(magD, mag(clipVelocity));
282 }
283 }
284
285 // Get internal field to insert values into
286 Field<vector>& iF = const_cast<Field<vector>&>(this->primitiveField());
287
288 setInInternalField(iF, (currentDisplacement+offset)());
291}
292
293
295(
296 Ostream& os
297) const
298{
300 os.writeEntry("file", featFileName_);
301 os.writeEntryIfDifferent<word>
302 (
303 "frozenPointsZone",
305 frozenPointsZone_
306 );
307 os.writeEntry("velocity", velocity_);
308}
309
310
311// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312
313namespace Foam
314{
315
317(
320);
321
322} // End namespace Foam
323
324// ************************************************************************* //
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...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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
Random number generator.
Definition Random.H:56
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Area discretisation.
Definition edgeFaMesh.H:50
static void read(const objectRegistry &obr, const dictionary &dict)
Read (& store) geometry. Exposed so point attraction can reuse it.
A class for handling file names.
Definition fileName.H:75
Non-pointer based hierarchical recursive searching.
A class for handling keywords in dictionaries.
Definition keyType.H:69
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const word & name() const noexcept
The patch name.
Displacement by attraction to nearest point. Use in a displacementMotionSolver as a bc on the pointDi...
pointAttractionDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const pointMesh & mesh() const noexcept
Return the mesh reference.
const pointPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
void patchInternalField(const UList< Type1 > &internalData, const labelUList &addressing, UList< Type1 > &pfld) const
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelUList &meshPoints) const
const Field< vector > & primitiveField() const noexcept
virtual void write(Ostream &os) const
Write.
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
virtual const labelList & meshPoints() const =0
Return mesh points.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition pointPatch.H:250
Standard boundBox with extra functionality for use in octree.
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
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)
const pointField & points
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)
messageStream Info
Information stream (stdout output on master, null elsewhere).
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
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.
pointPatchField< vector > pointPatchVectorField
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
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.
Random rndGen
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))