pointHistory.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) 2019 Zeljko Tukovic, FSB Zagreb.
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 "pointHistory.H"
30 #include "IOmanip.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(pointHistory, 0);
37 
39  (
40  functionObject,
41  pointHistory,
42  dictionary
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 bool Foam::pointHistory::writeData()
50 {
51  const fvMesh& mesh = time_.lookupObject<fvMesh>(polyMesh::defaultRegion);
52 
53  vector position(Zero);
54 
55  if (processor_ == Pstream::myProcNo())
56  {
57  position = mesh.points()[historyPointID_];
58  }
59 
60  reduce(position, sumOp<vector>());
61 
62  if (Pstream::master())
63  {
64  historyFilePtr_() << setprecision(12);
65 
66  historyFilePtr_()
67  << time_.time().value() << tab
68  << position.x() << tab
69  << position.y() << tab
70  << position.z() << endl;
71  }
72 
73  return true;
74 }
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
78 Foam::pointHistory::pointHistory
79 (
80  const word& name,
81  const Time& runTime,
82  const dictionary& dict
83 )
84 :
86  name_(name),
87  time_(runTime),
88  regionName_(polyMesh::defaultRegion),
89  historyPointID_(-1),
90  refHistoryPoint_(dict.lookup("refHistoryPoint")),
91  processor_(-1),
92  fileName_(dict.get<word>("fileName")),
93  historyFilePtr_(nullptr)
94 {
95  Info<< "Creating " << this->name() << " function object." << endl;
96 
97  dict.readIfPresent("region", regionName_);
98  dict.readIfPresent("historyPointID", historyPointID_);
99 
100  const fvMesh& mesh = time_.lookupObject<fvMesh>(regionName_);
101 
102  const vectorField& points = mesh.points();
103 
104  List<scalar> minDist(Pstream::nProcs(), GREAT);
105 
106  if (historyPointID_ == -1)
107  {
108  forAll(points, pointI)
109  {
110  scalar dist = mag(refHistoryPoint_ - points[pointI]);
111 
112  if (dist < minDist[Pstream::myProcNo()])
113  {
114  minDist[Pstream::myProcNo()] = dist;
115  historyPointID_ = pointI;
116  }
117  }
118  }
119 
120  Pstream::allGatherList(minDist);
121 
122  processor_ = -1;
123  scalar min = GREAT;
124 
125  forAll(minDist, procI)
126  {
127  if (minDist[procI] < min)
128  {
129  min = minDist[procI];
130  processor_ = procI;
131  }
132  }
133 
134  if (processor_ == Pstream::myProcNo())
135  {
136  Pout<< "History point ID: " << historyPointID_ << nl
137  << "History point coordinates: "
138  << points[historyPointID_] << nl
139  << "Reference point coordinates: " << refHistoryPoint_
140  << endl;
141  }
142 
143  // Create history file if not already created
144  if (!historyFilePtr_)
145  {
146  // File update
147  if (Pstream::master())
148  {
149  const word startTimeName =
150  time_.timeName(mesh.time().startTime().value());
151 
152  const fileName historyDir =
153  time_.globalPath()/"history"/startTimeName;
154 
155  // Create directory if does not exist.
156  mkDir(historyDir);
157 
158 
159  // Open new file at start up
160 
161  // OStringStream FileName;
162  // FileName() << "point_" << historyPointID_ << ".dat";
163 
164  historyFilePtr_.reset
165  (
166  new OFstream(historyDir/fileName_)
167  );
168 
169  // Add headers to output data
170  if (historyFilePtr_)
171  {
172  historyFilePtr_()
173  << "# Time" << tab << "X" << tab << "Y" << tab << "Z"
174  << endl;
175  }
176  }
177  }
178 
179  // Write start time data
180  writeData();
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 {
188  return writeData();
189 }
190 
191 
193 {
194  dict.readIfPresent("region", regionName_);
195 
196  return true;
197 }
198 
199 
200 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
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...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static void writeData(Ostream &os, const Type &val)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
Abstract base-class for Time/database function objects.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Omanip< int > setprecision(const int i)
Definition: IOmanip.H:205
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1077
virtual bool execute()
execute is called at each ++ or += of the time-loop
Definition: pointHistory.C:179
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Istream and Ostream manipulators taking arguments.
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
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
virtual bool read(const dictionary &dict)
Read and set the function object if its data has changed.
Definition: pointHistory.C:185
defineTypeNameAndDebug(combustionModel, 0)
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual linear/tree communicat...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127