sph_harm.hpp 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208
  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. static void Grid2SHC(const Vector<Real>& X_in, Long Nt_in, Long Np_in, Long p_out, Vector<Real>& S_out, SHCArrange arrange_out);
  27. static void SHC2Grid(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Long Nt_out, Long Np_out, Vector<Real>* X_out, Vector<Real>* X_theta_out=nullptr, Vector<Real>* X_phi_out=nullptr);
  28. static void SHCEval(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& theta_phi_in, Vector<Real>& X_out);
  29. static void SHC2Pole(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Vector<Real>& P_out);
  30. 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());
  31. // Vector Spherical Harmonics
  32. static void Grid2VecSHC(const Vector<Real>& X_in, Long Nt_in, Long Np_in, Long p_out, Vector<Real>& S_out, SHCArrange arrange_out);
  33. static void VecSHC2Grid(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Long Nt_out, Long Np_out, Vector<Real>* X_out);
  34. static void VecSHCEval(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& theta_phi_in, Vector<Real>& X_out);
  35. static void StokesEvalSL(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& theta_phi_in, Vector<Real>& X_out);
  36. static void test() {
  37. int p = 4;
  38. int dof = 6;
  39. int Nt = p+1, Np = 2*p+1;
  40. auto print_coeff = [&](Vector<Real> S) {
  41. Long idx=0;
  42. for (Long k=0;k<dof;k++) {
  43. for (Long n=0;n<=p;n++) {
  44. std::cout<<Vector<Real>(2*n+2, S.begin()+idx);
  45. idx+=2*n+2;
  46. }
  47. }
  48. std::cout<<'\n';
  49. };
  50. Vector<Real> r_theta_phi, theta_phi;
  51. { // Set r_theta_phi, theta_phi
  52. Vector<Real> leg_nodes = LegendreNodes(Nt-1);
  53. for (Long i=0;i<Nt;i++) {
  54. for (Long j=0;j<Np;j++) {
  55. r_theta_phi.PushBack(1);
  56. r_theta_phi.PushBack(leg_nodes[i]);
  57. r_theta_phi.PushBack(j * 2 * const_pi<Real>() / Np);
  58. theta_phi.PushBack(leg_nodes[i]);
  59. theta_phi.PushBack(j * 2 * const_pi<Real>() / Np);
  60. }
  61. }
  62. }
  63. int Ncoeff = (p + 1) * (p + 1);
  64. Vector<Real> Xcoeff(dof * Ncoeff), Xgrid;
  65. for (int i=0;i<Xcoeff.Dim();i++) Xcoeff[i]=i+1;
  66. VecSHC2Grid(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, Nt, Np, &Xgrid);
  67. std::cout<<Matrix<Real>(Nt*dof, Np, Xgrid.begin())<<'\n';
  68. {
  69. Vector<Real> val;
  70. VecSHCEval(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, theta_phi, val);
  71. //StokesEvalSL(Xcoeff, sctl::SHCArrange::COL_MAJOR_NONZERO, p, r_theta_phi, val);
  72. Matrix<Real>(dof, val.Dim()/dof, val.begin(), false) = Matrix<Real>(val.Dim()/dof, dof, val.begin()).Transpose();
  73. std::cout<<Matrix<Real>(val.Dim()/Np, Np, val.begin()) - Matrix<Real>(Nt*dof, Np, Xgrid.begin())<<'\n';
  74. }
  75. Grid2VecSHC(Xgrid, Nt, Np, p, Xcoeff, sctl::SHCArrange::ROW_MAJOR);
  76. print_coeff(Xcoeff);
  77. //SphericalHarmonics<Real>::WriteVTK("test", nullptr, &Xcoeff, sctl::SHCArrange::ROW_MAJOR, p, 32);
  78. Clear();
  79. }
  80. static void Clear() { MatrixStore().Resize(0); }
  81. private:
  82. // Probably don't work anymore, need to be updated :(
  83. static void SHC2GridTranspose(const Vector<Real>& X, Long p0, Long p1, Vector<Real>& S);
  84. static void RotateAll(const Vector<Real>& S, Long p0, Long dof, Vector<Real>& S_);
  85. static void RotateTranspose(const Vector<Real>& S_, Long p0, Long dof, Vector<Real>& S);
  86. static void StokesSingularInteg(const Vector<Real>& S, Long p0, Long p1, Vector<Real>* SLMatrix=nullptr, Vector<Real>* DLMatrix=nullptr);
  87. static void Grid2SHC_(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& B1);
  88. static void SHCArrange0(const Vector<Real>& B1, Long p, Vector<Real>& S, SHCArrange arrange);
  89. 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);
  90. static void SHCArrange1(const Vector<Real>& S_in, SHCArrange arrange_out, Long p, Vector<Real>& S_out);
  91. /**
  92. * \brief Computes all the Associated Legendre Polynomials (normalized) up to the specified degree.
  93. * \param[in] degree The degree up to which the Legendre polynomials have to be computed.
  94. * \param[in] X The input values for which the polynomials have to be computed.
  95. * \param[in] N The number of input points.
  96. * \param[out] poly_val The output array of size (degree+1)*(degree+2)*N/2 containing the computed polynomial values.
  97. * The output values are in the order:
  98. * 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],
  99. * P(2,0)[0], ..., P(degree,0)[N-1], P(1,1)[0], ...,P(2,1)[0], ..., P(degree,degree)[N-1]}
  100. */
  101. static void LegPoly(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  102. static void LegPolyDeriv(Vector<Real>& poly_val, const Vector<Real>& X, Long degree);
  103. static const Vector<Real>& LegendreNodes(Long p1);
  104. static const Vector<Real>& LegendreWeights(Long p1);
  105. static const Vector<Real>& SingularWeights(Long p1);
  106. static const Matrix<Real>& MatFourier(Long p0, Long p1);
  107. static const Matrix<Real>& MatFourierInv(Long p0, Long p1);
  108. static const Matrix<Real>& MatFourierGrad(Long p0, Long p1);
  109. static const FFT<Real>& OpFourier(Long Np);
  110. static const FFT<Real>& OpFourierInv(Long Np);
  111. static const std::vector<Matrix<Real>>& MatLegendre(Long p0, Long p1);
  112. static const std::vector<Matrix<Real>>& MatLegendreInv(Long p0, Long p1);
  113. static const std::vector<Matrix<Real>>& MatLegendreGrad(Long p0, Long p1);
  114. // Evaluate all Spherical Harmonic basis functions up to order p at (theta, phi) coordinates.
  115. static void SHBasisEval(Long p, const Vector<Real>& cos_theta_phi, Matrix<Real>& M);
  116. static void VecSHBasisEval(Long p, const Vector<Real>& cos_theta_phi, Matrix<Real>& M);
  117. static const std::vector<Matrix<Real>>& MatRotate(Long p0);
  118. template <bool SLayer, bool DLayer> static void StokesSingularInteg_(const Vector<Real>& X0, Long p0, Long p1, Vector<Real>& SL, Vector<Real>& DL);
  119. struct MatrixStorage{
  120. MatrixStorage(){
  121. const Long size = SCTL_SHMAXDEG;
  122. Resize(size);
  123. }
  124. void Resize(Long size){
  125. Qx_ .resize(size);
  126. Qw_ .resize(size);
  127. Sw_ .resize(size);
  128. Mf_ .resize(size*size);
  129. Mdf_.resize(size*size);
  130. Ml_ .resize(size*size);
  131. Mdl_.resize(size*size);
  132. Mr_ .resize(size);
  133. Mfinv_ .resize(size*size);
  134. Mlinv_ .resize(size*size);
  135. Mfft_.resize(size);
  136. Mfftinv_.resize(size);
  137. }
  138. std::vector<Vector<Real>> Qx_;
  139. std::vector<Vector<Real>> Qw_;
  140. std::vector<Vector<Real>> Sw_;
  141. std::vector<Matrix<Real>> Mf_ ;
  142. std::vector<Matrix<Real>> Mdf_;
  143. std::vector<std::vector<Matrix<Real>>> Ml_ ;
  144. std::vector<std::vector<Matrix<Real>>> Mdl_;
  145. std::vector<std::vector<Matrix<Real>>> Mr_;
  146. std::vector<Matrix<Real>> Mfinv_ ;
  147. std::vector<std::vector<Matrix<Real>>> Mlinv_ ;
  148. std::vector<FFT<Real>> Mfft_;
  149. std::vector<FFT<Real>> Mfftinv_;
  150. };
  151. static MatrixStorage& MatrixStore(){
  152. static MatrixStorage storage;
  153. return storage;
  154. }
  155. };
  156. template class SphericalHarmonics<double>;
  157. } // end namespace
  158. #include SCTL_INCLUDE(sph_harm.txx)
  159. #endif // _SCTL_SPH_HARM_HPP_