sph_harm.hpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499
  1. #ifndef _SCTL_SPH_HARM_HPP_
  2. #define _SCTL_SPH_HARM_HPP_
  3. #define SCTL_SHMAXDEG 1024
  4. #include SCTL_INCLUDE(matrix.hpp)
  5. #include SCTL_INCLUDE(fft_wrapper.hpp)
  6. #include SCTL_INCLUDE(common.hpp)
  7. namespace SCTL_NAMESPACE {
  8. enum class SHCArrange {
  9. // (p+1) x (p+1) complex elements in row-major order.
  10. // A : { A(0,0), A(0,1), ... A(0,p), A(1,0), ... A(p,p) }
  11. // where, A(n,m) = { Ar(n,m), Ai(n,m) } (real and imaginary parts)
  12. ALL,
  13. // (p+1)(p+2)/2 complex elements in row-major order (lower triangular part)
  14. // A : { A(0,0), A(1,0), A(1,1), A(2,0), A(2,1), A(2,2), ... A(p,p) }
  15. // where, A(n,m) = { Ar(n,m), Ai(n,m) } (real and imaginary parts)
  16. ROW_MAJOR,
  17. // (p+1)(p+1) real elements in col-major order (non-zero lower triangular part)
  18. // A : { Ar(0,0), Ar(1,0), ... Ar(p,0), Ar(1,1), ... Ar(p,1), Ai(1,1), ... Ai(p,1), ..., Ar(p,p), Ai(p,p)
  19. // where, A(n,m) = { Ar(n,m), Ai(n,m) } (real and imaginary parts)
  20. COL_MAJOR_NONZERO
  21. };
  22. template <class Real> class SphericalHarmonics{
  23. static constexpr Integer COORD_DIM = 3;
  24. public:
  25. // Scalar Spherical Harmonics
  26. /**
  27. * \brief Compute spherical harmonic coefficients from grid values.
  28. * \param[in] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
  29. * \param[in] Nt Number of grid points \theta \in (1,pi).
  30. * \param[in] Np Number of grid points \phi \in (1,2*pi).
  31. * \param[in] p Order of spherical harmonic expansion.
  32. * \param[in] arrange Arrangement of the coefficients.
  33. * \param[out] S Spherical harmonic coefficients.
  34. */
  35. static void Grid2SHC(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& S, SHCArrange arrange);
  36. /**
  37. * \brief Evaluate grid values from spherical harmonic coefficients.
  38. * \param[in] S Spherical harmonic coefficients.
  39. * \param[in] arrange Arrangement of the coefficients.
  40. * \param[in] p Order of spherical harmonic expansion.
  41. * \param[in] Nt Number of grid points \theta \in (1,pi).
  42. * \param[in] Np Number of grid points \phi \in (1,2*pi).
  43. * \param[out] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
  44. * \param[out] X_theta \theta derivative of X evaluated at grid points.
  45. * \param[out] X_phi \phi derivative of X evaluated at grid points.
  46. */
  47. static void SHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p, Long Nt, Long Np, Vector<Real>* X, Vector<Real>* X_theta=nullptr, Vector<Real>* X_phi=nullptr);
  48. /**
  49. * \brief Evaluate point values from spherical harmonic coefficients.
  50. * \param[in] S Spherical harmonic coefficients.
  51. * \param[in] arrange Arrangement of the coefficients.
  52. * \param[in] p Order of spherical harmonic expansion.
  53. * \param[in] theta_phi Evaluation coordinates given as {t0,p0, t1,p1, ... }.
  54. * \param[out] X Evaluated values {X0, X1, ... }.
  55. */
  56. static void SHCEval(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& theta_phi, Vector<Real>& X);
  57. static void SHC2Pole(const Vector<Real>& S, SHCArrange arrange, Long p, Vector<Real>& P);
  58. static void WriteVTK(const char* fname, const Vector<Real>* S, const Vector<Real>* f_val, SHCArrange arrange, Long p_in, Long p_out, Real period=0, const Comm& comm = Comm::World());
  59. // Vector Spherical Harmonics
  60. /**
  61. * \brief Compute vector spherical harmonic coefficients from grid values.
  62. * \param[in] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), ... , Y(t0,p0), ... , Z(t0,p0), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
  63. * \param[in] Nt Number of grid points \theta \in (1,pi).
  64. * \param[in] Np Number of grid points \phi \in (1,2*pi).
  65. * \param[in] p Order of spherical harmonic expansion.
  66. * \param[in] arrange Arrangement of the coefficients.
  67. * \param[out] S Vector spherical harmonic coefficients.
  68. */
  69. static void Grid2VecSHC(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& S, SHCArrange arrange);
  70. /**
  71. * \brief Evaluate grid values from vector spherical harmonic coefficients.
  72. * \param[in] S Vector spherical harmonic coefficients.
  73. * \param[in] arrange Arrangement of the coefficients.
  74. * \param[in] p Order of spherical harmonic expansion.
  75. * \param[in] Nt Number of grid points \theta \in (1,pi).
  76. * \param[in] Np Number of grid points \phi \in (1,2*pi).
  77. * \param[out] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... , Y(t0,p0), ... , Z(t0,p0), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
  78. */
  79. static void VecSHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p, Long Nt, Long Np, Vector<Real>& X);
  80. /**
  81. * \brief Evaluate point values from vector spherical harmonic coefficients.
  82. * \param[in] S Vector spherical harmonic coefficients.
  83. * \param[in] arrange Arrangement of the coefficients.
  84. * \param[in] p Order of spherical harmonic expansion.
  85. * \param[in] theta_phi Evaluation coordinates given as {t0,p0, t1,p1, ... }.
  86. * \param[out] X Evaluated values {X0,Y0,Z0, X1,Y1,Z1, ... }.
  87. */
  88. static void VecSHCEval(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& theta_phi, Vector<Real>& X);
  89. /**
  90. * \brief Evaluate Stokes single-layer operator at point values from the vector spherical harmonic coefficients for the density.
  91. * \param[in] S Vector spherical harmonic coefficients.
  92. * \param[in] arrange Arrangement of the coefficients.
  93. * \param[in] p Order of spherical harmonic expansion.
  94. * \param[in] Evaluation coordinates given as {x0,y0,z0, x1,y1,z1, ... }.
  95. * \param[out] U Evaluated values {Ux0,Uy0,Uz0, Ux1,Uy1,Uz1, ... }.
  96. */
  97. static void StokesEvalSL(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, bool interior, Vector<Real>& U);
  98. /**
  99. * \brief Evaluate Stokes double-layer operator at point values from the vector spherical harmonic coefficients for the density.
  100. * \param[in] S Vector spherical harmonic coefficients.
  101. * \param[in] arrange Arrangement of the coefficients.
  102. * \param[in] p Order of spherical harmonic expansion.
  103. * \param[in] Evaluation coordinates given as {x0,y0,z0, x1,y1,z1, ... }.
  104. * \param[out] U Evaluated values {Ux0,Uy0,Uz0, Ux1,Uy1,Uz1, ... }.
  105. */
  106. static void StokesEvalDL(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, bool interior, Vector<Real>& U);
  107. static void StokesEvalKL(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, const Vector<Real>& norm, bool interior, Vector<Real>& U);
  108. static void StokesEvalKSelf(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, bool interior, Vector<Real>& U);
  109. static void test_stokes() {
  110. int p = 6;
  111. int dof = 3;
  112. int Nt = 100, Np = 200;
  113. auto print_coeff = [&](Vector<Real> S) {
  114. Long idx=0;
  115. for (Long k=0;k<dof;k++) {
  116. for (Long n=0;n<=p;n++) {
  117. std::cout<<Vector<Real>(2*n+2, S.begin()+idx);
  118. idx+=2*n+2;
  119. }
  120. }
  121. std::cout<<'\n';
  122. };
  123. Vector<Real> Fcoeff(dof*(p+1)*(p+2));
  124. for (Long i=0;i<Fcoeff.Dim();i++) Fcoeff[i]=i+1;
  125. //Fcoeff = 0; Fcoeff[2] = 1;
  126. print_coeff(Fcoeff);
  127. Vector<Real> Fgrid;
  128. VecSHC2Grid(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, Nt, Np, Fgrid);
  129. Matrix<Real>(Fgrid.Dim()/3,3,Fgrid.begin(),false) = Matrix<Real>(3,Fgrid.Dim()/3,Fgrid.begin(),false).Transpose();
  130. const Vector<Real> CosTheta = LegendreNodes(Nt-1);
  131. const Vector<Real> LWeights = LegendreWeights(Nt-1);
  132. auto stokes_evalSL = [&](const Vector<Real>& trg, Vector<Real>& Sf) {
  133. Sf.ReInit(3);
  134. Sf=0;
  135. Real s = 1/(8*const_pi<Real>());
  136. for (Long i=0;i<Nt;i++) {
  137. Real cos_theta = CosTheta[i];
  138. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  139. for (Long j=0;j<Np;j++) {
  140. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  141. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  142. Real qw = LWeights[i]*2*const_pi<Real>()/Np;
  143. Real x[3], dr[3], f[3];
  144. f[0] = Fgrid[(i*Np+j)*3+0];
  145. f[1] = Fgrid[(i*Np+j)*3+1];
  146. f[2] = Fgrid[(i*Np+j)*3+2];
  147. x[0] = sin_theta*cos_phi;
  148. x[1] = sin_theta*sin_phi;
  149. x[2] = cos_theta;
  150. dr[0] = x[0]-trg[0];
  151. dr[1] = x[1]-trg[1];
  152. dr[2] = x[2]-trg[2];
  153. Real oor2 = 1/(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  154. Real oor1 = sqrt(oor2);
  155. Real oor3 = oor2*oor1;
  156. Real rdotf = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  157. Sf[0] += s*(f[0]*oor1 + dr[0]*rdotf*oor3) * qw;
  158. Sf[1] += s*(f[1]*oor1 + dr[1]*rdotf*oor3) * qw;
  159. Sf[2] += s*(f[2]*oor1 + dr[2]*rdotf*oor3) * qw;
  160. }
  161. }
  162. };
  163. auto stokes_evalDL = [&](const Vector<Real>& trg, Vector<Real>& Sf) {
  164. Sf.ReInit(3);
  165. Sf=0;
  166. Real s = 6/(8*const_pi<Real>());
  167. for (Long i=0;i<Nt;i++) {
  168. Real cos_theta = CosTheta[i];
  169. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  170. for (Long j=0;j<Np;j++) {
  171. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  172. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  173. Real qw = LWeights[i]*2*const_pi<Real>()/Np;
  174. Real x[3], dr[3], f[3], n[3];
  175. f[0] = Fgrid[(i*Np+j)*3+0];
  176. f[1] = Fgrid[(i*Np+j)*3+1];
  177. f[2] = Fgrid[(i*Np+j)*3+2];
  178. x[0] = sin_theta*cos_phi;
  179. x[1] = sin_theta*sin_phi;
  180. x[2] = cos_theta;
  181. dr[0] = x[0]-trg[0];
  182. dr[1] = x[1]-trg[1];
  183. dr[2] = x[2]-trg[2];
  184. n[0] = x[0];
  185. n[1] = x[1];
  186. n[2] = x[2];
  187. Real oor2 = 1/(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  188. Real oor5 = oor2*oor2*sqrt(oor2);
  189. Real rdotn = dr[0]*n[0]+dr[1]*n[1]+dr[2]*n[2];
  190. Real rdotf = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  191. Sf[0] += -s*dr[0]*rdotn*rdotf*oor5 * qw;
  192. Sf[1] += -s*dr[1]*rdotn*rdotf*oor5 * qw;
  193. Sf[2] += -s*dr[2]*rdotn*rdotf*oor5 * qw;
  194. }
  195. }
  196. };
  197. auto stokes_evalKL = [&](const Vector<Real>& trg, const Vector<Real>& nor, Vector<Real>& Sf) {
  198. Sf.ReInit(3);
  199. Sf=0;
  200. Real scal = 1/(8*const_pi<Real>());
  201. for (Long i=0;i<Nt;i++) {
  202. Real cos_theta = CosTheta[i];
  203. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  204. for (Long j=0;j<Np;j++) {
  205. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  206. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  207. Real qw = LWeights[i]*2*const_pi<Real>()/Np; // quadrature weights * area-element
  208. Real f[3]; // source density
  209. f[0] = Fgrid[(i*Np+j)*3+0];
  210. f[1] = Fgrid[(i*Np+j)*3+1];
  211. f[2] = Fgrid[(i*Np+j)*3+2];
  212. Real x[3]; // source coordinates
  213. x[0] = sin_theta*cos_phi;
  214. x[1] = sin_theta*sin_phi;
  215. x[2] = cos_theta;
  216. Real dr[3];
  217. dr[0] = trg[0] - x[0];
  218. dr[1] = trg[1] - x[1];
  219. dr[2] = trg[2] - x[2];
  220. Real invr = 1 / sqrt(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  221. Real invr2 = invr*invr;
  222. Real invr3 = invr2*invr;
  223. Real invr5 = invr2*invr3;
  224. Real fdotr = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  225. Real du[9];
  226. du[0] = ( fdotr*invr3 - 3*dr[0]*dr[0]*fdotr*invr5) * scal;
  227. du[1] = ((dr[0]*f[1]-dr[1]*f[0])*invr3 - 3*dr[0]*dr[1]*fdotr*invr5) * scal;
  228. du[2] = ((dr[0]*f[2]-dr[2]*f[0])*invr3 - 3*dr[0]*dr[2]*fdotr*invr5) * scal;
  229. du[3] = ((dr[1]*f[0]-dr[0]*f[1])*invr3 - 3*dr[1]*dr[0]*fdotr*invr5) * scal;
  230. du[4] = ( fdotr*invr3 - 3*dr[1]*dr[1]*fdotr*invr5) * scal;
  231. du[5] = ((dr[1]*f[2]-dr[2]*f[1])*invr3 - 3*dr[1]*dr[2]*fdotr*invr5) * scal;
  232. du[6] = ((dr[2]*f[0]-dr[0]*f[2])*invr3 - 3*dr[2]*dr[0]*fdotr*invr5) * scal;
  233. du[7] = ((dr[2]*f[1]-dr[1]*f[2])*invr3 - 3*dr[2]*dr[1]*fdotr*invr5) * scal;
  234. du[8] = ( fdotr*invr3 - 3*dr[2]*dr[2]*fdotr*invr5) * scal;
  235. Real p = (2*fdotr*invr3) * scal;
  236. Real K[9];
  237. K[0] = du[0] + du[0] - p; K[1] = du[1] + du[3] - 0; K[2] = du[2] + du[6] - 0;
  238. K[3] = du[3] + du[1] - 0; K[4] = du[4] + du[4] - p; K[5] = du[5] + du[7] - 0;
  239. K[6] = du[6] + du[2] - 0; K[7] = du[7] + du[5] - 0; K[8] = du[8] + du[8] - p;
  240. Sf[0] += (K[0]*nor[0] + K[1]*nor[1] + K[2]*nor[2]) * qw;
  241. Sf[1] += (K[3]*nor[0] + K[4]*nor[1] + K[5]*nor[2]) * qw;
  242. Sf[2] += (K[6]*nor[0] + K[7]*nor[1] + K[8]*nor[2]) * qw;
  243. }
  244. }
  245. };
  246. for (Long i = 0; i < 40; i++) { // Evaluate
  247. Real R0 = (0.01 + i/20.0);
  248. Vector<Real> x(3), n(3);
  249. x[0] = drand48()-0.5;
  250. x[1] = drand48()-0.5;
  251. x[2] = drand48()-0.5;
  252. n[0] = drand48()-0.5;
  253. n[1] = drand48()-0.5;
  254. n[2] = drand48()-0.5;
  255. Real R = sqrt<Real>(x[0]*x[0]+x[1]*x[1]+x[2]*x[2]);
  256. x[0] *= R0 / R;
  257. x[1] *= R0 / R;
  258. x[2] *= R0 / R;
  259. Vector<Real> Sf, Sf_;
  260. Vector<Real> Df, Df_;
  261. Vector<Real> Kf, Kf_;
  262. StokesEvalSL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, R0<1, Sf);
  263. StokesEvalDL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, R0<1, Df);
  264. StokesEvalKL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, n, R0<1, Kf);
  265. stokes_evalSL(x, Sf_);
  266. stokes_evalDL(x, Df_);
  267. stokes_evalKL(x, n, Kf_);
  268. auto max_val = [](const Vector<Real>& v) {
  269. Real max_v = 0;
  270. for (auto& x : v) max_v = std::max(max_v, fabs(x));
  271. return max_v;
  272. };
  273. auto errSL = (Sf-Sf_)/max_val(Sf+0.01);
  274. auto errDL = (Df-Df_)/max_val(Df+0.01);
  275. auto errKL = (Kf-Kf_)/max_val(Kf+0.01);
  276. for (auto& x:errSL) x=log(fabs(x))/log(10);
  277. for (auto& x:errDL) x=log(fabs(x))/log(10);
  278. for (auto& x:errKL) x=log(fabs(x))/log(10);
  279. std::cout<<"R = "<<(0.01 + i/20.0)<<"; SL-error = ";
  280. std::cout<<errSL;
  281. std::cout<<"R = "<<(0.01 + i/20.0)<<"; DL-error = ";
  282. std::cout<<errDL;
  283. std::cout<<"R = "<<(0.01 + i/20.0)<<"; KL-error = ";
  284. std::cout<<errKL;
  285. }
  286. Clear();
  287. }
  288. static void test() {
  289. int p = 3;
  290. int dof = 1;
  291. int Nt = p+1, Np = 2*p+1;
  292. auto print_coeff = [&](Vector<Real> S) {
  293. Long idx=0;
  294. for (Long k=0;k<dof;k++) {
  295. for (Long n=0;n<=p;n++) {
  296. std::cout<<Vector<Real>(2*n+2, S.begin()+idx);
  297. idx+=2*n+2;
  298. }
  299. }
  300. std::cout<<'\n';
  301. };
  302. Vector<Real> theta_phi;
  303. { // Set theta_phi
  304. Vector<Real> leg_nodes = LegendreNodes(Nt-1);
  305. for (Long i=0;i<Nt;i++) {
  306. for (Long j=0;j<Np;j++) {
  307. theta_phi.PushBack(acos(leg_nodes[i]));
  308. theta_phi.PushBack(j * 2 * const_pi<Real>() / Np);
  309. }
  310. }
  311. }
  312. int Ncoeff = (p + 1) * (p + 1);
  313. Vector<Real> Xcoeff(dof * Ncoeff), Xgrid;
  314. for (int i=0;i<Xcoeff.Dim();i++) Xcoeff[i]=i+1;
  315. SHC2Grid(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, Nt, Np, &Xgrid);
  316. std::cout<<Matrix<Real>(Nt*dof, Np, Xgrid.begin())<<'\n';
  317. {
  318. Vector<Real> val;
  319. SHCEval(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, theta_phi, val);
  320. Matrix<Real>(dof, val.Dim()/dof, val.begin(), false) = Matrix<Real>(val.Dim()/dof, dof, val.begin()).Transpose();
  321. std::cout<<Matrix<Real>(val.Dim()/Np, Np, val.begin()) - Matrix<Real>(Nt*dof, Np, Xgrid.begin())+1e-10<<'\n';
  322. }
  323. Grid2SHC(Xgrid, Nt, Np, p, Xcoeff, sctl::SHCArrange::ROW_MAJOR);
  324. print_coeff(Xcoeff);
  325. //SphericalHarmonics<Real>::WriteVTK("test", nullptr, &Xcoeff, sctl::SHCArrange::ROW_MAJOR, p, 32);
  326. Clear();
  327. }
  328. /**
  329. * \brief Clear all precomputed data. This must be done before the program exits to avoid memory leaks.
  330. */
  331. static void Clear() { MatrixStore().Resize(0); }
  332. private:
  333. // Probably don't work anymore, need to be updated :(
  334. static void SHC2GridTranspose(const Vector<Real>& X, Long p0, Long p1, Vector<Real>& S);
  335. static void RotateAll(const Vector<Real>& S, Long p0, Long dof, Vector<Real>& S_);
  336. static void RotateTranspose(const Vector<Real>& S_, Long p0, Long dof, Vector<Real>& S);
  337. static void StokesSingularInteg(const Vector<Real>& S, Long p0, Long p1, Vector<Real>* SLMatrix=nullptr, Vector<Real>* DLMatrix=nullptr);
  338. static void Grid2SHC_(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& B1);
  339. static void SHCArrange0(const Vector<Real>& B1, Long p, Vector<Real>& S, SHCArrange arrange);
  340. static void SHC2Grid_(const Vector<Real>& S, Long p, Long Nt, Long Np, Vector<Real>* X, Vector<Real>* X_theta=nullptr, Vector<Real>* X_phi=nullptr);
  341. static void SHCArrange1(const Vector<Real>& S_in, SHCArrange arrange_out, Long p, Vector<Real>& S_out);
  342. /**
  343. * \brief Computes all the Associated Legendre Polynomials (normalized) up to the specified degree.
  344. * \param[in] degree The degree up to which the Legendre polynomials have to be computed.
  345. * \param[in] X The input values for which the polynomials have to be computed.
  346. * \param[in] N The number of input points.
  347. * \param[out] poly_val The output array of size (degree+1)*(degree+2)*N/2 containing the computed polynomial values.
  348. * The output values are in the order:
  349. * P(n,m)[i] => {P(0,0)[0], P(0,0)[1], ..., P(0,0)[N-1], P(1,0)[0], ..., P(1,0)[N-1],
  350. * P(2,0)[0], ..., P(degree,0)[N-1], P(1,1)[0], ...,P(2,1)[0], ..., P(degree,degree)[N-1]}
  351. */
  352. static void LegPoly(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  353. static void LegPoly_(Vector<Real>& poly_val, const Vector<Real>& theta, Long degree);
  354. static void LegPolyDeriv(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  355. static void LegPolyDeriv_(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  356. static const Vector<Real>& LegendreNodes(Long p1);
  357. static const Vector<Real>& LegendreWeights(Long p1);
  358. static const Vector<Real>& SingularWeights(Long p1);
  359. static const Matrix<Real>& MatFourier(Long p0, Long p1);
  360. static const Matrix<Real>& MatFourierInv(Long p0, Long p1);
  361. static const Matrix<Real>& MatFourierGrad(Long p0, Long p1);
  362. static const FFT<Real>& OpFourier(Long Np);
  363. static const FFT<Real>& OpFourierInv(Long Np);
  364. static const std::vector<Matrix<Real>>& MatLegendre(Long p0, Long p1);
  365. static const std::vector<Matrix<Real>>& MatLegendreInv(Long p0, Long p1);
  366. static const std::vector<Matrix<Real>>& MatLegendreGrad(Long p0, Long p1);
  367. // Evaluate all Spherical Harmonic basis functions up to order p at (theta, phi) coordinates.
  368. static void SHBasisEval(Long p, const Vector<Real>& theta_phi, Matrix<Real>& M);
  369. static void VecSHBasisEval(Long p, const Vector<Real>& theta_phi, Matrix<Real>& M);
  370. static const std::vector<Matrix<Real>>& MatRotate(Long p0);
  371. template <bool SLayer, bool DLayer> static void StokesSingularInteg_(const Vector<Real>& X0, Long p0, Long p1, Vector<Real>& SL, Vector<Real>& DL);
  372. struct MatrixStorage{
  373. MatrixStorage(){
  374. const Long size = SCTL_SHMAXDEG;
  375. Resize(size);
  376. }
  377. void Resize(Long size){
  378. Qx_ .resize(size);
  379. Qw_ .resize(size);
  380. Sw_ .resize(size);
  381. Mf_ .resize(size*size);
  382. Mdf_.resize(size*size);
  383. Ml_ .resize(size*size);
  384. Mdl_.resize(size*size);
  385. Mr_ .resize(size);
  386. Mfinv_ .resize(size*size);
  387. Mlinv_ .resize(size*size);
  388. Mfft_.resize(size);
  389. Mfftinv_.resize(size);
  390. }
  391. std::vector<Vector<Real>> Qx_;
  392. std::vector<Vector<Real>> Qw_;
  393. std::vector<Vector<Real>> Sw_;
  394. std::vector<Matrix<Real>> Mf_ ;
  395. std::vector<Matrix<Real>> Mdf_;
  396. std::vector<std::vector<Matrix<Real>>> Ml_ ;
  397. std::vector<std::vector<Matrix<Real>>> Mdl_;
  398. std::vector<std::vector<Matrix<Real>>> Mr_;
  399. std::vector<Matrix<Real>> Mfinv_ ;
  400. std::vector<std::vector<Matrix<Real>>> Mlinv_ ;
  401. std::vector<FFT<Real>> Mfft_;
  402. std::vector<FFT<Real>> Mfftinv_;
  403. };
  404. static MatrixStorage& MatrixStore(){
  405. static MatrixStorage storage;
  406. if (!storage.Qx_.size()) storage.Resize(SCTL_SHMAXDEG);
  407. return storage;
  408. }
  409. };
  410. template class SphericalHarmonics<double>;
  411. } // end namespace
  412. #include SCTL_INCLUDE(sph_harm.txx)
  413. #endif // _SCTL_SPH_HARM_HPP_