lumpedPointDisplacementPointPatchVectorField.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) 2016-2020 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
29#include "lumpedPointMovement.H"
32#include "pointFields.H"
33#include "surfaceFields.H"
34#include "volFields.H"
35#include "Time.H"
36#include "polyMesh.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44 (
47 );
48}
49
50
51// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
52
54(
55 const pointVectorField& pvf,
56 const pointField& points0
57)
58{
59 label count = 0;
60
63
64 forAll(bf, patchi)
65 {
66 // Patch of this type
67 const auto* p = isA<patchType>(bf[patchi]);
68
69 if (p)
70 {
71 // Patch controls (mapping) for calculating forces/moments
72 const_cast<lumpedPointMovement&>(p->movement())
73 .setPatchControl
74 (
75 patches[patchi],
76 p->controllers(),
78 );
79
80 ++count;
81 }
82 }
83
84 return count;
85}
86
87
89(
90 const pointVectorField& pvf,
91 const pointField& points0
92)
93{
94 label count = 0;
95
97
98 forAll(bf, patchi)
99 {
100 // Patch of this type
101 const auto* p = isA<patchType>(bf[patchi]);
102
103 if (p)
104 {
105 // Patch controls (mapping) for calculating forces/moments
106 const_cast<lumpedPointMovement&>(p->movement())
107 .setInterpolator
108 (
109 p->patch(),
110 points0
111 );
112
113 ++count;
114 }
115 }
116
117 return count;
118}
119
120
123(
124 const pointVectorField& pvf
125)
126{
128
129 DynamicList<label> patchLst(bf.size());
130 forAll(bf, patchi)
131 {
132 // Patch of this type
133 if (isA<patchType>(bf[patchi]))
134 {
135 patchLst.append(patchi);
136 // or patchLst.append(bf[patchi].patch().index());
137 }
138 }
139
140 return patchLst.shrink();
141}
142
143
144// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
145
146const Foam::pointField&
148{
149 const objectRegistry& obr = this->patch().boundaryMesh().mesh().db();
150
151 // Obtain starting locations from the motionSolver (when possible)
152 const auto* solver =
153 obr.cfindObject<displacementMotionSolver>("dynamicMeshDict");
154
155 if (solver)
156 {
157 if (points0Ptr_)
158 {
159 points0Ptr_.reset(nullptr);
160 }
161 return solver->points0();
162 }
163 else if (!points0Ptr_)
164 {
165 points0Ptr_.reset
166 (
167 new pointIOField
168 (
170 (
171 this->patch().boundaryMesh().mesh().mesh()
172 )
173 )
174 );
175 }
176
177 return *points0Ptr_;
178}
179
180
183{
184 const objectRegistry& obr = this->patch().boundaryMesh().mesh().db();
185
188
189 if (ptr)
190 {
191 return *ptr; // Already exists
192 }
193
194 // Create and register with this patch as the owner
195 ptr = lumpedPointIOMovement::New(obr, this->patch().index()).ptr();
196
197 return objectRegistry::store(ptr);
198}
199
200
201// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202
203Foam::lumpedPointDisplacementPointPatchVectorField::
204lumpedPointDisplacementPointPatchVectorField
205(
206 const pointPatch& p,
208)
209:
211 controllers_(),
212 dataWritten_(0, 0),
213 points0Ptr_(nullptr)
214{}
215
216
217Foam::lumpedPointDisplacementPointPatchVectorField::
218lumpedPointDisplacementPointPatchVectorField
219(
220 const pointPatch& p,
222 const dictionary& dict
223)
224:
226 controllers_(),
227 dataWritten_(0, 0),
228 points0Ptr_(nullptr)
229{
230 dict.readIfPresent("controllers", controllers_);
231
232 dict.readIfPresent("dataWritten", dataWritten_);
233
234 if (controllers_.empty())
235 {
237 << "No controllers specified, using all lumped points for patch: "
238 << this->patch().name() << nl << nl;
239 }
240
241 // controllers_ : check? remove duplicates?
242}
243
244
245Foam::lumpedPointDisplacementPointPatchVectorField::
246lumpedPointDisplacementPointPatchVectorField
247(
249 const pointPatch& p,
251 const pointPatchFieldMapper& mapper
252)
253:
254 fixedValuePointPatchField<vector>(rhs, p, iF, mapper),
255 controllers_(rhs.controllers_),
256 dataWritten_(rhs.dataWritten_),
257 points0Ptr_(nullptr)
258{}
259
260
261Foam::lumpedPointDisplacementPointPatchVectorField::
262lumpedPointDisplacementPointPatchVectorField
263(
266)
267:
269 controllers_(rhs.controllers_),
270 dataWritten_(rhs.dataWritten_),
271 points0Ptr_(nullptr)
272{}
273
274
275// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
276
279{
280 // de-register movement if in use and managed by this patch
283 (
284 this->patch().boundaryMesh().mesh().db()
285 );
286
287 if (ptr && ptr->ownerId() == this->patch().index())
288 {
289 movement().coupler().shutdown();
290
291 ptr->checkOut();
292 }
293}
294
295
296// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
297
299{
300 if (this->updated())
301 {
302 return;
303 }
304
305 const label timeIndex = this->db().time().timeIndex();
306
308
309 if (movement().ownerId() == this->patch().index())
310 {
311 // The ownerId is always the lowest patch number,
312 // thus will always be triggered first
313
315 {
316 Pout<<"masterPatch: " << this->patch().index() << endl;
317 }
318
319 const polyMesh& mesh = this->patch().boundaryMesh().mesh().mesh();
320
321 // Mapping for calculating forces
322 // as well as face point interpolation
323 if (!movement().hasMapping())
324 {
325 // Add mapping for calculating forces/moments
326 setPatchControls
327 (
328 static_cast<const pointVectorField&>
329 (
330 this->internalField()
331 ),
332 this->points0()
333 );
334 }
335
336 int triggered = 0;
337
338 if
339 (
340 movement().coupler().slaveFirst()
341 && !movement().coupler().initialized()
342 )
343 {
344 // Master does nothing yet, but slave will
345 triggered = -1;
346 }
347 else if (movement().couplingPending(timeIndex))
348 {
349 // Trigger is pending, or coupling not yet initialized
350 triggered = 1;
351 }
352
353 if (triggered > 0)
354 {
355 // Synchronized for all processes
356 List<vector> forces, moments;
357 movement().forcesAndMoments(mesh, forces, moments);
358
360 {
361 Pout<<"gatherForces: " << forces << " called from patch "
362 << this->patch().index() << endl;
363
364 Info<< "output forces to file: called from patch "
365 << this->patch().index() << nl
366 << "# " << forces.size() << " force entries" << nl
367 << "# fx fy fz" << nl
368 << "output forces to file: "
369 << forces << " called from patch "
370 << this->patch().index() << endl;
371 }
372
373 // Update times when data (forces) were written
374 // With first=time, second=prevTime
375
376 dataWritten_.second() = dataWritten_.first();
377 dataWritten_.first() = this->db().time().timeOutputValue();
378
379 if (Pstream::master())
380 {
381 movement().writeData(forces, moments, &dataWritten_);
382
383 // Signal external source to execute
384 movement().coupler().useSlave();
385 }
386 }
387
388 if (triggered)
389 {
390 // Wait for slave to provide data (includes MPI barrier)
391 // and catch any abort information sent from slave
392 action = movement().coupler().waitForSlave();
393
394 const_cast<lumpedPointMovement&>(movement()).readState();
395
396 movement().couplingCompleted(timeIndex);
397 }
398 }
399
400 if (!movement().hasInterpolator(this->patch()))
401 {
402 const_cast<lumpedPointMovement&>(movement()).setInterpolator
403 (
404 this->patch(),
405 this->points0()
406 );
407 }
408
409 tmp<pointField> tdisp =
410 movement().pointsDisplacement
411 (
412 this->patch(),
413 this->points0()
414 );
415
416 this->operator==(tdisp);
417
418
420
421 // Process any abort information sent from slave
422 if
423 (
424 action != this->db().time().stopAt()
426 )
427 {
428 this->db().time().stopAt(action);
429 }
430}
431
432
434const
435{
437
438 if (controllers_.size())
439 {
440 os.writeEntry("controllers", controllers_);
441 }
442
443 // Times when data were written is only meaningful on the owner patch
444 if (movement().ownerId() == this->patch().index())
445 {
446 os.writeEntry("dataWritten", dataWritten_);
447 }
448
449 writeEntry("value", os);
450}
451
452
453// ************************************************************************* //
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...
const Mesh & mesh() const
Return mesh.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
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
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
stopAtControls
Definition: Time.H:98
@ saUnknown
Dummy no-op. Do not change current value.
Definition: Time.H:103
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Virtual base class for displacement motion solver.
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
A FixedValue boundary condition for pointField.
virtual bool write()
Write the output fields.
This is the point-patch responsible for managing the force integration on a 'lumped-point' basis,...
static label setPatchControls(const pointVectorField &pvf, const pointField &points0)
Set all patch controls for patches of this type.
virtual ~lumpedPointDisplacementPointPatchVectorField()
Destructor. De-register movement if in use and managed by this patch.
const pointField & points0() const
The starting locations (obtained from the motionSolver).
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static label setInterpolators(const pointVectorField &pvf, const pointField &points0)
Set all patch controls for patches of this type.
const lumpedPointMovement & movement() const
The auto-vivifying singleton for movement.
IO-registered version of lumpedPointMovement.
static lumpedPointIOMovement * getMovementObject(const objectRegistry &obr)
Find the movement object or nullptr if not found.
The movement driver that describes initial point locations, the current state of the points/rotations...
void couplingCompleted(const label timeIndex) const
Register that coupling is completed at this calcFrequency.
label ownerId() const
An owner Id, if needed for bookkeeping purposes.
static int debug
Debug switch.
Registry of regIOobjects.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Foam::pointPatchFieldMapper.
const pointPatch & patch() const
Return patch.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
virtual const word & name() const =0
Return name.
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
int count() const noexcept
Return the current reference count.
Definition: refCount.H:74
bool checkOut()
Remove all file watches and remove object from registry.
Definition: regIOobject.C:224
Base class for solution control classes.
Definition: solver.H:54
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
labelList patchIds
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))