xref: /petsc/src/ts/utils/dmplexlandau/land_tensors.h (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 #define LANDAU_INVSQRT(q) (1. / PetscSqrtReal(q))
2 
3 #if defined(__CUDA_ARCH__)
4 #define PETSC_DEVICE_FUNC_DECL __device__
5 #elif defined(KOKKOS_INLINE_FUNCTION)
6 #define PETSC_DEVICE_FUNC_DECL KOKKOS_INLINE_FUNCTION
7 #else
8 #define PETSC_DEVICE_FUNC_DECL static
9 #endif
10 
11 #if LANDAU_DIM == 2
12 /* elliptic functions
13  */
14 PETSC_DEVICE_FUNC_DECL PetscReal polevl_10(PetscReal x, const PetscReal coef[]) {
15   PetscReal ans;
16   PetscInt  i;
17   ans = coef[0];
18   for (i = 1; i < 11; i++) ans = ans * x + coef[i];
19   return (ans);
20 }
21 PETSC_DEVICE_FUNC_DECL PetscReal polevl_9(PetscReal x, const PetscReal coef[]) {
22   PetscReal ans;
23   PetscInt  i;
24   ans = coef[0];
25   for (i = 1; i < 10; i++) ans = ans * x + coef[i];
26   return (ans);
27 }
28 /*
29  *      Complete elliptic integral of the second kind
30  */
31 PETSC_DEVICE_FUNC_DECL void ellipticE(PetscReal x, PetscReal *ret) {
32 #if defined(PETSC_USE_REAL_SINGLE)
33   static const PetscReal P2[] = {1.53552577301013293365E-4F, 2.50888492163602060990E-3F, 8.68786816565889628429E-3F, 1.07350949056076193403E-2F, 7.77395492516787092951E-3F, 7.58395289413514708519E-3F,
34                                  1.15688436810574127319E-2F, 2.18317996015557253103E-2F, 5.68051945617860553470E-2F, 4.43147180560990850618E-1F, 1.00000000000000000299E0F};
35   static const PetscReal Q2[] = {3.27954898576485872656E-5F, 1.00962792679356715133E-3F, 6.50609489976927491433E-3F, 1.68862163993311317300E-2F, 2.61769742454493659583E-2F,
36                                  3.34833904888224918614E-2F, 4.27180926518931511717E-2F, 5.85936634471101055642E-2F, 9.37499997197644278445E-2F, 2.49999999999888314361E-1F};
37 #else
38   static const PetscReal P2[] = {1.53552577301013293365E-4, 2.50888492163602060990E-3, 8.68786816565889628429E-3, 1.07350949056076193403E-2, 7.77395492516787092951E-3, 7.58395289413514708519E-3,
39                                  1.15688436810574127319E-2, 2.18317996015557253103E-2, 5.68051945617860553470E-2, 4.43147180560990850618E-1, 1.00000000000000000299E0};
40   static const PetscReal Q2[] = {3.27954898576485872656E-5, 1.00962792679356715133E-3, 6.50609489976927491433E-3, 1.68862163993311317300E-2, 2.61769742454493659583E-2,
41                                  3.34833904888224918614E-2, 4.27180926518931511717E-2, 5.85936634471101055642E-2, 9.37499997197644278445E-2, 2.49999999999888314361E-1};
42 #endif
43   x    = 1 - x; /* where m = 1 - m1 */
44   *ret = polevl_10(x, P2) - PetscLogReal(x) * (x * polevl_9(x, Q2));
45 }
46 /*
47  *      Complete elliptic integral of the first kind
48  */
49 PETSC_DEVICE_FUNC_DECL void ellipticK(PetscReal x, PetscReal *ret) {
50 #if defined(PETSC_USE_REAL_SINGLE)
51   static const PetscReal P1[] = {1.37982864606273237150E-4F, 2.28025724005875567385E-3F, 7.97404013220415179367E-3F, 9.85821379021226008714E-3F, 6.87489687449949877925E-3F, 6.18901033637687613229E-3F,
52                                  8.79078273952743772254E-3F, 1.49380448916805252718E-2F, 3.08851465246711995998E-2F, 9.65735902811690126535E-2F, 1.38629436111989062502E0F};
53   static const PetscReal Q1[] = {2.94078955048598507511E-5F, 9.14184723865917226571E-4F, 5.94058303753167793257E-3F, 1.54850516649762399335E-2F, 2.39089602715924892727E-2F, 3.01204715227604046988E-2F,
54                                  3.73774314173823228969E-2F, 4.88280347570998239232E-2F, 7.03124996963957469739E-2F, 1.24999999999870820058E-1F, 4.99999999999999999821E-1F};
55 #else
56   static const PetscReal P1[] = {1.37982864606273237150E-4, 2.28025724005875567385E-3, 7.97404013220415179367E-3, 9.85821379021226008714E-3, 6.87489687449949877925E-3, 6.18901033637687613229E-3,
57                                  8.79078273952743772254E-3, 1.49380448916805252718E-2, 3.08851465246711995998E-2, 9.65735902811690126535E-2, 1.38629436111989062502E0};
58   static const PetscReal Q1[] = {2.94078955048598507511E-5, 9.14184723865917226571E-4, 5.94058303753167793257E-3, 1.54850516649762399335E-2, 2.39089602715924892727E-2, 3.01204715227604046988E-2,
59                                  3.73774314173823228969E-2, 4.88280347570998239232E-2, 7.03124996963957469739E-2, 1.24999999999870820058E-1, 4.99999999999999999821E-1};
60 #endif
61   x    = 1 - x; /* where m = 1 - m1 */
62   *ret = polevl_10(x, P1) - PetscLogReal(x) * polevl_10(x, Q1);
63 }
64 /* flip sign. papers use du/dt = C, PETSc uses form G(u) = du/dt - C(u) = 0 */
65 PETSC_DEVICE_FUNC_DECL void LandauTensor2D(const PetscReal x[], const PetscReal rp, const PetscReal zp, PetscReal Ud[][2], PetscReal Uk[][2], const PetscReal mask) {
66   PetscReal l, s, r = x[0], z = x[1], i1func, i2func, i3func, ks, es, pi4pow, sqrt_1s, r2, rp2, r2prp2, zmzp, zmzp2, tt;
67   //PetscReal mask /* = !!(r!=rp || z!=zp) */;
68   /* !!(zmzp2 > 1.e-12 || (r-rp) >  1.e-12 || (r-rp) < -1.e-12); */
69   r2      = PetscSqr(r);
70   zmzp    = z - zp;
71   rp2     = PetscSqr(rp);
72   zmzp2   = PetscSqr(zmzp);
73   r2prp2  = r2 + rp2;
74   l       = r2 + rp2 + zmzp2;
75   /* if      (zmzp2 >  PETSC_SMALL) mask = 1; */
76   /* else if ((tt=(r-rp)) >  PETSC_SMALL) mask = 1; */
77   /* else if  (tt         < -PETSC_SMALL) mask = 1; */
78   /* else mask = 0; */
79   s       = mask * 2 * r * rp / l; /* mask for vectorization */
80   tt      = 1. / (1 + s);
81   pi4pow  = 4 * PETSC_PI * LANDAU_INVSQRT(PetscSqr(l) * l);
82   sqrt_1s = PetscSqrtReal(1. + s);
83   /* sp.ellipe(2.*s/(1.+s)) */
84   ellipticE(2 * s * tt, &es); /* 44 flops * 2 + 75 = 163 flops including 2 logs, 1 sqrt, 1 pow, 21 mult */
85   /* sp.ellipk(2.*s/(1.+s)) */
86   ellipticK(2 * s * tt, &ks); /* 44 flops + 75 in rest, 21 mult */
87   /* mask is needed here just for single precision */
88   i2func   = 2. / ((1 - s) * sqrt_1s) * es;
89   i1func   = 4. / (PetscSqr(s) * sqrt_1s + PETSC_MACHINE_EPSILON) * mask * (ks - (1. + s) * es);
90   i3func   = 2. / ((1 - s) * (s)*sqrt_1s + PETSC_MACHINE_EPSILON) * (es - (1 - s) * ks);
91   Ud[0][0] = -pi4pow * (rp2 * i1func + PetscSqr(zmzp) * i2func);
92   Ud[0][1] = Ud[1][0] = Uk[0][1] = pi4pow * (zmzp) * (r * i2func - rp * i3func);
93   Uk[1][1] = Ud[1][1] = -pi4pow * ((r2prp2)*i2func - 2 * r * rp * i3func) * mask;
94   Uk[0][0]            = -pi4pow * (zmzp2 * i3func + r * rp * i1func);
95   Uk[1][0]            = pi4pow * (zmzp) * (r * i3func - rp * i2func); /* 48 mults + 21 + 21 = 90 mults and divs */
96 }
97 #else
98 /* integration point functions */
99 /* Evaluates the tensor U=(I-(x-y)(x-y)/(x-y)^2)/|x-y| at point x,y */
100 /* if x==y we will return zero. This is not the correct result */
101 /* since the tensor diverges for x==y but when integrated */
102 /* the divergent part is antisymmetric and vanishes. This is not  */
103 /* trivial, but can be proven. */
104 PETSC_DEVICE_FUNC_DECL void LandauTensor3D(const PetscReal x1[], const PetscReal xp, const PetscReal yp, const PetscReal zp, PetscReal U[][3], PetscReal mask) {
105   PetscReal dx[3], inorm3, inorm, inorm2, norm2, x2[] = {xp, yp, zp};
106   PetscInt  d;
107   for (d = 0, norm2 = PETSC_MACHINE_EPSILON; d < 3; ++d) {
108     dx[d] = x2[d] - x1[d];
109     norm2 += dx[d] * dx[d];
110   }
111   inorm2 = mask / norm2;
112   inorm  = PetscSqrtReal(inorm2);
113   inorm3 = inorm2 * inorm;
114   for (d = 0; d < 3; ++d) U[d][d] = -(inorm - inorm3 * dx[d] * dx[d]);
115   U[1][0] = U[0][1] = inorm3 * dx[0] * dx[1];
116   U[1][2] = U[2][1] = inorm3 * dx[2] * dx[1];
117   U[2][0] = U[0][2] = inorm3 * dx[0] * dx[2];
118 }
119 /* Relativistic form */
120 #define GAMMA3(_x, _c02) PetscSqrtReal(1.0 + ((_x[0] * _x[0]) + (_x[1] * _x[1]) + (_x[2] * _x[2])) / (_c02))
121 PETSC_DEVICE_FUNC_DECL void LandauTensor3DRelativistic(const PetscReal a_x1[], const PetscReal xp, const PetscReal yp, const PetscReal zp, PetscReal U[][3], PetscReal mask, PetscReal c0) {
122   const PetscReal x2[3] = {xp, yp, zp}, x1[3] = {a_x1[0], a_x1[1], a_x1[2]}, c02 = c0 * c0, g1 = GAMMA3(x1, c02), g2 = GAMMA3(x2, c02);
123   PetscReal       fact, u1u2, diff[3], udiff2, u12, u22, wsq, rsq, tt;
124   PetscInt        i, j;
125 
126   if (mask == 0.0) {
127     for (i = 0; i < 3; ++i) {
128       for (j = 0; j < 3; ++j) { U[i][j] = 0; }
129     }
130   } else {
131     for (i = 0, u1u2 = u12 = u22 = udiff2 = 0; i < 3; ++i) {
132       diff[i] = x1[i] - x2[i];
133       udiff2 += diff[i] * diff[i];
134       u12 += x1[i] * x1[i];
135       u22 += x2[i] * x2[i];
136       u1u2 += x1[i] * x2[i];
137     }
138     tt   = 2. * u1u2 * (1. - g1 * g2) + (u12 * u22 + u1u2 * u1u2) / c02; // these two terms are about the same with opposite sign
139     wsq  = udiff2 + tt;
140     //wsq = udiff2 + 2.*u1u2*(1.-g1*g2) + (u12*u22 + u1u2*u1u2)/c02;
141     rsq  = 1. + wsq / c02;
142     fact = -rsq / (g1 * g2 * PetscSqrtReal(wsq)); /* flip sign. papers use du/dt = C, PETSc uses form G(u) = du/dt - C(u) = 0 */
143     for (i = 0; i < 3; ++i) {
144       for (j = 0; j < 3; ++j) { U[i][j] = fact * (-diff[i] * diff[j] / wsq + (PetscSqrtReal(rsq) - 1.) * (x1[i] * x2[j] + x1[j] * x2[i]) / wsq); }
145       U[i][i] += fact;
146     }
147   }
148 }
149 
150 #endif
151