Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace 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_()
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::listGather(dataNormalised, sumOp<scalar>());
159 Pstream::listGather(dataCount, sumOp<label>());
160
161
162 // Write histogram data
163 histogramModel::write(dataNormalised, dataCount, mag(binMidPoints));
164
165 return true;
166}
167
168
169// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const T & max() const noexcept
The max value.
Definition MinMax.H:217
const T & min() const noexcept
The min value.
Definition MinMax.H:207
bool good() const
Range is non-inverted.
Definition MinMaxI.H:153
void reset()
Reset to an inverted (invalid) range.
Definition MinMaxI.H:160
T span() const
The min to max span. Zero for invalid range.
Definition MinMaxI.H:131
static void listGather(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather (reduce) list elements, applying bop to each list element.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A base class for histogram models.
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
volScalarField & getOrReadField(const word &fieldName) const
Return requested field from the object registry or read+register the field to the object registry.
void write(scalarField &dataNormalised, const labelField &dataCount, const scalarField &magMidBin)
Write histogram data.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
const word & fieldName() const noexcept
Return const reference to the operand field name.
histogramModel(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Histogram model which groups data into bins of equal width.
equalBinWidth(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual bool write(const bool log)
Write data to stream and files.
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
volScalarField & p
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const expr V(m.psi().mesh().V())
A namespace for various histogram model implementations.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299