cutPolyIntegralI.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 \*---------------------------------------------------------------------------*/
25 
26 #include "cutPolyIntegral.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace cutPoly
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 struct OpIndex
38 {
39  const label i_;
40 
41  inline OpIndex(const label i)
42  :
43  i_(i)
44  {}
45 
46  template<class Type>
47  inline const Type& operator()(const List<Type>& xs) const
48  {
49  return xs[i_];
50  }
51 };
52 
53 
54 struct OpBegin
55 {
56  template<class Type>
57  inline auto operator()(const Type& x) const
58  {
59  return x.begin();
60  }
61 };
62 
63 
65 {
66  template<class Type>
67  inline auto operator()(const Type& x) const
68  {
69  return *x;
70  }
71 };
72 
73 
74 struct OpNext
75 {
76  template<class Type>
77  inline auto operator()(const Type& x) const
78  {
79  return x.next();
80  }
81 };
82 
83 
84 template<class ScaleType>
85 struct OpScaled
86 {
87  const ScaleType s_;
88 
89  inline OpScaled(const ScaleType& s)
90  :
91  s_(s)
92  {}
93 
94  template<class Type>
95  inline auto operator()(const Type& x) const
96  {
97  return s_*x;
98  }
99 };
100 
101 
103 {
104  const vector& v_;
105 
106  inline OpPreInner(const vector& v)
107  :
108  v_(v)
109  {}
110 
111  template<class Type>
112  inline auto operator()(const Type& x) const
113  {
114  return v_ & x;
115  }
116 };
117 
118 
120 {
121  const labelUList& is_;
122 
123  inline OpIndirectAverage(const labelUList& is)
124  :
125  is_(is)
126  {}
127 
128  template<class Container>
129  inline auto operator()(const Container& xs) const
130  {
131  typename Container::value_type nResult =
133 
134  forAll(is_, i)
135  {
136  nResult += xs[is_[i]];
137  }
138 
139  return nResult/is_.size();
140  }
141 };
142 
143 
145 {
146  template<class Container>
147  inline auto operator()(const Container& xs) const
148  {
149  label n = 0;
150 
151  typename Container::value_type nResult =
153 
154  forAllConstIter(typename Container, xs, iter)
155  {
156  ++ n;
157  nResult += *iter;
158  }
159 
160  return nResult/n;
161  }
162 };
163 
164 
166 {
167  const face& f_;
170  const scalar isoAlpha_;
171  const bool below_;
172 
174  (
175  const face& f,
176  const List<labelPair>& fCuts,
177  const scalarField& pAlphas,
178  const scalar isoAlpha,
179  const bool below
180  )
181  :
182  f_(f),
183  fCuts_(fCuts),
184  pAlphas_(pAlphas),
185  isoAlpha_(isoAlpha),
186  below_(below)
187  {}
188 
189  template<class Type>
190  inline auto operator()(const Field<Type>& pPsis) const
191  {
192  return
194  (
195  f_,
196  fCuts_,
197  pPsis,
198  pAlphas_,
199  isoAlpha_,
200  below_
201  );
202  }
203 };
204 
205 
207 {
208  const cell& c_;
211  const faceList& fs_;
213  const scalar isoAlpha_;
214 
216  (
217  const cell& c,
218  const cellEdgeAddressing& cAddr,
219  const labelListList& cCuts,
220  const faceList& fs,
221  const scalarField& pAlphas,
222  const scalar isoAlpha
223  )
224  :
225  c_(c),
226  cAddr_(cAddr),
227  cCuts_(cCuts),
228  fs_(fs),
229  pAlphas_(pAlphas),
230  isoAlpha_(isoAlpha)
231  {}
232 
233  template<class Type>
234  inline auto operator()(const Field<Type>& pPsis) const
235  {
236  return
238  (
239  c_,
240  cAddr_,
241  cCuts_,
242  fs_,
243  pPsis,
244  pAlphas_,
245  isoAlpha_
246  );
247  }
248 };
249 
250 
252 {
253  template<class Type>
254  inline void operator()(Type& x) const
255  {
256  ++ x;
257  }
258 };
259 
260 
262 {
263  template<class Type>
264  inline auto operator()(const Type& a, const Type& b) const
265  {
266  return a + b;
267  }
268 };
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 } // End namespace cutPoly
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
278 (
279  const face& f,
280  const point& fPAvg,
281  const pointField& ps
282 )
283 {
284  return
286  (
287  FaceValues<point>(f, ps),
288  fPAvg,
289  std::make_tuple(),
290  std::make_tuple()
291  ).first();
292 }
293 
294 
296 (
297  const face& f,
298  const pointField& ps
299 )
300 {
301  return
302  faceArea
303  (
304  f,
305  OpIndirectAverage(f)(ps),
306  ps
307  );
308 }
309 
310 
311 template<class Type>
314 (
315  const face& f,
316  const point& fPAvg,
317  const Type& fPsiAvg,
318  const pointField& ps,
319  const Field<Type>& pPsis
320 )
321 {
322  auto result =
324  (
325  FaceValues<point>(f, ps),
326  fPAvg,
327  std::make_tuple(FaceValues<Type>(f, pPsis)),
328  std::make_tuple(fPsiAvg)
329  );
330 
331  return
333  (
334  result.first(),
335  std::get<0>(result.second())
336  );
337 }
338 
339 
340 template<class Type>
343 (
344  const face& f,
345  const pointField& ps,
346  const Field<Type>& pPsis
347 )
348 {
349  return
351  (
352  f,
353  OpIndirectAverage(f)(ps),
354  OpIndirectAverage(f)(pPsis),
355  ps,
356  pPsis
357  );
358 }
359 
360 
361 template<class Type>
363 (
364  const face& f,
365  const point& fPAvg,
366  const Type& fPsiAvg,
367  const pointField& ps,
368  const Field<Type>& pPsis
369 )
370 {
371  auto result =
373  (
374  FaceValues<point>(f, ps),
375  fPAvg,
376  std::make_tuple(FaceValues<Type>(f, pPsis)),
377  std::make_tuple(fPsiAvg)
378  );
379 
380  return
381  Tuple2<vector, Type>
382  (
383  result.first(),
384  std::get<0>(result.second())
385  );
386 }
387 
388 
389 template<class Type>
391 (
392  const face& f,
393  const pointField& ps,
394  const Field<Type>& pPsis
395 )
396 {
397  return
399  (
400  f,
401  OpIndirectAverage(f)(ps),
402  OpIndirectAverage(f)(pPsis),
403  ps,
404  pPsis
405  );
406 }
407 
408 
410 (
411  const face& f,
412  const vector& fArea,
413  const List<labelPair>& fCuts,
414  const pointField& ps,
415  const scalarField& pAlphas,
416  const scalar isoAlpha,
417  const bool below
418 )
419 {
420  return
422  (
423  f,
424  fArea,
425  std::make_tuple(),
426  fCuts,
427  ps,
428  std::make_tuple(),
429  pAlphas,
430  isoAlpha,
431  below
432  ).first();
433 }
434 
435 
436 template<class Type>
439 (
440  const face& f,
441  const vector& fArea,
442  const Type& fPsi,
443  const List<labelPair>& fCuts,
444  const pointField& ps,
445  const Field<Type>& pPsis,
446  const scalarField& pAlphas,
447  const scalar isoAlpha,
448  const bool below
449 )
450 {
451  auto result =
452  faceCutAreaIntegral<Type>
453  (
454  f,
455  fArea,
456  std::make_tuple(fPsi),
457  fCuts,
458  ps,
459  std::forward_as_tuple(pPsis),
460  pAlphas,
461  isoAlpha,
462  below
463  );
464 
465  return
467  (
468  result.first(),
469  std::get<0>(result.second())
470  );
471 }
472 
473 
474 inline Foam::scalar Foam::cutPoly::cellVolume
475 (
476  const cell& c,
477  const cellEdgeAddressing& cAddr,
478  const point& cPAvg,
479  const vectorField& fAreas,
480  const pointField& fCentres
481 )
482 {
483  return
485  (
486  c,
487  cAddr,
488  cPAvg,
489  std::make_tuple(),
490  fAreas,
491  fCentres,
492  std::make_tuple()
493  ).first();
494 }
495 
496 
497 inline Foam::scalar Foam::cutPoly::cellVolume
498 (
499  const cell& c,
500  const cellEdgeAddressing& cAddr,
501  const vectorField& fAreas,
502  const pointField& fCentres
503 )
504 {
505  return
506  cellVolume
507  (
508  c,
509  cAddr,
510  OpIndirectAverage(c)(fCentres),
511  fAreas,
512  fCentres
513  );
514 }
515 
516 
517 template<class Type>
519 (
520  const cell& c,
521  const cellEdgeAddressing& cAddr,
522  const point& cPAvg,
523  const Type& cPsiAvg,
524  const vectorField& fAreas,
525  const pointField& fCentres,
526  const Field<Type>& fPsis
527 )
528 {
529  auto result =
530  cellVolumeIntegral<Type>
531  (
532  c,
533  cAddr,
534  cPAvg,
535  std::make_tuple(cPsiAvg),
536  fAreas,
537  fCentres,
538  std::forward_as_tuple(fPsis)
539  );
540 
541  return
543  (
544  result.first(),
545  std::get<0>(result.second())
546  );
547 }
548 
549 
550 template<class Type>
552 (
553  const cell& c,
554  const cellEdgeAddressing& cAddr,
555  const vectorField& fAreas,
556  const pointField& fCentres,
557  const Field<Type>& fPsis
558 )
559 {
560  return
562  (
563  c,
564  cAddr,
565  OpIndirectAverage(c)(fCentres),
566  OpIndirectAverage(c)(fPsis),
567  fAreas,
568  fCentres,
569  fPsis
570  );
571 }
572 
573 
575 (
576  const cell& c,
577  const cellEdgeAddressing& cAddr,
578  const scalar cVolume,
579  const labelListList& cCuts,
580  const faceUList& fs,
581  const vectorField& fAreas,
582  const pointField& fCentres,
583  const vectorField& fCutAreas,
584  const pointField& ps,
585  const scalarField& pAlphas,
586  const scalar isoAlpha,
587  const bool below
588 )
589 {
590  return
592  (
593  c,
594  cAddr,
595  cVolume,
596  std::make_tuple(),
597  cCuts,
598  fs,
599  fAreas,
600  fCentres,
601  std::make_tuple(),
602  fCutAreas,
603  std::make_tuple(),
604  ps,
605  std::make_tuple(),
606  pAlphas,
607  isoAlpha,
608  below
609  ).first();
610 }
611 
612 
613 template<class Type>
615 (
616  const cell& c,
617  const cellEdgeAddressing& cAddr,
618  const scalar cVolume,
619  const Type& cPsi,
620  const labelListList& cCuts,
621  const faceUList& fs,
622  const vectorField& fAreas,
623  const pointField& fCentres,
624  const Field<Type>& fPsis,
625  const vectorField& fCutAreas,
626  const Field<Type>& fCutPsis,
627  const pointField& ps,
628  const Field<Type>& pPsis,
629  const scalarField& pAlphas,
630  const scalar isoAlpha,
631  const bool below
632 )
633 {
634  auto result =
636  (
637  c,
638  cAddr,
639  cVolume,
640  std::make_tuple(cPsi),
641  cCuts,
642  fs,
643  fAreas,
644  fCentres,
645  std::forward_as_tuple(fPsis),
646  fCutAreas,
647  std::forward_as_tuple(fCutPsis),
648  ps,
649  std::forward_as_tuple(pPsis),
650  pAlphas,
651  isoAlpha,
652  below
653  );
654 
655  return Tuple2<scalar, Type>(result.first(), std::get<0>(result.second()));
656 }
657 
658 
659 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
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
Traits class for primitives.
Definition: pTraits.H:53
volVectorField vectorField(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
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.
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.
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
vector point
Point is a vector.
Definition: point.H:41
labelList f(nPoints)
auto operator()(const Type &a, const Type &b) const
auto operator()(const Type &x) const
auto operator()(const Field< Type > &pPsis) const
OpCellCutValues(const cell &c, const cellEdgeAddressing &cAddr, const labelListList &cCuts, const faceList &fs, const scalarField &pAlphas, const scalar isoAlpha)
const labelListList & cCuts_
const cellEdgeAddressing & cAddr_
auto operator()(const Type &x) const
auto operator()(const Field< Type > &pPsis) const
const List< labelPair > & fCuts_
OpFaceCutValues(const face &f, const List< labelPair > &fCuts, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
const Type & operator()(const List< Type > &xs) const
OpIndirectAverage(const labelUList &is)
auto operator()(const Container &xs) const
auto operator()(const Container &xs) const
auto operator()(const Type &x) const
auto operator()(const Type &x) const
OpScaled(const ScaleType &s)
auto operator()(const Type &x) const