layerAverage.H
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) 2011-2022 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 Class
25  Foam::layerAverage
26 
27 Description
28  Generates plots of fields averaged over the layers in the mesh
29 
30  Example of function object specification:
31  \verbatim
32  layerAverage1
33  {
34  type layerAverage;
35  libs ("libfieldFunctionObjects.so");
36 
37  writeControl writeTime;
38 
39  setFormat raw;
40 
41  patches (bottom);
42  zones (quarterPlane threeQuartersPlane);
43 
44  axis y;
45 
46  symmetric true;
47 
48  fields (pMean pPrime2Mean UMean UPrime2Mean k);
49  }
50  \endverbatim
51 
52 Usage
53  \table
54  Property | Description | Required | Default value
55  type | Type name: layerAverage | yes |
56  setFormat | Output plotting format | yes |
57  patches | Patches that layers extrude from | no | ()
58  zones | Face zones that the layers extrude from | no | ()
59  axis | Component of the position to plot against | yes |
60  symmetric | Is the geometry symmetric around the centre layer? \
61  | no | false
62  field | Fields to average and plot | yes |
63  \endtable
64 
65 SourceFiles
66  layerAverage.C
67  layerAverageTemplates.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef layerAverage_H
72 #define layerAverage_H
73 
74 #include "fvMeshFunctionObject.H"
75 #include "setWriter.H"
76 #include "boolList.H"
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 namespace functionObjects
83 {
84 
85 /*---------------------------------------------------------------------------*\
86  Class layerAverage Declaration
87 \*---------------------------------------------------------------------------*/
88 
89 class layerAverage
90 :
91  public fvMeshFunctionObject
92 {
93  // Private Data
94 
95  //- Patches which form the start of the layers
96  labelList patchIDs_;
97 
98  //- Zones which form the start of the layers
99  labelList zoneIDs_;
100 
101  //- Zones on which the layers are considered to end
102  labelList endZoneIDs_;
103 
104  //- Is the case symmetric?
105  bool symmetric_;
106 
107  //- The direction over which to plot the results
108  coordSet::axisType axis_;
109 
110  //- Per cell the global layer
111  label nLayers_;
112 
113  //- Per cell the global layer
114  labelList cellLayer_;
115 
116  //- Per global layer the number of cells
117  scalarField layerCount_;
118 
119  //- The average centre of each layer
120  pointField layerCentre_;
121 
122  //- Fields to sample
123  wordList fields_;
124 
125  //- Set formatter
126  autoPtr<setWriter> formatter_;
127 
128 
129  // Private Member Functions
130 
131  //- Create the layer information, the sort map, and the scalar axis
132  void calcLayers();
133 
134  //- Return the coefficient to multiply onto symmetric values
135  template<class T>
136  T symmetricCoeff() const;
137 
138  //- Sum field per layer
139  template<class T>
140  Field<T> sum(const Field<T>& cellField) const;
141 
142  //- Average a field per layer
143  template<class T>
144  Field<T> average(const Field<T>& cellField) const;
145 
146 
147 public:
148 
149  //- Runtime type information
150  TypeName("layerAverage");
151 
152 
153  // Constructors
154 
155  //- Construct from Time and dictionary
157  (
158  const word& name,
159  const Time& runTime,
160  const dictionary& dict
161  );
162 
163  //- Disallow default bitwise copy construction
164  layerAverage(const layerAverage&) = delete;
165 
166 
167  //- Destructor
168  virtual ~layerAverage();
169 
170 
171  // Member Functions
172 
173  //- Read the field average data
174  virtual bool read(const dictionary&);
175 
176  //- Return the list of fields required
177  virtual wordList fields() const;
178 
179  //- Do nothing
180  virtual bool execute();
181 
182  //- Calculate and write the graphs
183  virtual bool write();
184 
185  //- Update for mesh point-motion
186  virtual void movePoints(const polyMesh&);
187 
188  //- Update topology using the given map
189  virtual void topoChange(const polyTopoChangeMap&);
190 
191  //- Update from another mesh using the given map
192  virtual void mapMesh(const polyMeshMap&);
193 
194  //- Redistribute or update using the given distribution map
195  virtual void distribute(const polyDistributionMap&);
196 
197 
198  // Member Operators
199 
200  //- Disallow default bitwise assignment
201  void operator=(const layerAverage&) = delete;
202 };
203 
204 
205 template<>
206 vector layerAverage::symmetricCoeff<vector>() const;
207 
208 template<>
209 symmTensor layerAverage::symmetricCoeff<symmTensor>() const;
210 
211 template<>
212 tensor layerAverage::symmetricCoeff<tensor>() const;
213 
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace functionObjects
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #ifdef NoRepository
223  #include "layerAverageTemplates.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
axisType
Enumeration defining the output format for coordinates.
Definition: coordSet.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return the name of this functionObject.
virtual wordList fields() const
Return the list of fields required.
Definition: layerAverage.C:256
void operator=(const layerAverage &)=delete
Disallow default bitwise assignment.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: layerAverage.C:368
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: layerAverage.C:391
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: layerAverage.C:380
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: layerAverage.C:357
layerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: layerAverage.C:191
TypeName("layerAverage")
Runtime type information.
virtual bool execute()
Do nothing.
Definition: layerAverage.C:262
virtual bool write()
Calculate and write the graphs.
Definition: layerAverage.C:268
virtual ~layerAverage()
Destructor.
Definition: layerAverage.C:205
virtual bool read(const dictionary &)
Read the field average data.
Definition: layerAverage.C:211
Generates plots of fields averaged over the layers in the mesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict