cell_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  Basic Lagrangian averaging process in which values are averaged in the
29  cells and then interpolated assuming a constant value within the cell
30 
31 SourceFiles
32  cell_LagrangianAverage.C
33  cell_LagrangianAverages.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef cell_LagrangianAverage_H
38 #define cell_LagrangianAverage_H
39 
40 #include "LagrangianAverage.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace LagrangianAverages
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class LagrangianAverage Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class Type>
54 class cell
55 :
56  public LagrangianAverage<Type>
57 {
58 private:
59 
60  // Private Classes
61 
62  //- Data structure for an average
63  struct data
64  {
65  //- Map from the cell to the cell average
66  List<label> cellCellAvg_;
67 
68  //- Map from the cell average to the cell
69  DynamicList<label> cellAvgCell_;
70 
71  //- Number of samples in each cell average
72  DynamicList<label> cellAvgCount_;
73 
74  //- Weight sums for each cell average
75  autoPtr<DynamicList<scalar>> cellAvgWeightSumPtr_;
76 
77  //- Value sums for each cell average
78  DynamicList<Type> cellAvgSum_;
79 
80  //- Construct given a number of cells and a flag to specify whether
81  // or not weight sums need to be stored
82  data(const label nCells, const bool hasWeightSum);
83  };
84 
85 
86  // Private Data
87 
88  //- Constant cell weight sums. E.g., cell volume. Might be null.
89  const Field<scalar>& weightSum_;
90 
91  //- Base data
92  data data_;
93 
94  //- Flag to indicate whether or not difference data is available
95  bool dDataIsValid_;
96 
97  //- Difference data for correction later
98  data dData_;
99 
100 
101  // Private Member Functions
102 
103  //- Remove all values from the average data
104  static void clear(data& d);
105 
106  //- Remove values from the average data
107  static void remove
108  (
109  const LagrangianSubSubField<scalar>& weight,
111  data& d
112  );
113 
114  //- Add values to the average data
115  static void add
116  (
117  const LagrangianSubSubField<scalar>& weight,
119  data& d
120  );
121 
122  //- Interpolate into a sub-field, where possible. Calling code deals
123  // with default values where no elements are available to interpolate.
124  virtual void interpolate(LagrangianSubField<Type>& result) const;
125 
126 
127 public:
128 
129  //- Runtime type information
130  TypeName("cell");
131 
132 
133  // Constructors
134 
135  //- Construct with a name, for a mesh and with given dimensions
136  cell
137  (
138  const word& name,
139  const LagrangianMesh& mesh,
140  const dimensionSet& dimensions,
141  const Field<scalar>& weightSum
142  );
143 
144 
145  //- Destructor
146  virtual ~cell();
147 
148 
149  // Member Functions
150 
151  //- Remove weighted values from the average
152  virtual void remove
153  (
154  const LagrangianSubSubField<scalar>& weight,
156  );
157 
158  //- Add weighted values to the average
159  virtual void add
160  (
161  const LagrangianSubSubField<scalar>& weight,
163  const bool cache
164  );
165 
166  //- Correct weighted values in the average
167  virtual void correct
168  (
169  const LagrangianSubSubField<scalar>& weight,
171  const bool cache
172  );
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace LagrangianAverages
179 } // End namespace Foam
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 #ifdef NoRepository
184  #include "cell_LagrangianAverage.C"
185 #endif
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
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...
cell(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions, const Field< scalar > &weightSum)
Construct with a name, for a mesh and with given dimensions.
TypeName("cell")
Runtime type information.
virtual void correct(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)
Correct weighted values in the average.
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
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)
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.
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)