Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
28#include "pointHistory.H"
30#include "IOmanip.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37
39 (
43 );
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49bool 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
78Foam::pointHistory::pointHistory
79(
80 const word& name,
81 const Time& runTime,
82 const dictionary& dict
83)
84:
85 functionObject(name),
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// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
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:1697
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
functionObject(const word &name, const bool withNamePrefix=defaultUseNamePrefix)
Construct from components.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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...
virtual bool read(const dictionary &dict)
Read and set the function object if its data has changed.
virtual bool execute()
execute is called at each ++ or += of the time-loop
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
Lookup type of boundary radiation properties.
Definition lookup.H:60
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
engineTime & runTime
auto & name
const pointField & points
Namespace for OpenFOAM.
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:616
messageStream Info
Information stream (stdout output on master, null elsewhere).
Omanip< int > setprecision(const int i)
Definition IOmanip.H:205
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.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static void writeData(Ostream &os, const Type &val)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299