erfInv.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "MathFunctions.H"
30 #include "error.H"
31 #include <cmath>
32 #include <limits>
33 
34 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
35 
36 Foam::scalar Foam::Math::erfInv(const scalar y)
37 {
38  #ifdef FULLDEBUG
39  if (mag(y) >= scalar(1))
40  {
42  << "The domain of inverse error function argument "
43  << "(i.e. y) should be limited to (-1, 1):" << nl
44  << " y = " << y
45  << endl;
46 
47  return std::numeric_limits<scalar>::infinity();
48  }
49  #endif
50 
51  // (W:p. 2) to reduce the max relative error to O(1e-4)
52  constexpr scalar a = 0.147;
53 
54  const scalar k =
55  scalar(2)/(a*constant::mathematical::pi) + 0.5*log(scalar(1) - sqr(y));
56 
57  const scalar h = log(scalar(1) - sqr(y))/a;
58 
59  // (W:Eq. 7)
60  const scalar x = sqrt(-k + sqrt(sqr(k) - h));
61 
62  if (y < 0)
63  {
64  return -x;
65  }
66 
67  return x;
68 }
69 
70 
71 // ************************************************************************* //
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
label k
Boltzmann constant.
scalar erfInv(const scalar y)
Inverse error function of a real-number argument.
Definition: erfInv.C:29
scalar y
constexpr scalar pi(M_PI)
const dimensionedScalar h
Planck constant.
#define WarningInFunction
Report a warning using Foam::Warning.