unintegrableForNonZeroQ.C
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) 2023-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 \*---------------------------------------------------------------------------*/
25 
27 
28 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
29 
31 (
32  const label q,
33  const scalarField& x
34 ) const
35 {
36  if (q == 0)
37  {
38  return PhiForZeroQ(x);
39  }
40  else
41  {
42  return unintegrable::Phi(q, x);
43  }
44 }
45 
46 
48 (
49  const label q
50 ) const
51 {
52  if (q == 0)
53  {
54  const scalarField x(scalarList({min(), max()}));
55  const scalarField Phi(this->PhiForZeroQ(x));
56  return Pair<scalar>(Phi.first(), Phi.last());
57  }
58  else
59  {
60  return unintegrable::Phi01(q);
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
68 {
69  if (q() == 0)
70  {
71  return sampleForZeroQ();
72  }
73  else
74  {
75  return unintegrable::sample();
76  }
77 }
78 
79 
82 (
83  const scalarField& x,
84  const label e,
85  const bool consistent
86 ) const
87 {
88  if (!consistent && q() + e == 0)
89  {
90  const Pair<scalar>& Phi01 = this->Phi01();
91  const scalarField xClip(Foam::min(Foam::max(x, min()), max()));
92  return (PhiForZeroQ(xClip) - Phi01[0])/(Phi01[1] - Phi01[0]);
93  }
94  else
95  {
97  }
98 }
99 
100 
101 // ************************************************************************* //
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
T & first()
Return the first element of the list.
Definition: UListI.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:128
label q() const
Return the effective distribution size exponent.
Definition: distribution.H:107
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const
Return the integral of the PDF multiplied by an integer power of x.
virtual scalar sample() const
Sample the distribution.
virtual Pair< scalar > Phi01(const label q) const
Access cached values of the un-normalised CDF at the minimum and.
Definition: unintegrable.C:330
virtual tmp< scalarField > PhiForZeroQ(const scalarField &x) const =0
Return values of the un-normalised CDF for zero effective size.
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: unintegrable.C:320
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const
Return the integral of the PDF multiplied by an integer power of x.
Definition: unintegrable.C:445
virtual scalar sample() const
Sample the distribution.
Definition: unintegrable.C:416
const Pair< scalar > & Phi01() const
Access cached values of the un-normalised CDF at the minimum and.
Definition: unintegrable.C:340
A class for managing temporary objects.
Definition: tmp.H:55
const doubleScalar e
Definition: doubleScalar.H:106
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
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)