histogram.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "histogram.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
36  defineTypeNameAndDebug(histogram, 0);
37  addToRunTimeSelectionTable(functionObject, histogram, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::functionObjects::histogram::writeGraph
45 (
46  const coordSet& coords,
47  const word& fieldName,
48  const scalarField& values
49 ) const
50 {
51  const wordList fieldNames(1, fieldName);
52 
53  fileName outputPath = file_.baseTimeDir();
54  mkDir(outputPath);
55  OFstream graphFile
56  (
57  outputPath/formatterPtr_().getFileName(coords, fieldNames)
58  );
59 
60  Log << " Writing histogram of " << fieldName
61  << " to " << graphFile.name() << endl;
62 
63  List<const scalarField*> yPtrs(1);
64  yPtrs[0] = &values;
65  formatterPtr_().write(coords, fieldNames, yPtrs, graphFile);
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const word& name,
74  const Time& runTime,
75  const dictionary& dict
76 )
77 :
78  fvMeshFunctionObject(name, runTime, dict),
79  file_(obr_, name)
80 {
81  read(dict);
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86 
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
94 {
95  dict.lookup("field") >> fieldName_;
96  dict.lookup("max") >> max_;
97  min_ = dict.lookupOrDefault<scalar>("min", 0);
98  dict.lookup("nBins") >> nBins_;
99 
100  word format(dict.lookup("setFormat"));
101  formatterPtr_ = writer<scalar>::New(format);
102 
103  return true;
104 }
105 
106 
108 {
109  return true;
110 }
111 
112 
114 {
115  Log << type() << " " << name() << " write:" << nl;
116 
117  autoPtr<volScalarField> fieldPtr;
118  if (obr_.foundObject<volScalarField>(fieldName_))
119  {
120  Log << " Looking up field " << fieldName_ << endl;
121  }
122  else
123  {
124  Log << " Reading field " << fieldName_ << endl;
125  fieldPtr.reset
126  (
127  new volScalarField
128  (
129  IOobject
130  (
131  fieldName_,
132  mesh_.time().timeName(),
133  mesh_,
136  ),
137  mesh_
138  )
139  );
140  }
141 
142  const volScalarField& field =
143  (
144  fieldPtr.valid()
145  ? fieldPtr()
146  : obr_.lookupObject<volScalarField>(fieldName_)
147  );
148 
149  // Calculate the mid-points of bins for the graph axis
150  pointField xBin(nBins_);
151  const scalar delta = (max_- min_)/nBins_;
152 
153  scalar x = min_ + 0.5*delta;
154  forAll(xBin, i)
155  {
156  xBin[i] = point(x, 0, 0);
157  x += delta;
158  }
159 
160  scalarField volFrac(nBins_, 0);
161  const scalarField& V = mesh_.V();
162 
163  forAll(field, celli)
164  {
165  const label bini = (field[celli] - min_)/delta;
166  if (bini >= 0 && bini < nBins_)
167  {
168  volFrac[bini] += V[celli];
169  }
170  }
171 
173 
174  if (Pstream::master())
175  {
176  const scalar sumVol = sum(volFrac);
177 
178  if (sumVol > small)
179  {
180  volFrac /= sumVol;
181 
182  const coordSet coords
183  (
184  "Volume_Fraction",
185  "x",
186  xBin,
187  mag(xBin)
188  );
189 
190  writeGraph(coords, field.name(), volFrac);
191  }
192  }
193 
194  return true;
195 }
196 
197 
198 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Calculate the histogram and write.
Definition: histogram.C:113
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:59
const word & name() const
Return name.
Definition: IOobject.H:303
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
word format(conversionProperties.lookup("format"))
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
virtual bool execute()
Execute, currently does nothing.
Definition: histogram.C:107
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static List< word > fieldNames
Definition: globalFoam.H:46
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool read(const dictionary &)
Read the histogram data.
Definition: histogram.C:93
histogram(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: histogram.C:72
Holds list of sampling positions.
Definition: coordSet.H:49
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35
static const char nl
Definition: Ostream.H:260
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
rDeltaTY field()
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812