cutPolyIntegral.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) 2022-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 Namespace
25  Foam::cutPoly
26 
27 Description
28  Functions for computing integrals over cut faces and cells
29 
30 SourceFiles
31  cutPolyIntegralTemplates.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef cutPolyIntegral_H
36 #define cutPolyIntegral_H
37 
38 #include "cutPolyValue.H"
39 #include <tuple>
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace cutPoly
46 {
47 
48 //- Definition for the type of a face-area-integral
49 template<class Type>
52 
53 //- Compute the face-area and face-area-integrals of the given properties
54 template<template<class> class FaceValues, class ... Types>
56 (
57  const FaceValues<point>& fPs,
58  const point& fPAvg,
59  const std::tuple<FaceValues<Types> ...>& fPsis,
60  const std::tuple<Types ...>& fPsiAvgs
61 );
62 
63 //- Compute the face-area
64 inline vector faceArea
65 (
66  const face& f,
67  const point& fPAvg,
68  const pointField& ps
69 );
70 
71 //- Compute the face-area
72 inline vector faceArea
73 (
74  const face& f,
75  const pointField& ps
76 );
77 
78 //- Compute the face-area and face-area-integral of a given property
79 template<class Type>
81 (
82  const face& f,
83  const point& fPAvg,
84  const Type& fPsiAvg,
85  const pointField& ps,
86  const Field<Type>& pPsis
87 );
88 
89 //- Compute the face-area and face-area-integral of a given property
90 template<class Type>
92 (
93  const face& f,
94  const pointField& ps,
95  const Field<Type>& pPsis
96 );
97 
98 //- Compute the face-area and face-area-averages of the given properties
99 template<template<class> class FaceValues, class ... Types>
100 Tuple2<vector, std::tuple<Types ...>> faceAreaAverage
101 (
102  const FaceValues<point>& fPs,
103  const point& fPAvg,
104  const std::tuple<FaceValues<Types> ...>& fPsis,
105  const std::tuple<Types ...>& fPsiAvgs
106 );
107 
108 //- Compute the face-area and face-area-average of a given property
109 template<class Type>
111 (
112  const face& f,
113  const point& fPAvg,
114  const Type& fPsiAvg,
115  const pointField& ps,
116  const Field<Type>& pPsis
117 );
118 
119 //- Compute the face-area and face-area-average of a given property
120 template<class Type>
122 (
123  const face& f,
124  const pointField& ps,
125  const Field<Type>& pPsis
126 );
127 
128 //- Compute the face-cut-area and face-cut-area-integral of the given properties
129 template<class ... Types>
131 (
132  const face& f,
133  const vector& fArea,
134  const std::tuple<Types ...>& fPsis,
135  const List<labelPair>& fCuts,
136  const pointField& ps,
137  const std::tuple<const Field<Types>& ...>& pPsis,
138  const scalarField& pAlphas,
139  const scalar isoAlpha,
140  const bool below
141 );
142 
143 //- Compute the face-cut-area
144 inline vector faceCutArea
145 (
146  const face& f,
147  const vector& fArea,
148  const List<labelPair>& fCuts,
149  const pointField& ps,
150  const scalarField& pAlphas,
151  const scalar isoAlpha,
152  const bool below
153 );
154 
155 //- Compute the face-cut-area and face-cut-area-integral of the given property
156 template<class Type>
158 (
159  const face& f,
160  const vector& fArea,
161  const Type& fPsi,
162  const List<labelPair>& fCuts,
163  const pointField& ps,
164  const Field<Type>& pPsis,
165  const scalarField& pAlphas,
166  const scalar isoAlpha,
167  const bool below
168 );
169 
170 //- Compute the cell-volume and cell-volume-integrals of the given properties
171 template<class ... Types>
172 Tuple2<scalar, std::tuple<Types ...>> cellVolumeIntegral
173 (
174  const cell& c,
175  const cellEdgeAddressing& cAddr,
176  const point& cPAvg,
177  const std::tuple<Types ...>& cPsiAvgs,
178  const vectorField& fAreas,
179  const pointField& fCentres,
180  const std::tuple<const Field<Types>& ...>& fPsis
181 );
182 
183 //- Compute the cell-volume
184 inline scalar cellVolume
185 (
186  const cell& c,
187  const cellEdgeAddressing& cAddr,
188  const point& cPAvg,
189  const vectorField& fAreas,
190  const pointField& fCentres
191 );
192 
193 //- Compute the cell-volume
194 inline scalar cellVolume
195 (
196  const cell& c,
197  const cellEdgeAddressing& cAddr,
198  const vectorField& fAreas,
199  const pointField& fCentres
200 );
201 
202 //- Compute the cell-volume and cell-volume-integrals of the given property
203 template<class Type>
205 (
206  const cell& c,
207  const cellEdgeAddressing& cAddr,
208  const point& cPAvg,
209  const Type& cPsiAvg,
210  const vectorField& fAreas,
211  const pointField& fCentres,
212  const Field<Type>& fPsis
213 );
214 
215 //- Compute the cell-volume and cell-volume-integrals of the given property
216 template<class Type>
218 (
219  const cell& c,
220  const cellEdgeAddressing& cAddr,
221  const vectorField& fAreas,
222  const pointField& fCentres,
223  const Field<Type>& fPsis
224 );
225 
226 //- Compute the cell-cut-volume and cell-cut-volume-integral
227 // of the given properties
228 template<class ... Types>
229 Tuple2<scalar, std::tuple<Types ...>> cellCutVolumeIntegral
230 (
231  const cell& c,
232  const cellEdgeAddressing& cAddr,
233  const scalar cVolume,
234  const std::tuple<Types ...>& cPsis,
235  const labelListList& cCuts,
236  const faceUList& fs,
237  const vectorField& fAreas,
238  const pointField& fCentres,
239  const std::tuple<const Field<Types>& ...>& fPsis,
240  const vectorField& fCutAreas,
241  const std::tuple<const Field<Types>& ...>& fCutPsis,
242  const pointField& ps,
243  const std::tuple<const Field<Types>& ...>& pPsis,
244  const scalarField& pAlphas,
245  const scalar isoAlpha,
246  const bool below
247 );
248 
249 //- Compute the cell-cut-volume
250 inline scalar cellCutVolume
251 (
252  const cell& c,
253  const cellEdgeAddressing& cAddr,
254  const scalar cVolume,
255  const labelListList& cCuts,
256  const faceUList& fs,
257  const vectorField& fAreas,
258  const pointField& fCentres,
259  const vectorField& fCutAreas,
260  const pointField& ps,
261  const scalarField& pAlphas,
262  const scalar isoAlpha,
263  const bool below
264 );
265 
266 //- Compute the cell-cut-volume and cell-cut-volume-integral
267 // of the given property
268 template<class Type>
270 (
271  const cell& c,
272  const cellEdgeAddressing& cAddr,
273  const scalar cVolume,
274  const Type& cPsi,
275  const labelListList& cCuts,
276  const faceUList& fs,
277  const vectorField& fAreas,
278  const pointField& fCentres,
279  const Field<Type>& fPsis,
280  const vectorField& fCutAreas,
281  const Field<Type>& fCutPsis,
282  const pointField& ps,
283  const Field<Type>& pPsis,
284  const scalarField& pAlphas,
285  const scalar isoAlpha,
286  const bool below
287 );
288 
289 } // End namespace cutPoly
290 } // End namespace Foam
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 #include "cutPolyIntegralI.H"
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 #ifdef NoRepository
299  #include "cutPolyIntegralTemplates.C"
300 #endif
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #endif
305 
306 // ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
const dimensionedScalar c
Speed of light in a vacuum.
Tuple2< vector, std::tuple< AreaIntegralType< Types > ... > > faceAreaIntegral(const FaceValues< point > &fPs, const point &fPAvg, const std::tuple< FaceValues< Types > ... > &fPsis, const std::tuple< Types ... > &fPsiAvgs)
Compute the face-area and face-area-integrals of the given properties.
Tuple2< vector, std::tuple< AreaIntegralType< Types > ... > > faceCutAreaIntegral(const face &f, const vector &fArea, const std::tuple< Types ... > &fPsis, const List< labelPair > &fCuts, const pointField &ps, const std::tuple< const Field< Types > &... > &pPsis, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the face-cut-area and face-cut-area-integral of the given properties.
scalar cellVolume(const cell &c, const cellEdgeAddressing &cAddr, const point &cPAvg, const vectorField &fAreas, const pointField &fCentres)
Compute the cell-volume.
vector faceCutArea(const face &f, const vector &fArea, const List< labelPair > &fCuts, const pointField &ps, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the face-cut-area.
Tuple2< scalar, std::tuple< Types ... > > cellCutVolumeIntegral(const cell &c, const cellEdgeAddressing &cAddr, const scalar cVolume, const std::tuple< Types ... > &cPsis, const labelListList &cCuts, const faceUList &fs, const vectorField &fAreas, const pointField &fCentres, const std::tuple< const Field< Types > &... > &fPsis, const vectorField &fCutAreas, const std::tuple< const Field< Types > &... > &fCutPsis, const pointField &ps, const std::tuple< const Field< Types > &... > &pPsis, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the cell-cut-volume and cell-cut-volume-integral.
Tuple2< scalar, std::tuple< Types ... > > cellVolumeIntegral(const cell &c, const cellEdgeAddressing &cAddr, const point &cPAvg, const std::tuple< Types ... > &cPsiAvgs, const vectorField &fAreas, const pointField &fCentres, const std::tuple< const Field< Types > &... > &fPsis)
Compute the cell-volume and cell-volume-integrals of the given properties.
typename Foam::outerProduct< Foam::vector, Type >::type AreaIntegralType
Definition for the type of a face-area-integral.
Tuple2< vector, std::tuple< Types ... > > faceAreaAverage(const FaceValues< point > &fPs, const point &fPAvg, const std::tuple< FaceValues< Types > ... > &fPsis, const std::tuple< Types ... > &fPsiAvgs)
Compute the face-area and face-area-averages of the given properties.
vector faceArea(const face &f, const point &fPAvg, const pointField &ps)
Compute the face-area.
scalar cellCutVolume(const cell &c, const cellEdgeAddressing &cAddr, const scalar cVolume, const labelListList &cCuts, const faceUList &fs, const vectorField &fAreas, const pointField &fCentres, const vectorField &fCutAreas, const pointField &ps, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the cell-cut-volume.
Namespace for OpenFOAM.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
labelList f(nPoints)