lumpedPointForces.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 -------------------------------------------------------------------------------
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 Application
27  lumpedPointForces
28 
29 Description
30  Extract force/moment information from simulation results that
31  use the lumped points movement description.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "argList.H"
36 #include "Time.H"
37 #include "timeSelector.H"
38 #include "volFields.H"
39 #include "IOobjectList.H"
40 #include "foamVtkSeriesWriter.H"
41 #include "lumpedPointTools.H"
42 #include "lumpedPointIOMovement.H"
43 
44 using namespace Foam;
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 template<class GeoFieldType>
49 autoPtr<GeoFieldType> loadField
50 (
51  const volMesh::Mesh& mesh,
52  const IOobject* io
53 )
54 {
55  if (io && io->isHeaderClass<GeoFieldType>())
56  {
57  Info<< "Reading " << io->headerClassName()
58  << ' ' << io->name() << endl;
59 
61  (
62  IOobject
63  (
64  io->name(),
65  io->instance(),
66  io->local(),
67  io->db(),
71  ),
72  mesh
73  );
74  }
75 
76  return nullptr;
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 int main(int argc, char *argv[])
83 {
85  (
86  "Extract force/moment information from simulation results that"
87  " use the lumped points movement description."
88  );
89 
91  (
92  "vtk",
93  "Create visualization files of the forces"
94  );
95 
96  timeSelector::addOptions(true, false);
97  argList::noFunctionObjects(); // Never use function objects
98 
99  #include "addRegionOption.H"
100  #include "setRootCase.H"
101 
102  const bool withVTK = args.found("vtk");
103 
104  #include "createTime.H"
105 
107 
108  #include "createNamedMesh.H"
109 
111 
112  if (!movement)
113  {
114  Info<< "No valid movement found" << endl;
115  return 1;
116  }
117 
119  if (!nPatches)
120  {
121  Info<< "No point patches with lumped movement found" << endl;
122  return 2;
123  }
124 
125  Info<<"Lumped point patch controls set on " << nPatches
126  << " patches" << nl;
127 
128 
129  vtk::seriesWriter forceSeries;
130  List<vector> forces, moments;
131 
132  forAll(timeDirs, timei)
133  {
134  runTime.setTime(timeDirs[timei], timei);
135 
136  Info<< "Time = " << runTime.timeName() << endl;
137 
138  if (mesh.readUpdate())
139  {
140  Info<< " Read new mesh" << nl;
141  }
142 
143  // Search for list of objects for this time
144  IOobjectList objects(mesh, runTime.timeName());
145 
146  // Pressure field
148  = loadField<volScalarField>(mesh, objects.findObject("p"));
149 
150  // The forces per zone
151  if (movement().forcesAndMoments(mesh, forces, moments))
152  {
153  Info<<"forces per zone: " << forces << endl;
154  Info<<"moments per zone: " << moments << endl;
155 
156  if (withVTK && Pstream::master())
157  {
158  const word outputName =
159  word::printf("forces_%06d.vtp", runTime.timeIndex());
160 
161  Info<<" " << outputName << endl;
162 
163  movement().writeForcesAndMomentsVTP
164  (
165  outputName,
166  forces,
167  moments
168  );
169 
170  forceSeries.append(runTime.timeIndex(), outputName);
171  }
172  }
173  }
174 
175 
176  // Create file series
177 
178  if (forceSeries.size())
179  {
180  forceSeries.write("forces.vtp");
181  }
182 
183  Info<< "\nEnd\n" << endl;
184 
185  return 0;
186 }
187 
188 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:547
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
label size() const noexcept
The number of data sets.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Required Classes.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Required Classes.
word outputName("finiteArea-edges.obj")
dynamicFvMesh & mesh
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:704
Provides a means of accumulating and generating VTK file series.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:935
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:234
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
const word & headerClassName() const noexcept
Return name of the class name read from header.
Definition: IOobjectI.H:213
label setPatchControls(const pointVectorField &pvf, const pointField &points0)
Return the patch-ids associated with a "lumpedPointDisplacement" type.
static void write(const fileName &base, const UList< instant > &series, const char sep='_')
Write file series (JSON format) to disk, for specified instances.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
Automatically write from objectRegistry::writeObject()
static autoPtr< lumpedPointIOMovement > New(const objectRegistry &obr, label ownerId=-1)
Create a movement object in the registry by reading system dictionary.
messageStream Info
Information stream (stdout output on master, null elsewhere)
bool registerObject() const noexcept
Should objects created with this IOobject be registered?
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
const fileName & local() const noexcept
Read access to local path component.
Definition: IOobjectI.H:278
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
bool append(const fileNameInstant &inst)
Append the specified file instant.
bool isHeaderClass() const
Check if headerClassName() equals Type::typeName.
Definition: IOobjectI.H:258