histogram.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 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "histogram.H"
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(histogram, 0);
40  addToRunTimeSelectionTable(functionObject, histogram, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::functionObjects::histogram::writeGraph
48 (
49  const coordSet& coords,
50  const word& fieldName,
51  const scalarField& normalizedValues,
52  const scalarField& absoluteValues
53 ) const
54 {
55  fileName outputPath = baseTimeDir();
56  mkDir(outputPath);
57  OFstream graphFile
58  (
59  outputPath
60  /formatterPtr_().getFileName
61  (
62  coords,
63  wordList(1, fieldName)
64  )
65  );
66 
67  Log << " Writing histogram of " << fieldName
68  << " to " << graphFile.name() << endl;
69 
71  fieldNames[0] = fieldName;
72  fieldNames[1] = fieldName + "Count";
73  List<const scalarField*> yPtrs(2);
74  yPtrs[0] = &normalizedValues;
75  yPtrs[1] = &absoluteValues;
76  formatterPtr_().write(coords, fieldNames, yPtrs, graphFile);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const word& name,
85  const Time& runTime,
86  const dictionary& dict
87 )
88 :
90  writeFile(obr_, name),
91  max_(-GREAT),
92  min_(GREAT)
93 {
94  read(dict);
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
104 
105  dict.readEntry("field", fieldName_);
106 
107  max_ = dict.getOrDefault<scalar>("max", -GREAT);
108  min_ = dict.getOrDefault<scalar>("min", GREAT);
109  dict.readEntry("nBins", nBins_);
110 
111  if (nBins_ < 1)
112  {
114  << "Number of histogram bins = " << nBins_
115  << " cannot be negative or zero."
116  << abort(FatalError);
117  }
118 
119  const word format(dict.get<word>("setFormat"));
120  formatterPtr_ = writer<scalar>::New(format);
121 
122  return true;
123 }
124 
125 
127 {
128  return true;
129 }
130 
131 
133 {
134  Log << type() << " " << name() << " write:" << nl;
135 
136  autoPtr<volScalarField> fieldPtr;
137  if (obr_.foundObject<volScalarField>(fieldName_))
138  {
139  Log << " Looking up field " << fieldName_ << endl;
140  }
141  else
142  {
143  Log << " Reading field " << fieldName_ << endl;
144  fieldPtr.reset
145  (
146  new volScalarField
147  (
148  IOobject
149  (
150  fieldName_,
151  mesh_.time().timeName(),
152  mesh_,
155  ),
156  mesh_
157  )
158  );
159  }
160 
161  const volScalarField& field =
162  (
163  fieldPtr.valid()
164  ? fieldPtr()
165  : obr_.lookupObject<volScalarField>(fieldName_)
166  );
167 
168 
169  scalar histMax = max_;
170  scalar histMin = min_;
171 
172  if (max_ == -GREAT)
173  {
174  // Determine current min and max
175  histMax = max(field).value();
176 
177  if (min_ == GREAT)
178  {
179  histMin = min(field).value();
180  }
181  Log << " Determined histogram bounds from field"
182  << " min/max(" << fieldName_ << ") = "
183  << histMin << ' ' << histMax << endl;
184  }
185  else if (min_ == GREAT)
186  {
187  histMin = 0;
188  }
189 
190  // Calculate the mid-points of bins for the graph axis
191  pointField xBin(nBins_);
192  const scalar delta = (histMax- histMin)/nBins_;
193 
194  scalar x = histMin + 0.5*delta;
195  forAll(xBin, i)
196  {
197  xBin[i] = point(x, 0, 0);
198  x += delta;
199  }
200 
201  scalarField dataNormalized(nBins_, Zero);
202  labelField dataCount(nBins_, Zero);
203  const scalarField& V = mesh_.V();
204 
205  forAll(field, celli)
206  {
207  const label bini = (field[celli] - histMin)/delta;
208  if (bini >= 0 && bini < nBins_)
209  {
210  dataNormalized[bini] += V[celli];
211  dataCount[bini]++;
212  }
213  }
214 
215  Pstream::listCombineGather(dataNormalized, plusEqOp<scalar>());
217 
218  if (Pstream::master())
219  {
220  const scalar sumData = sum(dataNormalized);
221 
222  if (sumData > SMALL)
223  {
224  dataNormalized /= sumData;
225 
226  const coordSet coords
227  (
228  fieldName_,
229  "x",
230  xBin,
231  mag(xBin)
232  );
233 
234 
235  // Convert count field from labelField to scalarField
236  scalarField count(dataCount.size());
237  forAll(count, i)
238  {
239  count[i] = 1.0*dataCount[i];
240  }
241 
242  writeGraph(coords, fieldName_, dataNormalized, count);
243  }
244  }
245 
246  return true;
247 }
248 
249 
250 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::functionObjects::histogram::write
virtual bool write()
Calculate the histogram and write.
Definition: histogram.C:132
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
histogram.H
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::autoPtr::valid
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: autoPtr.H:272
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::histogram::execute
virtual bool execute()
Execute (effectively no-op)
Definition: histogram.C:126
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
format
word format(conversionProperties.get< word >("format"))
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:213
Foam::Field< vector >
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
field
rDeltaTY field()
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::functionObjects::histogram::read
virtual bool read(const dictionary &)
Read the histogram data.
Definition: histogram.C:100
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::functionObjects::histogram::histogram
histogram(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: histogram.C:83
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
Foam::plusEqOp
Definition: ops.H:72
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::functionObjects::writeFile::baseTimeDir
fileName baseTimeDir() const
Return the base directory for the current time value.
Definition: writeFile.C:76
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38