cellPoint_LagrangianAverage.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) 2025 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::LagrangianAverage
26 
27 Description
28  Lagrangian averaging process in which values are averaged onto the mesh
29  cell-centres and points using the basis functions associated with the
30  tetrahedral decomposition. Interpolation back to the Lagrangian points is
31  then done with the same basis functions.
32 
33 SourceFiles
34  cellPoint_LagrangianAverage.C
35  cellPoint_LagrangianAverages.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef cellPoint_LagrangianAverage_H
40 #define cellPoint_LagrangianAverage_H
41 
42 #include "LagrangianAverage.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 class polyMesh;
50 
51 namespace LagrangianAverages
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class LagrangianAverage Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 template<class Type>
59 class cellPoint
60 :
61  public LagrangianAverage<Type>
62 {
63 private:
64 
65  // Private Classes
66 
67  //- Data structure for an average
68  struct data
69  {
70  //- Map from the cell to the cell average
71  List<label> cellCellAvg_;
72 
73  //- Map from the point to the point average
74  List<label> pointPointAvg_;
75 
76  //- Map from the cell average to the cell
77  DynamicList<label> cellAvgCell_;
78 
79  //- Map from the point average to the point
80  DynamicList<label> pointAvgPoint_;
81 
82  //- Number of samples in each cell average
83  DynamicList<label> cellAvgCount_;
84 
85  //- Number of samples in each point average
86  DynamicList<label> pointAvgCount_;
87 
88  //- Weight sums for each cell average
89  autoPtr<DynamicList<scalar>> cellAvgWeightSumPtr_;
90 
91  //- Weight sums for each point average
92  autoPtr<DynamicList<scalar>> pointAvgWeightSumPtr_;
93 
94  //- Value sums for each cell average
95  DynamicList<Type> cellAvgSum_;
96 
97  //- Value sums for each point average
98  DynamicList<Type> pointAvgSum_;
99 
100  //- Construct given a number of cells and points and a flag to
101  // specify whether or not weight sums need to be stored
102  data
103  (
104  const label nCells,
105  const label nPoints,
106  const bool hasWeightSum
107  );
108  };
109 
110 
111  // Private Data
112 
113  //- Constant cell weight sums. E.g., the volume of the cell basis
114  // function. Might be empty.
115  const autoPtr<Field<scalar>> cellWeightSumPtr_;
116 
117  //- Constant point weight sums. E.g., the volume of the point basis
118  // function. Might be empty.
119  const autoPtr<Field<scalar>> pointWeightSumPtr_;
120 
121  //- Base data
122  data data_;
123 
124  //- Flag to indicate whether or not difference data is available
125  bool dDataIsValid_;
126 
127  //- Difference data for correction later
128  data dData_;
129 
130 
131  // Private Member Functions
132 
133  //- Calculate a volume-weighted sum over the cells
134  template<class CellWeight>
135  static tmp<Field<scalar>> cellVolumeWeightedSum
136  (
137  const polyMesh& pMesh,
138  const CellWeight& cellWeight
139  );
140 
141  //- Calculate a volume-weighted sum over the points
142  static tmp<Field<scalar>> cellWeightSum
143  (
144  const polyMesh& pMesh,
145  const Field<scalar>& weightSum
146  );
147 
148  //- Calculate the cell constant weight sum
149  template<class CellWeight>
150  static tmp<Field<scalar>> pointVolumeWeightedSum
151  (
152  const polyMesh& pMesh,
153  const CellWeight& cellWeight
154  );
155 
156  //- Calculate the point constant weight sum
157  static tmp<Field<scalar>> pointWeightSum
158  (
159  const polyMesh& pMesh,
160  const Field<scalar>& weightSum
161  );
162 
163  //- Remove all values from the average data
164  static void clear(data& d);
165 
166  //- Remove values from the cell average data
167  static void removeFromCells
168  (
169  const LagrangianSubSubField<scalar>& weight,
171  data& d
172  );
173 
174  //- Add values to the cell average data
175  static void addToCells
176  (
177  const LagrangianSubSubField<scalar>& weight,
179  data& d
180  );
181 
182  //- Add values to the point average data
183  static void addToPoints
184  (
185  const LagrangianSubSubField<scalar>& weight,
187  data& d
188  );
189 
190  //- Remove values from the point average data
191  static void removeFromPoints(const data& dd, data& d);
192 
193  //- Add values to the point average data
194  static void addToPoints(const data& dd, data& d);
195 
196  //- Interpolate into a sub-field, where possible. Calling code deals
197  // with default values where no elements are available to interpolate.
198  virtual void interpolate(LagrangianSubField<Type>& result) const;
199 
200 
201 public:
202 
203  //- Runtime type information
204  TypeName("cellPoint");
205 
206 
207  // Constructors
208 
209  //- Construct with a name, for a mesh and with given dimensions
210  cellPoint
211  (
212  const word& name,
213  const LagrangianMesh& mesh,
214  const dimensionSet& dimensions,
215  const Field<scalar>& weightSum
216  );
217 
218 
219  //- Destructor
220  virtual ~cellPoint();
221 
222 
223  // Member Functions
224 
225  //- Remove weighted values from the average
226  virtual void remove
227  (
228  const LagrangianSubSubField<scalar>& weight,
230  );
231 
232  //- Add weighted values to the average
233  virtual void add
234  (
235  const LagrangianSubSubField<scalar>& weight,
237  const bool cache
238  );
239 
240  //- Correct weighted values in the average
241  virtual void correct
242  (
243  const LagrangianSubSubField<scalar>& weight,
245  const bool cache
246  );
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace LagrangianAverages
253 } // End namespace Foam
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 #ifdef NoRepository
259 #endif
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 #endif
264 
265 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated as...
virtual void remove(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi)
Remove weighted values from the average.
virtual void add(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)
Add weighted values to the average.
cellPoint(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions, const Field< scalar > &weightSum)
Construct with a name, for a mesh and with given dimensions.
virtual void correct(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)
Correct weighted values in the average.
TypeName("cellPoint")
Runtime type information.
Class containing Lagrangian geometry and topology.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Dimension set for the base types.
Definition: dimensionSet.H:125
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label nPoints
const volScalarField & psi
Namespace for OpenFOAM.
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.