Loading...
Searching...
No Matches
unequalBinWidth.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) 2022 OpenCFD Ltd.
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 "unequalBinWidth.H"
29#include "histogramModel.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace histogramModels
37{
40}
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
47(
48 const word& name,
49 const fvMesh& mesh,
50 const dictionary& dict
51)
52:
54 nBins_(-1),
55 ranges_()
57 read(dict);
58}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
64{
66 {
67 return false;
68 }
69
70 ranges_ = dict.get<List<scalarMinMax>>("ranges");
71 nBins_ = ranges_.size();
72
73 forAll(ranges_, bini)
74 {
75 const auto& range = ranges_[bini];
76
77 if (!range.good())
78 {
80 << "Histogram bin-" << bini
81 << " has invalid range: " << range
83 }
84 }
85
86 if (nBins_ < 1)
87 {
89 << "Invalid number of histogram bins: " << nBins_
91 }
92
93 return true;
94}
95
96
98{
99 // Retrieve operand field
101
102 // Calculate the mid-points of bins for the graph axis
103 pointField midBin(nBins_, Zero);
104
105 forAll(ranges_, bini)
106 {
107 midBin[bini].x() = ranges_[bini].centre();
108 }
109
110
111 // Calculate the histogram data
112 scalarField dataNormalised(nBins_, Zero);
113 labelField dataCount(nBins_, Zero);
114 const scalarField& V = mesh().V();
115
116 forAll(field, celli)
117 {
118 forAll(ranges_, bini)
119 {
120 const auto& range = ranges_[bini];
121
122 // Like range.contains(field[celli]) but exclusive on max()
123
124 if (field[celli] >= range.min() && field[celli] < range.max())
125 {
126 dataNormalised[bini] += V[celli];
127 dataCount[bini]++;
128 break;
129 }
130 }
131 }
132 Pstream::listGather(dataNormalised, sumOp<scalar>());
133 Pstream::listGather(dataCount, sumOp<label>());
134
135
136 // Write histogram data
137 histogramModel::write(dataNormalised, dataCount, mag(midBin));
138
139 return true;
140}
141
142
143// ************************************************************************* //
scalar range
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
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.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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 unequal widths.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
unequalBinWidth(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
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
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
A namespace for various histogram model implementations.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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 ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299