equalBinWidth.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-2022 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 "equalBinWidth.H"
30 #include "histogramModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace histogramModels
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const word& name,
50  const fvMesh& mesh,
51  const dictionary& dict
52 )
53 :
55  nBins_(0),
56  range_()
57 {
58  read(dict);
59 }
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
65 {
67  {
68  return false;
69  }
70 
71  range_.reset
72  (
73  dict.getOrDefault<scalar>("min", GREAT),
74  dict.getOrDefault<scalar>("max", -GREAT)
75  );
76 
77  nBins_ = dict.get<scalar>("nBins");
78 
79  if (nBins_ < 1)
80  {
82  << "Number of histogram bins = " << nBins_
83  << " cannot be negative or zero."
85  }
86 
87  return true;
88 }
89 
90 
92 {
93  // Retrieve operand field
95 
96  // Determine min and max from the operand field
97  // if the user did not provide any min or max
98 
99  scalarMinMax histRange(range_);
100 
101  if (histRange.max() == -GREAT)
102  {
103  histRange.max() = max(field).value();
104 
105  if (histRange.min() == GREAT)
106  {
107  histRange.min() = min(field).value();
108  }
109 
110  if (log)
111  {
112  Info<< " Determined histogram bounds from field"
113  << " min/max(" << fieldName() << ") = "
114  << histRange << endl;
115  }
116  }
117  else if (histRange.min() == GREAT)
118  {
119  histRange.min() = Zero;
120  }
121 
122  if (!histRange.good())
123  {
125  << "Invalid histogram range: " << histRange
126  << exit(FatalError);
127  }
128 
129 
130  // Calculate the mid-points of bins for the graph axis
131  pointField binMidPoints(nBins_, Zero);
132  const scalar delta = histRange.span()/nBins_;
133 
134  {
135  scalar x = histRange.min() + 0.5*delta;
136  for (point& p : binMidPoints)
137  {
138  p.x() = x;
139  x += delta;
140  }
141  }
142 
143 
144  // Calculate the histogram data
145  scalarField dataNormalised(nBins_, Zero);
146  labelField dataCount(nBins_, Zero);
147  const scalarField& V = mesh().V();
148 
149  forAll(field, celli)
150  {
151  const label bini = (field[celli] - histRange.min())/delta;
152  if (bini >= 0 && bini < nBins_)
153  {
154  dataNormalised[bini] += V[celli];
155  dataCount[bini]++;
156  }
157  }
158  Pstream::listCombineGather(dataNormalised, plusEqOp<scalar>());
159  Pstream::listCombineGather(dataCount, plusEqOp<label>());
160 
161 
162  // Write histogram data
163  histogramModel::write(dataNormalised, dataCount, mag(binMidPoints));
164 
165  return true;
166 }
167 
168 
169 // ************************************************************************* //
void write(scalarField &dataNormalised, const labelField &dataCount, const scalarField &magMidBin)
Write histogram data.
scalar delta
dictionary dict
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
defineTypeNameAndDebug(equalBinWidth, 0)
Histogram model which groups data into bins of equal width.
addToRunTimeSelectionTable(histogramModel, equalBinWidth, dictionary)
rDeltaTY field()
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual bool read(const dictionary &dict)
Read top-level dictionary.
Definition: equalBinWidth.C:57
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual bool read(const dictionary &dict)
Read top-level dictionary.
equalBinWidth(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: equalBinWidth.C:41
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void reset()
Reset to an inverted (invalid) range.
Definition: MinMaxI.H:172
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A base class for histogram models.
virtual bool write(const bool log)
Write data to stream and files.
Definition: equalBinWidth.C:84
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
vector point
Point is a vector.
Definition: point.H:37
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volScalarField & getOrReadField(const word &fieldName) const
Return requested field from the object registry or read+register the field to the object registry...
messageStream Info
Information stream (stdout output on master, null elsewhere)
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
volScalarField & p
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127