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-2024 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  patch | Patch that layers extrude from | no |
58  patches | Patches that layers extrude from | no | ()
59  zones | Face zones that the layers extrude from | no | ()
60  axis | Component of the position to plot against | yes |
61  symmetric | Is the geometry symmetric around the centre layer? \
62  | no | false
63  fields | Fields to average and plot | yes |
64  weightField | Field with which to weight the average | no | none
65  weightFields | Fields with which to weight the average | no | ()
66  \endtable
67 
68 SourceFiles
69  layerAverage.C
70  layerAverageTemplates.C
71 
72 \*---------------------------------------------------------------------------*/
73 
74 #ifndef layerAverage_H
75 #define layerAverage_H
76 
77 #include "fvMeshFunctionObject.H"
78 #include "setWriter.H"
79 #include "boolList.H"
80 #include "volFields.H"
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 namespace functionObjects
87 {
88 
89 /*---------------------------------------------------------------------------*\
90  Class layerAverage Declaration
91 \*---------------------------------------------------------------------------*/
92 
93 class layerAverage
94 :
95  public fvMeshFunctionObject
96 {
97  // Private Data
98 
99  //- Patches which form the start of the layers
100  labelList patchIndices_;
101 
102  //- Zones which form the start of the layers
103  labelList zoneIndices_;
104 
105  //- Is the case symmetric?
106  bool symmetric_;
107 
108  //- The direction over which to plot the results
109  coordSet::axisType axis_;
110 
111  //- Per cell the global layer
112  label nLayers_;
113 
114  //- Per cell the global layer
115  labelList cellLayer_;
116 
117  //- Per global layer the volume
118  scalarField layerVolume_;
119 
120  //- The average centre of each layer
121  pointField layerCentre_;
122 
123  //- Fields to average
124  wordList fields_;
125 
126  //- Fields with which to weight the averages
127  wordList weightFields_;
128 
129  //- Set formatter
130  autoPtr<setWriter> formatter_;
131 
132 
133  // Private Member Functions
134 
135  //- Create the layer information, the sort map, and the scalar axis
136  void calcLayers();
137 
138  //- Calculate and return the weight field, or a null pointer if there
139  // are no weight fields
140  tmp<VolInternalField<scalar>> weight() const;
141 
142  //- Return the coefficient to multiply onto symmetric values
143  template<class T>
144  T symmetricCoeff() const;
145 
146  //- Sum field per layer
147  template<class T>
148  tmp<Field<T>> sum(const VolInternalField<T>& cellField) const;
149 
150  //- Average a field per layer
151  template<class T>
152  tmp<Field<T>> average
153  (
154  const tmp<VolInternalField<scalar>>& cellWeight,
155  const tmp<Field<scalar>>& layerWeight,
156  const VolInternalField<T>& cellField
157  ) const;
158 
159 
160 public:
161 
162  //- Runtime type information
163  TypeName("layerAverage");
164 
165 
166  // Constructors
167 
168  //- Construct from Time and dictionary
170  (
171  const word& name,
172  const Time& runTime,
173  const dictionary& dict
174  );
175 
176  //- Disallow default bitwise copy construction
177  layerAverage(const layerAverage&) = delete;
178 
179 
180  //- Destructor
181  virtual ~layerAverage();
182 
183 
184  // Member Functions
185 
186  //- Read the field average data
187  virtual bool read(const dictionary&);
188 
189  //- Return the list of fields required
190  virtual wordList fields() const;
191 
192  //- Do nothing
193  virtual bool execute();
194 
195  //- Calculate and write the graphs
196  virtual bool write();
197 
198  //- Update for mesh point-motion
199  virtual void movePoints(const polyMesh&);
200 
201  //- Update topology using the given map
202  virtual void topoChange(const polyTopoChangeMap&);
203 
204  //- Update from another mesh using the given map
205  virtual void mapMesh(const polyMeshMap&);
206 
207  //- Redistribute or update using the given distribution map
208  virtual void distribute(const polyDistributionMap&);
209 
210 
211  // Member Operators
212 
213  //- Disallow default bitwise assignment
214  void operator=(const layerAverage&) = delete;
215 };
216 
217 
218 template<>
219 vector layerAverage::symmetricCoeff<vector>() const;
220 
221 template<>
222 symmTensor layerAverage::symmetricCoeff<symmTensor>() const;
223 
224 template<>
225 tensor layerAverage::symmetricCoeff<tensor>() const;
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace functionObjects
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #ifdef NoRepository
236  #include "layerAverageTemplates.C"
237 #endif
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
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:162
const word & name() const
Return the name of this functionObject.
virtual wordList fields() const
Return the list of fields required.
Definition: layerAverage.C:287
void operator=(const layerAverage &)=delete
Disallow default bitwise assignment.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: layerAverage.C:410
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: layerAverage.C:433
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: layerAverage.C:422
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: layerAverage.C:399
layerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: layerAverage.C:219
TypeName("layerAverage")
Runtime type information.
virtual bool execute()
Do nothing.
Definition: layerAverage.C:295
virtual bool write()
Calculate and write the graphs.
Definition: layerAverage.C:301
virtual ~layerAverage()
Destructor.
Definition: layerAverage.C:233
virtual bool read(const dictionary &)
Read the field average data.
Definition: layerAverage.C:239
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 managing temporary objects.
Definition: tmp.H:55
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.
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict