sph_harm.hpp 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493
  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] cos_theta_phi Evaluation coordinates given as {cos(t0),p0, cos(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>& cos_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] cos_theta_phi Evaluation coordinates given as {cos(t0),p0, cos(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>& cos_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 test_stokes() {
  109. int p = 6;
  110. int dof = 3;
  111. int Nt = 100, Np = 200;
  112. auto print_coeff = [&](Vector<Real> S) {
  113. Long idx=0;
  114. for (Long k=0;k<dof;k++) {
  115. for (Long n=0;n<=p;n++) {
  116. std::cout<<Vector<Real>(2*n+2, S.begin()+idx);
  117. idx+=2*n+2;
  118. }
  119. }
  120. std::cout<<'\n';
  121. };
  122. Vector<Real> Fcoeff(dof*(p+1)*(p+2));
  123. for (Long i=0;i<Fcoeff.Dim();i++) Fcoeff[i]=i+1;
  124. //Fcoeff = 0; Fcoeff[2] = 1;
  125. print_coeff(Fcoeff);
  126. Vector<Real> Fgrid;
  127. VecSHC2Grid(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, Nt, Np, Fgrid);
  128. Matrix<Real>(Fgrid.Dim()/3,3,Fgrid.begin(),false) = Matrix<Real>(3,Fgrid.Dim()/3,Fgrid.begin(),false).Transpose();
  129. const Vector<Real> CosTheta = LegendreNodes(Nt-1);
  130. const Vector<Real> LWeights = LegendreWeights(Nt-1);
  131. auto stokes_evalSL = [&](const Vector<Real>& trg, Vector<Real>& Sf) {
  132. Sf.ReInit(3);
  133. Sf=0;
  134. Real s = 1/(8*const_pi<Real>());
  135. for (Long i=0;i<Nt;i++) {
  136. Real cos_theta = CosTheta[i];
  137. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  138. for (Long j=0;j<Np;j++) {
  139. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  140. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  141. Real qw = LWeights[i]*2*const_pi<Real>()/Np;
  142. Real x[3], dr[3], f[3];
  143. f[0] = Fgrid[(i*Np+j)*3+0];
  144. f[1] = Fgrid[(i*Np+j)*3+1];
  145. f[2] = Fgrid[(i*Np+j)*3+2];
  146. x[0] = sin_theta*cos_phi;
  147. x[1] = sin_theta*sin_phi;
  148. x[2] = cos_theta;
  149. dr[0] = x[0]-trg[0];
  150. dr[1] = x[1]-trg[1];
  151. dr[2] = x[2]-trg[2];
  152. Real oor2 = 1/(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  153. Real oor1 = sqrt(oor2);
  154. Real oor3 = oor2*oor1;
  155. Real rdotf = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  156. Sf[0] += s*(f[0]*oor1 + dr[0]*rdotf*oor3) * qw;
  157. Sf[1] += s*(f[1]*oor1 + dr[1]*rdotf*oor3) * qw;
  158. Sf[2] += s*(f[2]*oor1 + dr[2]*rdotf*oor3) * qw;
  159. }
  160. }
  161. };
  162. auto stokes_evalDL = [&](const Vector<Real>& trg, Vector<Real>& Sf) {
  163. Sf.ReInit(3);
  164. Sf=0;
  165. Real s = 6/(8*const_pi<Real>());
  166. for (Long i=0;i<Nt;i++) {
  167. Real cos_theta = CosTheta[i];
  168. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  169. for (Long j=0;j<Np;j++) {
  170. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  171. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  172. Real qw = LWeights[i]*2*const_pi<Real>()/Np;
  173. Real x[3], dr[3], f[3], n[3];
  174. f[0] = Fgrid[(i*Np+j)*3+0];
  175. f[1] = Fgrid[(i*Np+j)*3+1];
  176. f[2] = Fgrid[(i*Np+j)*3+2];
  177. x[0] = sin_theta*cos_phi;
  178. x[1] = sin_theta*sin_phi;
  179. x[2] = cos_theta;
  180. dr[0] = x[0]-trg[0];
  181. dr[1] = x[1]-trg[1];
  182. dr[2] = x[2]-trg[2];
  183. n[0] = x[0];
  184. n[1] = x[1];
  185. n[2] = x[2];
  186. Real oor2 = 1/(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  187. Real oor5 = oor2*oor2*sqrt(oor2);
  188. Real rdotn = dr[0]*n[0]+dr[1]*n[1]+dr[2]*n[2];
  189. Real rdotf = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  190. Sf[0] += -s*dr[0]*rdotn*rdotf*oor5 * qw;
  191. Sf[1] += -s*dr[1]*rdotn*rdotf*oor5 * qw;
  192. Sf[2] += -s*dr[2]*rdotn*rdotf*oor5 * qw;
  193. }
  194. }
  195. };
  196. auto stokes_evalKL = [&](const Vector<Real>& trg, const Vector<Real>& nor, Vector<Real>& Sf) {
  197. Sf.ReInit(3);
  198. Sf=0;
  199. Real scal = 1/(8*const_pi<Real>());
  200. for (Long i=0;i<Nt;i++) {
  201. Real cos_theta = CosTheta[i];
  202. Real sin_theta = sqrt(1-cos_theta*cos_theta);
  203. for (Long j=0;j<Np;j++) {
  204. Real cos_phi = cos(2*const_pi<Real>()*j/Np);
  205. Real sin_phi = sin(2*const_pi<Real>()*j/Np);
  206. Real qw = LWeights[i]*2*const_pi<Real>()/Np; // quadrature weights * area-element
  207. Real f[3]; // source density
  208. f[0] = Fgrid[(i*Np+j)*3+0];
  209. f[1] = Fgrid[(i*Np+j)*3+1];
  210. f[2] = Fgrid[(i*Np+j)*3+2];
  211. Real x[3]; // source coordinates
  212. x[0] = sin_theta*cos_phi;
  213. x[1] = sin_theta*sin_phi;
  214. x[2] = cos_theta;
  215. Real dr[3];
  216. dr[0] = trg[0] - x[0];
  217. dr[1] = trg[1] - x[1];
  218. dr[2] = trg[2] - x[2];
  219. Real invr = 1 / sqrt(dr[0]*dr[0] + dr[1]*dr[1] + dr[2]*dr[2]);
  220. Real invr2 = invr*invr;
  221. Real invr3 = invr2*invr;
  222. Real invr5 = invr2*invr3;
  223. Real fdotr = dr[0]*f[0]+dr[1]*f[1]+dr[2]*f[2];
  224. Real du[9];
  225. du[0] = ( fdotr*invr3 - 3*dr[0]*dr[0]*fdotr*invr5) * scal;
  226. du[1] = ((dr[0]*f[1]-dr[1]*f[0])*invr3 - 3*dr[0]*dr[1]*fdotr*invr5) * scal;
  227. du[2] = ((dr[0]*f[2]-dr[2]*f[0])*invr3 - 3*dr[0]*dr[2]*fdotr*invr5) * scal;
  228. du[3] = ((dr[1]*f[0]-dr[0]*f[1])*invr3 - 3*dr[1]*dr[0]*fdotr*invr5) * scal;
  229. du[4] = ( fdotr*invr3 - 3*dr[1]*dr[1]*fdotr*invr5) * scal;
  230. du[5] = ((dr[1]*f[2]-dr[2]*f[1])*invr3 - 3*dr[1]*dr[2]*fdotr*invr5) * scal;
  231. du[6] = ((dr[2]*f[0]-dr[0]*f[2])*invr3 - 3*dr[2]*dr[0]*fdotr*invr5) * scal;
  232. du[7] = ((dr[2]*f[1]-dr[1]*f[2])*invr3 - 3*dr[2]*dr[1]*fdotr*invr5) * scal;
  233. du[8] = ( fdotr*invr3 - 3*dr[2]*dr[2]*fdotr*invr5) * scal;
  234. Real p = 0;//(2*fdotr*invr3) * scal;
  235. Real K[9];
  236. K[0] = du[0] + du[0] - p; K[1] = du[1] + du[3] - 0; K[2] = du[2] + du[6] - 0;
  237. K[3] = du[3] + du[1] - 0; K[4] = du[4] + du[4] - p; K[5] = du[5] + du[7] - 0;
  238. K[6] = du[6] + du[2] - 0; K[7] = du[7] + du[5] - 0; K[8] = du[8] + du[8] - p;
  239. Sf[0] += (K[0]*nor[0] + K[1]*nor[1] + K[2]*nor[2]) * qw;
  240. Sf[1] += (K[3]*nor[0] + K[4]*nor[1] + K[5]*nor[2]) * qw;
  241. Sf[2] += (K[6]*nor[0] + K[7]*nor[1] + K[8]*nor[2]) * qw;
  242. }
  243. }
  244. };
  245. for (Long i = 0; i < 40; i++) { // Evaluate
  246. Real R0 = (0.01 + i/20.0);
  247. Vector<Real> x(3), n(3);
  248. x[0] = drand48()-0.5;
  249. x[1] = drand48()-0.5;
  250. x[2] = drand48()-0.5;
  251. n[0] = drand48()-0.5;
  252. n[1] = drand48()-0.5;
  253. n[2] = drand48()-0.5;
  254. Real R = sqrt<Real>(x[0]*x[0]+x[1]*x[1]+x[2]*x[2]);
  255. x[0] *= R0 / R;
  256. x[1] *= R0 / R;
  257. x[2] *= R0 / R;
  258. Vector<Real> Sf, Sf_;
  259. Vector<Real> Df, Df_;
  260. Vector<Real> Kf, Kf_;
  261. StokesEvalSL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, R0<1, Sf);
  262. StokesEvalDL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, R0<1, Df);
  263. StokesEvalKL(Fcoeff, sctl::SHCArrange::ROW_MAJOR, p, x, n, R0<1, Kf);
  264. stokes_evalSL(x, Sf_);
  265. stokes_evalDL(x, Df_);
  266. stokes_evalKL(x, n, Kf_);
  267. auto errSL = (Sf-Sf_)/(Sf+0.01);
  268. auto errDL = (Df-Df_)/(Df+0.01);
  269. auto errKL = (Kf-Kf_)/(Kf+0.01);
  270. for (auto& x:errSL) x=log(fabs(x))/log(10);
  271. for (auto& x:errDL) x=log(fabs(x))/log(10);
  272. for (auto& x:errKL) x=log(fabs(x))/log(10);
  273. std::cout<<"R = "<<(0.01 + i/20.0)<<"; SL-error = ";
  274. std::cout<<errSL;
  275. std::cout<<"R = "<<(0.01 + i/20.0)<<"; DL-error = ";
  276. std::cout<<errDL;
  277. std::cout<<"R = "<<(0.01 + i/20.0)<<"; KL-error = ";
  278. std::cout<<errKL;
  279. }
  280. Clear();
  281. }
  282. static void test() {
  283. int p = 3;
  284. int dof = 1;
  285. int Nt = p+1, Np = 2*p+1;
  286. auto print_coeff = [&](Vector<Real> S) {
  287. Long idx=0;
  288. for (Long k=0;k<dof;k++) {
  289. for (Long n=0;n<=p;n++) {
  290. std::cout<<Vector<Real>(2*n+2, S.begin()+idx);
  291. idx+=2*n+2;
  292. }
  293. }
  294. std::cout<<'\n';
  295. };
  296. Vector<Real> r_theta_phi, theta_phi;
  297. { // Set r_theta_phi, theta_phi
  298. Vector<Real> leg_nodes = LegendreNodes(Nt-1);
  299. for (Long i=0;i<Nt;i++) {
  300. for (Long j=0;j<Np;j++) {
  301. r_theta_phi.PushBack(1);
  302. r_theta_phi.PushBack(leg_nodes[i]);
  303. r_theta_phi.PushBack(j * 2 * const_pi<Real>() / Np);
  304. theta_phi.PushBack(leg_nodes[i]);
  305. theta_phi.PushBack(j * 2 * const_pi<Real>() / Np);
  306. }
  307. }
  308. }
  309. int Ncoeff = (p + 1) * (p + 1);
  310. Vector<Real> Xcoeff(dof * Ncoeff), Xgrid;
  311. for (int i=0;i<Xcoeff.Dim();i++) Xcoeff[i]=i+1;
  312. SHC2Grid(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, Nt, Np, &Xgrid);
  313. std::cout<<Matrix<Real>(Nt*dof, Np, Xgrid.begin())<<'\n';
  314. {
  315. Vector<Real> val;
  316. SHCEval(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, theta_phi, val);
  317. Matrix<Real>(dof, val.Dim()/dof, val.begin(), false) = Matrix<Real>(val.Dim()/dof, dof, val.begin()).Transpose();
  318. std::cout<<Matrix<Real>(val.Dim()/Np, Np, val.begin()) - Matrix<Real>(Nt*dof, Np, Xgrid.begin())+1e-10<<'\n';
  319. }
  320. Grid2SHC(Xgrid, Nt, Np, p, Xcoeff, sctl::SHCArrange::ROW_MAJOR);
  321. print_coeff(Xcoeff);
  322. //SphericalHarmonics<Real>::WriteVTK("test", nullptr, &Xcoeff, sctl::SHCArrange::ROW_MAJOR, p, 32);
  323. Clear();
  324. }
  325. /**
  326. * \brief Clear all precomputed data. This must be done before the program exits to avoid memory leaks.
  327. */
  328. static void Clear() { MatrixStore().Resize(0); }
  329. private:
  330. // Probably don't work anymore, need to be updated :(
  331. static void SHC2GridTranspose(const Vector<Real>& X, Long p0, Long p1, Vector<Real>& S);
  332. static void RotateAll(const Vector<Real>& S, Long p0, Long dof, Vector<Real>& S_);
  333. static void RotateTranspose(const Vector<Real>& S_, Long p0, Long dof, Vector<Real>& S);
  334. static void StokesSingularInteg(const Vector<Real>& S, Long p0, Long p1, Vector<Real>* SLMatrix=nullptr, Vector<Real>* DLMatrix=nullptr);
  335. static void Grid2SHC_(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& B1);
  336. static void SHCArrange0(const Vector<Real>& B1, Long p, Vector<Real>& S, SHCArrange arrange);
  337. 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);
  338. static void SHCArrange1(const Vector<Real>& S_in, SHCArrange arrange_out, Long p, Vector<Real>& S_out);
  339. /**
  340. * \brief Computes all the Associated Legendre Polynomials (normalized) up to the specified degree.
  341. * \param[in] degree The degree up to which the Legendre polynomials have to be computed.
  342. * \param[in] X The input values for which the polynomials have to be computed.
  343. * \param[in] N The number of input points.
  344. * \param[out] poly_val The output array of size (degree+1)*(degree+2)*N/2 containing the computed polynomial values.
  345. * The output values are in the order:
  346. * 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],
  347. * P(2,0)[0], ..., P(degree,0)[N-1], P(1,1)[0], ...,P(2,1)[0], ..., P(degree,degree)[N-1]}
  348. */
  349. static void LegPoly(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  350. static void LegPolyDeriv(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  351. static const Vector<Real>& LegendreNodes(Long p1);
  352. static const Vector<Real>& LegendreWeights(Long p1);
  353. static const Vector<Real>& SingularWeights(Long p1);
  354. static const Matrix<Real>& MatFourier(Long p0, Long p1);
  355. static const Matrix<Real>& MatFourierInv(Long p0, Long p1);
  356. static const Matrix<Real>& MatFourierGrad(Long p0, Long p1);
  357. static const FFT<Real>& OpFourier(Long Np);
  358. static const FFT<Real>& OpFourierInv(Long Np);
  359. static const std::vector<Matrix<Real>>& MatLegendre(Long p0, Long p1);
  360. static const std::vector<Matrix<Real>>& MatLegendreInv(Long p0, Long p1);
  361. static const std::vector<Matrix<Real>>& MatLegendreGrad(Long p0, Long p1);
  362. // Evaluate all Spherical Harmonic basis functions up to order p at (theta, phi) coordinates.
  363. static void SHBasisEval(Long p, const Vector<Real>& cos_theta_phi, Matrix<Real>& M);
  364. static void VecSHBasisEval(Long p, const Vector<Real>& cos_theta_phi, Matrix<Real>& M);
  365. static const std::vector<Matrix<Real>>& MatRotate(Long p0);
  366. template <bool SLayer, bool DLayer> static void StokesSingularInteg_(const Vector<Real>& X0, Long p0, Long p1, Vector<Real>& SL, Vector<Real>& DL);
  367. struct MatrixStorage{
  368. MatrixStorage(){
  369. const Long size = SCTL_SHMAXDEG;
  370. Resize(size);
  371. }
  372. void Resize(Long size){
  373. Qx_ .resize(size);
  374. Qw_ .resize(size);
  375. Sw_ .resize(size);
  376. Mf_ .resize(size*size);
  377. Mdf_.resize(size*size);
  378. Ml_ .resize(size*size);
  379. Mdl_.resize(size*size);
  380. Mr_ .resize(size);
  381. Mfinv_ .resize(size*size);
  382. Mlinv_ .resize(size*size);
  383. Mfft_.resize(size);
  384. Mfftinv_.resize(size);
  385. }
  386. std::vector<Vector<Real>> Qx_;
  387. std::vector<Vector<Real>> Qw_;
  388. std::vector<Vector<Real>> Sw_;
  389. std::vector<Matrix<Real>> Mf_ ;
  390. std::vector<Matrix<Real>> Mdf_;
  391. std::vector<std::vector<Matrix<Real>>> Ml_ ;
  392. std::vector<std::vector<Matrix<Real>>> Mdl_;
  393. std::vector<std::vector<Matrix<Real>>> Mr_;
  394. std::vector<Matrix<Real>> Mfinv_ ;
  395. std::vector<std::vector<Matrix<Real>>> Mlinv_ ;
  396. std::vector<FFT<Real>> Mfft_;
  397. std::vector<FFT<Real>> Mfftinv_;
  398. };
  399. static MatrixStorage& MatrixStore(){
  400. static MatrixStorage storage;
  401. return storage;
  402. }
  403. };
  404. template class SphericalHarmonics<double>;
  405. } // end namespace
  406. #include SCTL_INCLUDE(sph_harm.txx)
  407. #endif // _SCTL_SPH_HARM_HPP_