lumpedPointState.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-2019 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 \*---------------------------------------------------------------------------*/
27 
28 #include "lumpedPointState.H"
29 #include "demandDrivenData.H"
30 #include "unitConversion.H"
32 #include "IFstream.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 const Foam::Enum
37 <
39 >
41 ({
42  { inputFormatType::PLAIN, "plain" },
43  { inputFormatType::DICTIONARY, "dictionary" },
44 });
45 
46 
47 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
48 
50 static Foam::string getLineNoComment
51 (
52  Foam::ISstream& is,
53  const char comment = '#'
54 )
55 {
56  Foam::string line;
57  do
58  {
59  is.getLine(line);
60  }
61  while ((line.empty() || line[0] == comment) && is.good());
62 
63  return line;
64 }
66 
67 
68 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69 
70 void Foam::lumpedPointState::calcRotations() const
71 {
72  rotationPtr_ = new tensorField(angles_.size());
73 
74  auto rotIter = rotationPtr_->begin();
75 
76  for (const vector& angles : angles_)
77  {
78  *rotIter =
80 
81  ++rotIter;
82  }
83 }
84 
85 
86 void Foam::lumpedPointState::readDict(const dictionary& dict)
87 {
88  dict.readEntry("points", points_);
89  dict.readEntry("angles", angles_);
90  order_ =
92  (
93  "order",
94  dict,
95  quaternion::eulerOrder::ZXZ
96  );
97  degrees_ = dict.lookupOrDefault("degrees", false);
98 
99  deleteDemandDrivenData(rotationPtr_);
100 }
101 
102 
103 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
104 
106 :
107  points_(),
108  angles_(),
109  order_(quaternion::eulerOrder::ZXZ),
110  degrees_(false),
111  rotationPtr_(nullptr)
112 {}
113 
114 
116 :
117  points_(rhs.points_),
118  angles_(rhs.angles_),
119  order_(rhs.order_),
120  degrees_(rhs.degrees_),
121  rotationPtr_(nullptr)
122 {}
123 
124 
126 :
127  points_(pts),
128  angles_(points_.size(), Zero),
129  order_(quaternion::eulerOrder::ZXZ),
130  degrees_(false),
131  rotationPtr_(nullptr)
132 {}
133 
134 
136 :
137  points_(pts),
138  angles_(points_.size(), Zero),
139  order_(quaternion::eulerOrder::ZXZ),
140  degrees_(false),
141  rotationPtr_(nullptr)
142 {}
143 
144 
146 :
147  points_(),
148  angles_(),
149  order_(quaternion::eulerOrder::ZXZ),
150  degrees_(false),
151  rotationPtr_(nullptr)
152 {
153  readDict(dict);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
158 
160 {
161  deleteDemandDrivenData(rotationPtr_);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 {
169  points_ = rhs.points_;
170  angles_ = rhs.angles_;
171  order_ = rhs.order_;
172  degrees_ = rhs.degrees_;
173 
174  deleteDemandDrivenData(rotationPtr_);
175 }
176 
177 
178 void Foam::lumpedPointState::scalePoints(const scalar scaleFactor)
179 {
180  if (scaleFactor > 0)
181  {
182  points_ *= scaleFactor;
183  }
184 }
185 
186 
188 (
189  const scalar alpha,
190  const lumpedPointState& prev
191 )
192 {
193  points_ = prev.points_ + alpha*(points_ - prev.points_);
194 
195  scalar convert = 1.0;
196  if (degrees_ != prev.degrees_)
197  {
198  if (prev.degrees_)
199  {
200  // Was degrees, now radians
201  convert = degToRad();
202  }
203  else
204  {
205  // Was radians, now degrees
206  convert = radToDeg();
207  }
208  }
209 
210  angles_ = convert*prev.angles_ + alpha*(angles_ - convert*prev.angles_);
211 
212  deleteDemandDrivenData(rotationPtr_);
213 }
214 
215 
217 {
218  // Assume generic input stream so we can do line-based input
219  ISstream& iss = dynamic_cast<ISstream&>(is);
220 
221  string line = getLineNoComment(iss);
222 
223  label count;
224  {
225  IStringStream isstr(line);
226  isstr >> count;
227  }
228 
229  points_.setSize(count);
230  angles_.setSize(count);
231 
232  count = 0;
233  forAll(points_, i)
234  {
235  iss.getLine(line);
236  IStringStream isstr(line);
237 
238  isstr
239  >> points_[count].x() >> points_[count].y() >> points_[count].z()
240  >> angles_[count].x() >> angles_[count].y() >> angles_[count].z();
241 
242  ++count;
243  }
244 
245  points_.setSize(count);
246  angles_.setSize(count);
247  order_ = quaternion::eulerOrder::ZXZ;
248  degrees_ = false;
249 
250  deleteDemandDrivenData(rotationPtr_);
251 
252  return count;
253 }
254 
255 
257 {
258  dictionary dict(is);
259  readDict(dict);
260 
261  return points_.size();
262 }
263 
264 
266 {
267  writeDict(os);
268  return true;
269 }
270 
271 
273 {
274  os.writeEntry("points", points_);
275  os.writeEntry("angles", angles_);
276  if (order_ != quaternion::eulerOrder::ZXZ)
277  {
278  os.writeEntry("order", quaternion::eulerOrderNames[order_]);
279  }
280  if (degrees_)
281  {
282  os.writeEntry("degrees", "true");
283  }
284 }
285 
286 
288 {
289  os <<"# input for OpenFOAM\n"
290  <<"# N, points, angles\n"
291  << points_.size() << "\n";
292 
293  forAll(points_, i)
294  {
295  const vector& pt = points_[i];
296 
297  os << pt.x() << ' ' << pt.y() << ' ' << pt.z();
298 
299  if (i < angles_.size())
300  {
301  os << ' ' << angles_[i].x()
302  << ' ' << angles_[i].y()
303  << ' ' << angles_[i].z() << '\n';
304  }
305  else
306  {
307  os << " 0 0 0\n";
308  }
309  }
310 }
311 
312 
314 (
315  const inputFormatType& fmt,
316  const fileName& file
317 )
318 {
319  bool ok = false;
320  if (Pstream::master())
321  {
322  IFstream is(file);
323 
324  if (fmt == inputFormatType::PLAIN)
325  {
326  ok = this->readPlain(is);
327  }
328  else
329  {
330  ok = this->readData(is);
331  }
332  }
333 
334  if (Pstream::parRun())
335  {
336  // Scatter master data using communication scheme
337 
338  const List<Pstream::commsStruct>& comms =
339  (
343  );
344 
345  // Get my communication order
346  const Pstream::commsStruct& myComm = comms[Pstream::myProcNo()];
347 
348  // Receive from up
349  if (myComm.above() != -1)
350  {
351  IPstream fromAbove
352  (
354  myComm.above(),
355  0,
358  );
359 
360  fromAbove >> points_ >> angles_ >> degrees_;
361  }
362 
363  // Send to downstairs neighbours
364  forAllReverse(myComm.below(), belowI)
365  {
366  OPstream toBelow
367  (
369  myComm.below()[belowI],
370  0,
373  );
374 
375  toBelow << points_ << angles_ << degrees_;
376  }
377 
378  deleteDemandDrivenData(rotationPtr_);
379 
380  // MPI barrier
381  Pstream::scatter(ok);
382  }
383 
384  return ok;
385 }
386 
387 
388 // ************************************************************************* //
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
EulerCoordinateRotation.H
Foam::lumpedPointState::writeDict
void writeDict(Ostream &os) const
Output as dictionary content.
Definition: lumpedPointState.C:272
Foam::ISstream::getLine
ISstream & getLine(std::string &str, char delim='\n')
Raw, low-level getline into a string function.
Definition: ISstreamI.H:78
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::IFstream
Input from file stream, using an ISstream.
Definition: IFstream.H:97
Foam::lumpedPointState::inputFormatType
inputFormatType
Input format types.
Definition: lumpedPointState.H:116
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::lumpedPointState::operator=
void operator=(const lumpedPointState &rhs)
Assignment operator.
Definition: lumpedPointState.C:167
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::lumpedPointState::relax
void relax(const scalar alpha, const lumpedPointState &prev)
Relax the state.
Definition: lumpedPointState.C:188
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::ISstream
Generic input stream using standard (STL) streams.
Definition: ISstream.H:54
Foam::lumpedPointState::~lumpedPointState
virtual ~lumpedPointState()
Destructor.
Definition: lumpedPointState.C:159
unitConversion.H
Unit conversion functions.
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:73
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
lumpedPointState.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::quaternion
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:56
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::lumpedPointState
The state of lumped points corresponds to positions and rotations.
Definition: lumpedPointState.H:111
Foam::lumpedPointState::scalePoints
void scalePoints(const scalar scaleFactor)
Scale points by given factor.
Definition: lumpedPointState.C:178
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::lumpedPointState::formatNames
static const Enum< inputFormatType > formatNames
Names for the input format types.
Definition: lumpedPointState.H:125
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::lumpedPointState::angles
const vectorField & angles() const
The orientation of the points (mass centres)
Definition: lumpedPointStateI.H:52
Foam::Enum::getOrDefault
EnumType getOrDefault(const word &key, const dictionary &dict, const EnumType deflt, const bool failsafe=false) const
Definition: Enum.C:190
Foam::lumpedPointState::writePlain
void writePlain(Ostream &os) const
Output as plain content.
Definition: lumpedPointState.C:287
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
IFstream.H
Foam::lumpedPointState::readData
bool readData(Istream &is)
Read input as dictionary content.
Definition: lumpedPointState.C:256
Foam::UPstream::commsTypes::scheduled
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::lumpedPointState::writeData
bool writeData(Ostream &os) const
Output as dictionary content.
Definition: lumpedPointState.C:265
Foam::IStringStream
Input from string buffer, using a ISstream.
Definition: StringStream.H:112
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::UPstream::commsStruct::above
label above() const
Definition: UPstream.H:129
Foam::UPstream::commsStruct
Structure for communicating between processors.
Definition: UPstream.H:80
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::UPstream::nProcsSimpleSum
static int nProcsSimpleSum
Definition: UPstream.H:270
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
Foam::quaternion::eulerOrderNames
static const Enum< eulerOrder > eulerOrderNames
The names for Euler-angle rotation order.
Definition: quaternion.H:112
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:74
Foam::Vector< scalar >
Foam::UPstream::worldComm
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:285
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::line
A line primitive.
Definition: line.H:59
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:303
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::coordinateRotations::euler::rotation
static tensor rotation(const vector &angles, bool degrees=false)
Definition: EulerCoordinateRotation.C:239
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::lumpedPointState::lumpedPointState
lumpedPointState()
Construct null.
Definition: lumpedPointState.C:105
Foam::UPstream::linearCommunication
static const List< commsStruct > & linearCommunication(const label communicator=0)
Communication schedule for linear all-to-master (proc 0)
Definition: UPstream.H:474
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:216
Foam::lumpedPointState::readPlain
bool readPlain(Istream &is)
Read input as plain content.
Definition: lumpedPointState.C:216
Foam::UPstream::commsStruct::below
const labelList & below() const
Definition: UPstream.H:134
Foam::UPstream::treeCommunication
static const List< commsStruct > & treeCommunication(const label communicator=0)
Communication schedule for tree all-to-master (proc 0)
Definition: UPstream.H:483