1 #define LANDAU_INVSQRT(q) (1./PetscSqrtReal(q)) 2 #define LANDAU_SQRT(q) PetscSqrtReal(q) 3 4 #if LANDAU_DIM==2 5 /* elliptic functions 6 */ 7 PETSC_DEVICE_FUNC_DECL PetscReal polevl_10(PetscReal x, const PetscReal coef[]) 8 { 9 PetscReal ans; 10 PetscInt i; 11 ans = coef[0]; 12 for (i=1; i<11; i++) ans = ans * x + coef[i]; 13 return(ans); 14 } 15 PETSC_DEVICE_FUNC_DECL PetscReal polevl_9(PetscReal x, const PetscReal coef[]) 16 { 17 PetscReal ans; 18 PetscInt i; 19 ans = coef[0]; 20 for (i=1; i<10; i++) ans = ans * x + coef[i]; 21 return(ans); 22 } 23 /* 24 * Complete elliptic integral of the second kind 25 */ 26 PETSC_DEVICE_FUNC_DECL void ellipticE(PetscReal x,PetscReal *ret) 27 { 28 #if defined(PETSC_USE_REAL_SINGLE) 29 static const PetscReal P2[] = { 30 1.53552577301013293365E-4F, 31 2.50888492163602060990E-3F, 32 8.68786816565889628429E-3F, 33 1.07350949056076193403E-2F, 34 7.77395492516787092951E-3F, 35 7.58395289413514708519E-3F, 36 1.15688436810574127319E-2F, 37 2.18317996015557253103E-2F, 38 5.68051945617860553470E-2F, 39 4.43147180560990850618E-1F, 40 1.00000000000000000299E0F 41 }; 42 static const PetscReal Q2[] = { 43 3.27954898576485872656E-5F, 44 1.00962792679356715133E-3F, 45 6.50609489976927491433E-3F, 46 1.68862163993311317300E-2F, 47 2.61769742454493659583E-2F, 48 3.34833904888224918614E-2F, 49 4.27180926518931511717E-2F, 50 5.85936634471101055642E-2F, 51 9.37499997197644278445E-2F, 52 2.49999999999888314361E-1F 53 }; 54 #else 55 static const PetscReal P2[] = { 56 1.53552577301013293365E-4, 57 2.50888492163602060990E-3, 58 8.68786816565889628429E-3, 59 1.07350949056076193403E-2, 60 7.77395492516787092951E-3, 61 7.58395289413514708519E-3, 62 1.15688436810574127319E-2, 63 2.18317996015557253103E-2, 64 5.68051945617860553470E-2, 65 4.43147180560990850618E-1, 66 1.00000000000000000299E0 67 }; 68 static const PetscReal Q2[] = { 69 3.27954898576485872656E-5, 70 1.00962792679356715133E-3, 71 6.50609489976927491433E-3, 72 1.68862163993311317300E-2, 73 2.61769742454493659583E-2, 74 3.34833904888224918614E-2, 75 4.27180926518931511717E-2, 76 5.85936634471101055642E-2, 77 9.37499997197644278445E-2, 78 2.49999999999888314361E-1 79 }; 80 #endif 81 x = 1 - x; /* where m = 1 - m1 */ 82 *ret = polevl_10(x,P2) - PetscLogReal(x) * (x * polevl_9(x,Q2)); 83 } 84 /* 85 * Complete elliptic integral of the first kind 86 */ 87 PETSC_DEVICE_FUNC_DECL void ellipticK(PetscReal x,PetscReal *ret) 88 { 89 #if defined(PETSC_USE_REAL_SINGLE) 90 static const PetscReal P1[] = 91 { 92 1.37982864606273237150E-4F, 93 2.28025724005875567385E-3F, 94 7.97404013220415179367E-3F, 95 9.85821379021226008714E-3F, 96 6.87489687449949877925E-3F, 97 6.18901033637687613229E-3F, 98 8.79078273952743772254E-3F, 99 1.49380448916805252718E-2F, 100 3.08851465246711995998E-2F, 101 9.65735902811690126535E-2F, 102 1.38629436111989062502E0F 103 }; 104 static const PetscReal Q1[] = 105 { 106 2.94078955048598507511E-5F, 107 9.14184723865917226571E-4F, 108 5.94058303753167793257E-3F, 109 1.54850516649762399335E-2F, 110 2.39089602715924892727E-2F, 111 3.01204715227604046988E-2F, 112 3.73774314173823228969E-2F, 113 4.88280347570998239232E-2F, 114 7.03124996963957469739E-2F, 115 1.24999999999870820058E-1F, 116 4.99999999999999999821E-1F 117 }; 118 #else 119 static const PetscReal P1[] = 120 { 121 1.37982864606273237150E-4, 122 2.28025724005875567385E-3, 123 7.97404013220415179367E-3, 124 9.85821379021226008714E-3, 125 6.87489687449949877925E-3, 126 6.18901033637687613229E-3, 127 8.79078273952743772254E-3, 128 1.49380448916805252718E-2, 129 3.08851465246711995998E-2, 130 9.65735902811690126535E-2, 131 1.38629436111989062502E0 132 }; 133 static const PetscReal Q1[] = 134 { 135 2.94078955048598507511E-5, 136 9.14184723865917226571E-4, 137 5.94058303753167793257E-3, 138 1.54850516649762399335E-2, 139 2.39089602715924892727E-2, 140 3.01204715227604046988E-2, 141 3.73774314173823228969E-2, 142 4.88280347570998239232E-2, 143 7.03124996963957469739E-2, 144 1.24999999999870820058E-1, 145 4.99999999999999999821E-1 146 }; 147 #endif 148 x = 1 - x; /* where m = 1 - m1 */ 149 *ret = polevl_10(x,P1) - PetscLogReal(x) * polevl_10(x,Q1); 150 } 151 PETSC_DEVICE_FUNC_DECL void LandauTensor2D(const PetscReal x[], const PetscReal rp, const PetscReal zp, PetscReal Ud[][2], PetscReal Uk[][2], const PetscReal mask) 152 { 153 PetscReal l,s,r=x[0],z=x[1],i1func,i2func,i3func,ks,es,pi4pow,sqrt_1s,r2,rp2,r2prp2,zmzp,zmzp2,tt; 154 //PetscReal mask /* = !!(r!=rp || z!=zp) */; 155 /* !!(zmzp2 > 1.e-12 || (r-rp) > 1.e-12 || (r-rp) < -1.e-12); */ 156 r2=PetscSqr(r); 157 zmzp=z-zp; 158 rp2=PetscSqr(rp); 159 zmzp2=PetscSqr(zmzp); 160 r2prp2=r2+rp2; 161 l = r2 + rp2 + zmzp2; 162 /* if (zmzp2 > PETSC_SMALL) mask = 1; */ 163 /* else if ((tt=(r-rp)) > PETSC_SMALL) mask = 1; */ 164 /* else if (tt < -PETSC_SMALL) mask = 1; */ 165 /* else mask = 0; */ 166 s = mask*2*r*rp/l; /* mask for vectorization */ 167 tt = 1./(1+s); 168 pi4pow = 4*PETSC_PI*LANDAU_INVSQRT(PetscSqr(l)*l); 169 sqrt_1s = LANDAU_SQRT(1.+s); 170 /* sp.ellipe(2.*s/(1.+s)) */ 171 ellipticE(2*s*tt,&es); /* 44 flops * 2 + 75 = 163 flops including 2 logs, 1 sqrt, 1 pow, 21 mult */ 172 /* sp.ellipk(2.*s/(1.+s)) */ 173 ellipticK(2*s*tt,&ks); /* 44 flops + 75 in rest, 21 mult */ 174 /* mask is needed here just for single precision */ 175 i2func = 2./((1-s)*sqrt_1s) * es; 176 i1func = 4./(PetscSqr(s)*sqrt_1s + PETSC_MACHINE_EPSILON) * mask * (ks - (1.+s) * es); 177 i3func = 2./((1-s)*(s)*sqrt_1s + PETSC_MACHINE_EPSILON) * (es - (1-s) * ks); 178 Ud[0][0]= pi4pow*(rp2*i1func+PetscSqr(zmzp)*i2func); 179 Ud[0][1]=Ud[1][0]=Uk[0][1]= -pi4pow*(zmzp)*(r*i2func-rp*i3func); 180 Uk[1][1]=Ud[1][1]= pi4pow*((r2prp2)*i2func-2*r*rp*i3func)*mask; 181 Uk[0][0]= pi4pow*(zmzp2*i3func+r*rp*i1func); 182 Uk[1][0]= -pi4pow*(zmzp)*(r*i3func-rp*i2func); /* 48 mults + 21 + 21 = 90 mults and divs */ 183 } 184 #else 185 /* integration point functions */ 186 /* Evaluates the tensor U=(I-(x-y)(x-y)/(x-y)^2)/|x-y| at point x,y */ 187 /* if x==y we will return zero. This is not the correct result */ 188 /* since the tensor diverges for x==y but when integrated */ 189 /* the divergent part is antisymmetric and vanishes. This is not */ 190 /* trivial, but can be proven. */ 191 PETSC_DEVICE_FUNC_DECL void LandauTensor3D(const PetscReal x1[], const PetscReal xp, const PetscReal yp, const PetscReal zp, PetscReal U[][3], PetscReal mask) 192 { 193 PetscReal dx[3],inorm3,inorm,inorm2,norm2,x2[] = {xp,yp,zp}; 194 PetscInt d; 195 for (d = 0, norm2 = PETSC_MACHINE_EPSILON; d < 3; ++d) { 196 dx[d] = x2[d] - x1[d]; 197 norm2 += dx[d] * dx[d]; 198 } 199 inorm2 = mask/norm2; 200 inorm = LANDAU_SQRT(inorm2); 201 inorm3 = inorm2*inorm; 202 for (d = 0; d < 3; ++d) U[d][d] = inorm - inorm3 * dx[d] * dx[d]; 203 U[1][0] = U[0][1] = -inorm3 * dx[0] * dx[1]; 204 U[1][2] = U[2][1] = -inorm3 * dx[2] * dx[1]; 205 U[2][0] = U[0][2] = -inorm3 * dx[0] * dx[2]; 206 } 207 #endif 208