sph_harm.txx 62 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789
  1. #include SCTL_INCLUDE(legendre_rule.hpp)
  2. // TODO: Replace work vectors with dynamic-arrays
  3. namespace SCTL_NAMESPACE {
  4. template <class Real> void SphericalHarmonics<Real>::Grid2SHC(const Vector<Real>& X, Long Nt, Long Np, Long p1, Vector<Real>& S, SHCArrange arrange){
  5. Long N = X.Dim() / (Np*Nt);
  6. assert(X.Dim() == N*Np*Nt);
  7. Vector<Real> B1(N*(p1+1)*(p1+1));
  8. Grid2SHC_(X, Nt, Np, p1, B1);
  9. B1 *= 1 / sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  10. SHCArrange0(B1, p1, S, arrange);
  11. }
  12. template <class Real> void SphericalHarmonics<Real>::Grid2SHC_(const Vector<Real>& X, Long Nt, Long Np, Long p1, Vector<Real>& B1){
  13. const auto& Mf = OpFourierInv(Np);
  14. assert(Mf.Dim(0) == Np);
  15. const std::vector<Matrix<Real>>& Ml = SphericalHarmonics<Real>::MatLegendreInv(Nt-1,p1);
  16. assert((Long)Ml.size() == p1+1);
  17. Long N = X.Dim() / (Np*Nt);
  18. assert(X.Dim() == N*Np*Nt);
  19. Vector<Real> B0((2*p1+1) * N*Nt);
  20. #pragma omp parallel
  21. { // B0 <-- Transpose(FFT(X))
  22. Integer tid=omp_get_thread_num();
  23. Integer omp_p=omp_get_num_threads();
  24. Long a=(tid+0)*N*Nt/omp_p;
  25. Long b=(tid+1)*N*Nt/omp_p;
  26. Vector<Real> buff(Mf.Dim(1));
  27. Long fft_coeff_len = std::min(buff.Dim(), 2*p1+2);
  28. Matrix<Real> B0_(2*p1+1, N*Nt, B0.begin(), false);
  29. const Matrix<Real> MX(N * Nt, Np, (Iterator<Real>)X.begin(), false);
  30. for (Long i = a; i < b; i++) {
  31. { // buff <-- FFT(Xi)
  32. const Vector<Real> Xi(Np, (Iterator<Real>)X.begin() + Np * i, false);
  33. Mf.Execute(Xi, buff);
  34. }
  35. { // B0 <-- Transpose(buff)
  36. B0_[0][i] = buff[0]; // skipping buff[1] == 0
  37. for (Long j = 2; j < fft_coeff_len; j++) B0_[j-1][i] = buff[j];
  38. for (Long j = fft_coeff_len; j < 2*p1+2; j++) B0_[j-1][i] = 0;
  39. }
  40. }
  41. }
  42. if (B1.Dim() != N*(p1+1)*(p1+1)) B1.ReInit(N*(p1+1)*(p1+1));
  43. #pragma omp parallel
  44. { // Evaluate Legendre polynomial
  45. Integer tid=omp_get_thread_num();
  46. Integer omp_p=omp_get_num_threads();
  47. Long offset0=0;
  48. Long offset1=0;
  49. for (Long i = 0; i < p1+1; i++) {
  50. Long N_ = (i==0 ? N : 2*N);
  51. Matrix<Real> Min (N_, Nt , B0.begin()+offset0, false);
  52. Matrix<Real> Mout(N_, p1+1-i, B1.begin()+offset1, false);
  53. { // Mout = Min * Ml[i] // split between threads
  54. Long a=(tid+0)*N_/omp_p;
  55. Long b=(tid+1)*N_/omp_p;
  56. if (a < b) {
  57. Matrix<Real> Min_ (b-a, Min .Dim(1), Min [a], false);
  58. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  59. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  60. }
  61. }
  62. offset0+=Min .Dim(0)*Min .Dim(1);
  63. offset1+=Mout.Dim(0)*Mout.Dim(1);
  64. }
  65. assert(offset0 == B0.Dim());
  66. assert(offset1 == B1.Dim());
  67. }
  68. }
  69. template <class Real> void SphericalHarmonics<Real>::SHCArrange0(const Vector<Real>& B1, Long p1, Vector<Real>& S, SHCArrange arrange){
  70. Long M = (p1+1)*(p1+1);
  71. Long N = B1.Dim() / M;
  72. assert(B1.Dim() == N*M);
  73. if (arrange == SHCArrange::ALL) { // S <-- Rearrange(B1)
  74. Long M = 2*(p1+1)*(p1+1);
  75. if(S.Dim() != N * M) S.ReInit(N * M);
  76. #pragma omp parallel
  77. { // S <-- Rearrange(B1)
  78. Integer tid=omp_get_thread_num();
  79. Integer omp_p=omp_get_num_threads();
  80. Long a=(tid+0)*N/omp_p;
  81. Long b=(tid+1)*N/omp_p;
  82. for (Long i = a; i < b; i++) {
  83. Long offset = 0;
  84. for (Long j = 0; j < p1+1; j++) {
  85. Long len = p1+1 - j;
  86. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  87. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  88. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 0;
  89. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = B_[k];
  90. offset += len;
  91. }
  92. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  93. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  94. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 1;
  95. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = B_[k];
  96. offset += len;
  97. } else {
  98. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 1;
  99. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = 0;
  100. }
  101. }
  102. }
  103. }
  104. }
  105. if (arrange == SHCArrange::ROW_MAJOR) { // S <-- Rearrange(B1)
  106. Long M = (p1+1)*(p1+2);
  107. if(S.Dim() != N * M) S.ReInit(N * M);
  108. #pragma omp parallel
  109. { // S <-- Rearrange(B1)
  110. Integer tid=omp_get_thread_num();
  111. Integer omp_p=omp_get_num_threads();
  112. Long a=(tid+0)*N/omp_p;
  113. Long b=(tid+1)*N/omp_p;
  114. for (Long i = a; i < b; i++) {
  115. Long offset = 0;
  116. for (Long j = 0; j < p1+1; j++) {
  117. Long len = p1+1 - j;
  118. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  119. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  120. Iterator<Real> S_ = S .begin() + i*M + 0;
  121. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = B_[k];
  122. offset += len;
  123. }
  124. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  125. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  126. Iterator<Real> S_ = S .begin() + i*M + 1;
  127. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = B_[k];
  128. offset += len;
  129. } else {
  130. Iterator<Real> S_ = S .begin() + i*M + 1;
  131. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = 0;
  132. }
  133. }
  134. }
  135. }
  136. }
  137. if (arrange == SHCArrange::COL_MAJOR_NONZERO) { // S <-- Rearrange(B1)
  138. Long M = (p1+1)*(p1+1);
  139. if(S.Dim() != N * M) S.ReInit(N * M);
  140. #pragma omp parallel
  141. { // S <-- Rearrange(B1)
  142. Integer tid=omp_get_thread_num();
  143. Integer omp_p=omp_get_num_threads();
  144. Long a=(tid+0)*N/omp_p;
  145. Long b=(tid+1)*N/omp_p;
  146. for (Long i = a; i < b; i++) {
  147. Long offset = 0;
  148. for (Long j = 0; j < p1+1; j++) {
  149. Long len = p1+1 - j;
  150. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  151. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  152. Iterator<Real> S_ = S .begin() + i*M + offset;
  153. for (Long k = 0; k < len; k++) S_[k] = B_[k];
  154. offset += len;
  155. }
  156. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  157. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  158. Iterator<Real> S_ = S .begin() + i*M + offset;
  159. for (Long k = 0; k < len; k++) S_[k] = B_[k];
  160. offset += len;
  161. }
  162. }
  163. }
  164. }
  165. }
  166. }
  167. template <class Real> void SphericalHarmonics<Real>::SHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p0, Long Nt, Long Np, Vector<Real>* X, Vector<Real>* X_phi, Vector<Real>* X_theta){
  168. Vector<Real> B0;
  169. SHCArrange1(S, arrange, p0, B0);
  170. B0 *= sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  171. SHC2Grid_(B0, p0, Nt, Np, X, X_phi, X_theta);
  172. }
  173. template <class Real> void SphericalHarmonics<Real>::SHCArrange1(const Vector<Real>& S, SHCArrange arrange, Long p0, Vector<Real>& B0){
  174. Long M, N;
  175. { // Set M, N
  176. M = 0;
  177. if (arrange == SHCArrange::ALL) M = 2*(p0+1)*(p0+1);
  178. if (arrange == SHCArrange::ROW_MAJOR) M = (p0+1)*(p0+2);
  179. if (arrange == SHCArrange::COL_MAJOR_NONZERO) M = (p0+1)*(p0+1);
  180. if (M == 0) return;
  181. N = S.Dim() / M;
  182. assert(S.Dim() == N * M);
  183. }
  184. if (B0.Dim() != N*(p0+1)*(p0+1)) B0.ReInit(N*(p0+1)*(p0+1));
  185. if (arrange == SHCArrange::ALL) { // B0 <-- Rearrange(S)
  186. #pragma omp parallel
  187. { // B0 <-- Rearrange(S)
  188. Integer tid=omp_get_thread_num();
  189. Integer omp_p=omp_get_num_threads();
  190. Long a=(tid+0)*N/omp_p;
  191. Long b=(tid+1)*N/omp_p;
  192. for (Long i = a; i < b; i++) {
  193. Long offset = 0;
  194. for (Long j = 0; j < p0+1; j++) {
  195. Long len = p0+1 - j;
  196. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  197. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  198. ConstIterator<Real> S_ = S .begin() + i*M + j*(p0+1)*2 + j*2 + 0;
  199. for (Long k = 0; k < len; k++) B_[k] = S_[k * (p0+1)*2];
  200. offset += len;
  201. }
  202. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  203. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  204. ConstIterator<Real> S_ = S .begin() + i*M + j*(p0+1)*2 + j*2 + 1;
  205. for (Long k = 0; k < len; k++) B_[k] = S_[k * (p0+1)*2];
  206. offset += len;
  207. }
  208. }
  209. }
  210. }
  211. }
  212. if (arrange == SHCArrange::ROW_MAJOR) { // B0 <-- Rearrange(S)
  213. #pragma omp parallel
  214. { // B0 <-- Rearrange(S)
  215. Integer tid=omp_get_thread_num();
  216. Integer omp_p=omp_get_num_threads();
  217. Long a=(tid+0)*N/omp_p;
  218. Long b=(tid+1)*N/omp_p;
  219. for (Long i = a; i < b; i++) {
  220. Long offset = 0;
  221. for (Long j = 0; j < p0+1; j++) {
  222. Long len = p0+1 - j;
  223. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  224. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  225. ConstIterator<Real> S_ = S .begin() + i*M + 0;
  226. for (Long k=0;k<len;k++) B_[k] = S_[(j+k)*(j+k+1) + 2*j];
  227. offset += len;
  228. }
  229. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  230. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  231. ConstIterator<Real> S_ = S .begin() + i*M + 1;
  232. for (Long k=0;k<len;k++) B_[k] = S_[(j+k)*(j+k+1) + 2*j];
  233. offset += len;
  234. }
  235. }
  236. }
  237. }
  238. }
  239. if (arrange == SHCArrange::COL_MAJOR_NONZERO) { // B0 <-- Rearrange(S)
  240. #pragma omp parallel
  241. { // B0 <-- Rearrange(S)
  242. Integer tid=omp_get_thread_num();
  243. Integer omp_p=omp_get_num_threads();
  244. Long a=(tid+0)*N/omp_p;
  245. Long b=(tid+1)*N/omp_p;
  246. for (Long i = a; i < b; i++) {
  247. Long offset = 0;
  248. for (Long j = 0; j < p0+1; j++) {
  249. Long len = p0+1 - j;
  250. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  251. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  252. ConstIterator<Real> S_ = S .begin() + i*M + offset;
  253. for (Long k = 0; k < len; k++) B_[k] = S_[k];
  254. offset += len;
  255. }
  256. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  257. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  258. ConstIterator<Real> S_ = S .begin() + i*M + offset;
  259. for (Long k = 0; k < len; k++) B_[k] = S_[k];
  260. offset += len;
  261. }
  262. }
  263. }
  264. }
  265. }
  266. }
  267. template <class Real> void SphericalHarmonics<Real>::SHC2Grid_(const Vector<Real>& B0, Long p0, Long Nt, Long Np, Vector<Real>* X, Vector<Real>* X_phi, Vector<Real>* X_theta){
  268. const auto& Mf = OpFourier(Np);
  269. assert(Mf.Dim(1) == Np);
  270. const std::vector<Matrix<Real>>& Ml =SphericalHarmonics<Real>::MatLegendre (p0,Nt-1);
  271. const std::vector<Matrix<Real>>& Mdl=SphericalHarmonics<Real>::MatLegendreGrad(p0,Nt-1);
  272. assert((Long)Ml .size() == p0+1);
  273. assert((Long)Mdl.size() == p0+1);
  274. Long N = B0.Dim() / ((p0+1)*(p0+1));
  275. assert(B0.Dim() == N*(p0+1)*(p0+1));
  276. if(X && X ->Dim()!=N*Np*Nt) X ->ReInit(N*Np*Nt);
  277. if(X_theta && X_theta->Dim()!=N*Np*Nt) X_theta->ReInit(N*Np*Nt);
  278. if(X_phi && X_phi ->Dim()!=N*Np*Nt) X_phi ->ReInit(N*Np*Nt);
  279. Vector<Real> B1(N*(2*p0+1)*Nt);
  280. if(X || X_phi){
  281. #pragma omp parallel
  282. { // Evaluate Legendre polynomial
  283. Integer tid=omp_get_thread_num();
  284. Integer omp_p=omp_get_num_threads();
  285. Long offset0=0;
  286. Long offset1=0;
  287. for(Long i=0;i<p0+1;i++){
  288. Long N_ = (i==0 ? N : 2*N);
  289. const Matrix<Real> Min (N_, p0+1-i, (Iterator<Real>)B0.begin()+offset0, false);
  290. Matrix<Real> Mout(N_, Nt , B1.begin()+offset1, false);
  291. { // Mout = Min * Ml[i] // split between threads
  292. Long a=(tid+0)*N_/omp_p;
  293. Long b=(tid+1)*N_/omp_p;
  294. if(a<b){
  295. const Matrix<Real> Min_ (b-a, Min .Dim(1), (Iterator<Real>)Min [a], false);
  296. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  297. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  298. }
  299. }
  300. offset0+=Min .Dim(0)*Min .Dim(1);
  301. offset1+=Mout.Dim(0)*Mout.Dim(1);
  302. }
  303. }
  304. #pragma omp parallel
  305. { // Transpose and evaluate Fourier
  306. Integer tid=omp_get_thread_num();
  307. Integer omp_p=omp_get_num_threads();
  308. Long a=(tid+0)*N*Nt/omp_p;
  309. Long b=(tid+1)*N*Nt/omp_p;
  310. Vector<Real> buff(Mf.Dim(0)); buff = 0;
  311. Long fft_coeff_len = std::min(buff.Dim(), 2*p0+2);
  312. Matrix<Real> B1_(2*p0+1, N*Nt, B1.begin(), false);
  313. for (Long i = a; i < b; i++) {
  314. { // buff <-- Transpose(B1)
  315. buff[0] = B1_[0][i];
  316. buff[1] = 0;
  317. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  318. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  319. }
  320. { // X <-- FFT(buff)
  321. Vector<Real> Xi(Np, X->begin() + Np * i, false);
  322. Mf.Execute(buff, Xi);
  323. }
  324. if(X_phi){ // Evaluate Fourier gradient
  325. { // buff <-- Transpose(B1)
  326. buff[0] = 0;
  327. buff[1] = 0;
  328. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  329. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  330. for (Long j = 1; j < buff.Dim()/2; j++) {
  331. Real x = buff[2*j+0];
  332. Real y = buff[2*j+1];
  333. buff[2*j+0] = -j*y;
  334. buff[2*j+1] = j*x;
  335. }
  336. }
  337. { // X_phi <-- FFT(buff)
  338. Vector<Real> Xi(Np, X_phi->begin() + Np * i, false);
  339. Mf.Execute(buff, Xi);
  340. }
  341. }
  342. }
  343. }
  344. }
  345. if(X_theta){
  346. #pragma omp parallel
  347. { // Evaluate Legendre gradient
  348. Integer tid=omp_get_thread_num();
  349. Integer omp_p=omp_get_num_threads();
  350. Long offset0=0;
  351. Long offset1=0;
  352. for(Long i=0;i<p0+1;i++){
  353. Long N_ = (i==0 ? N : 2*N);
  354. const Matrix<Real> Min (N_, p0+1-i, (Iterator<Real>)B0.begin()+offset0, false);
  355. Matrix<Real> Mout(N_, Nt , B1.begin()+offset1, false);
  356. { // Mout = Min * Mdl[i] // split between threads
  357. Long a=(tid+0)*N_/omp_p;
  358. Long b=(tid+1)*N_/omp_p;
  359. if(a<b){
  360. const Matrix<Real> Min_ (b-a, Min .Dim(1), (Iterator<Real>)Min [a], false);
  361. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  362. Matrix<Real>::GEMM(Mout_,Min_,Mdl[i]);
  363. }
  364. }
  365. offset0+=Min .Dim(0)*Min .Dim(1);
  366. offset1+=Mout.Dim(0)*Mout.Dim(1);
  367. }
  368. }
  369. #pragma omp parallel
  370. { // Transpose and evaluate Fourier
  371. Integer tid=omp_get_thread_num();
  372. Integer omp_p=omp_get_num_threads();
  373. Long a=(tid+0)*N*Nt/omp_p;
  374. Long b=(tid+1)*N*Nt/omp_p;
  375. Vector<Real> buff(Mf.Dim(0)); buff = 0;
  376. Long fft_coeff_len = std::min(buff.Dim(), 2*p0+2);
  377. Matrix<Real> B1_(2*p0+1, N*Nt, B1.begin(), false);
  378. for (Long i = a; i < b; i++) {
  379. { // buff <-- Transpose(B1)
  380. buff[0] = B1_[0][i];
  381. buff[1] = 0;
  382. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  383. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  384. }
  385. { // Xi <-- FFT(buff)
  386. Vector<Real> Xi(Np, X_theta->begin() + Np * i, false);
  387. Mf.Execute(buff, Xi);
  388. }
  389. }
  390. }
  391. }
  392. }
  393. template <class Real> void SphericalHarmonics<Real>::SHC2Pole(const Vector<Real>& S, SHCArrange arrange, Long p0, Vector<Real>& P){
  394. Vector<Real> QP[2];
  395. { // Set QP // TODO: store these weights
  396. Vector<Real> x(1), alp;
  397. const Real SQRT2PI = sqrt<Real>(4 * const_pi<Real>());
  398. for (Long i = 0; i < 2; i++) {
  399. x = (i ? -1 : 1);
  400. LegPoly(alp, x, p0);
  401. QP[i].ReInit(p0 + 1, alp.begin());
  402. QP[i] *= SQRT2PI;
  403. }
  404. }
  405. Long M, N;
  406. { // Set M, N
  407. M = 0;
  408. if (arrange == SHCArrange::ALL) M = 2*(p0+1)*(p0+1);
  409. if (arrange == SHCArrange::ROW_MAJOR) M = (p0+1)*(p0+2);
  410. if (arrange == SHCArrange::COL_MAJOR_NONZERO) M = (p0+1)*(p0+1);
  411. if (M == 0) return;
  412. N = S.Dim() / M;
  413. assert(S.Dim() == N * M);
  414. }
  415. if(P.Dim() != N * 2) P.ReInit(N * 2);
  416. if (arrange == SHCArrange::ALL) {
  417. #pragma omp parallel
  418. { // Compute pole
  419. Integer tid = omp_get_thread_num();
  420. Integer omp_p = omp_get_num_threads();
  421. Long a = (tid + 0) * N / omp_p;
  422. Long b = (tid + 1) * N / omp_p;
  423. for (Long i = a; i < b; i++) {
  424. Real P_[2] = {0, 0};
  425. for (Long j = 0; j < p0 + 1; j++) {
  426. P_[0] += S[i*M + j*(p0+1)*2] * QP[0][j];
  427. P_[1] += S[i*M + j*(p0+1)*2] * QP[1][j];
  428. }
  429. P[2*i+0] = P_[0];
  430. P[2*i+1] = P_[1];
  431. }
  432. }
  433. }
  434. if (arrange == SHCArrange::ROW_MAJOR) {
  435. #pragma omp parallel
  436. { // Compute pole
  437. Integer tid = omp_get_thread_num();
  438. Integer omp_p = omp_get_num_threads();
  439. Long a = (tid + 0) * N / omp_p;
  440. Long b = (tid + 1) * N / omp_p;
  441. for (Long i = a; i < b; i++) {
  442. Long idx = 0;
  443. Real P_[2] = {0, 0};
  444. for (Long j = 0; j < p0 + 1; j++) {
  445. P_[0] += S[i*M+idx] * QP[0][j];
  446. P_[1] += S[i*M+idx] * QP[1][j];
  447. idx += 2*(j+1);
  448. }
  449. P[2*i+0] = P_[0];
  450. P[2*i+1] = P_[1];
  451. }
  452. }
  453. }
  454. if (arrange == SHCArrange::COL_MAJOR_NONZERO) {
  455. #pragma omp parallel
  456. { // Compute pole
  457. Integer tid = omp_get_thread_num();
  458. Integer omp_p = omp_get_num_threads();
  459. Long a = (tid + 0) * N / omp_p;
  460. Long b = (tid + 1) * N / omp_p;
  461. for (Long i = a; i < b; i++) {
  462. Real P_[2] = {0, 0};
  463. for (Long j = 0; j < p0 + 1; j++) {
  464. P_[0] += S[i*M+j] * QP[0][j];
  465. P_[1] += S[i*M+j] * QP[1][j];
  466. }
  467. P[2*i+0] = P_[0];
  468. P[2*i+1] = P_[1];
  469. }
  470. }
  471. }
  472. }
  473. template <class Real> void SphericalHarmonics<Real>::WriteVTK(const char* fname, const Vector<Real>* S, const Vector<Real>* v_ptr, SHCArrange arrange, Long p0, Long p1, Real period, const Comm& comm){
  474. typedef double VTKReal;
  475. Vector<Real> SS;
  476. if (S == nullptr) {
  477. Integer p = 2;
  478. Integer Ncoeff = (p + 1) * (p + 1);
  479. Vector<Real> SSS(COORD_DIM * Ncoeff), SSS_grid;
  480. SSS.SetZero();
  481. SSS[1+0*p+0*Ncoeff] = sqrt<Real>(2.0)/sqrt<Real>(3.0);
  482. SSS[1+1*p+1*Ncoeff] = 1/sqrt<Real>(3.0);
  483. SSS[1+2*p+2*Ncoeff] = 1/sqrt<Real>(3.0);
  484. SphericalHarmonics<Real>::SHC2Grid(SSS, SHCArrange::COL_MAJOR_NONZERO, p, p+1, 2*p+2, &SSS_grid);
  485. SphericalHarmonics<Real>::Grid2SHC(SSS_grid, p+1, 2*p+2, p0, SS, arrange);
  486. S = &SS;
  487. }
  488. Vector<Real> X, Xp, V, Vp;
  489. { // Upsample X
  490. const Vector<Real>& X0=*S;
  491. SphericalHarmonics<Real>::SHC2Grid(X0, arrange, p0, p1+1, 2*p1, &X);
  492. SphericalHarmonics<Real>::SHC2Pole(X0, arrange, p0, Xp);
  493. }
  494. if(v_ptr){ // Upsample V
  495. const Vector<Real>& X0=*v_ptr;
  496. SphericalHarmonics<Real>::SHC2Grid(X0, arrange, p0, p1+1, 2*p1, &V);
  497. SphericalHarmonics<Real>::SHC2Pole(X0, arrange, p0, Vp);
  498. }
  499. std::vector<VTKReal> point_coord;
  500. std::vector<VTKReal> point_value;
  501. std::vector<int32_t> poly_connect;
  502. std::vector<int32_t> poly_offset;
  503. { // Set point_coord, point_value, poly_connect
  504. Long N_ves = X.Dim()/(2*p1*(p1+1)*COORD_DIM); // Number of vesicles
  505. assert(Xp.Dim() == N_ves*2*COORD_DIM);
  506. for(Long k=0;k<N_ves;k++){ // Set point_coord
  507. Real C[COORD_DIM]={0,0,0};
  508. if(period>0){
  509. for(Integer l=0;l<COORD_DIM;l++) C[l]=0;
  510. for(Long i=0;i<p1+1;i++){
  511. for(Long j=0;j<2*p1;j++){
  512. for(Integer l=0;l<COORD_DIM;l++){
  513. C[l]+=X[j+2*p1*(i+(p1+1)*(l+k*COORD_DIM))];
  514. }
  515. }
  516. }
  517. for(Integer l=0;l<COORD_DIM;l++) C[l]+=Xp[0+2*(l+k*COORD_DIM)];
  518. for(Integer l=0;l<COORD_DIM;l++) C[l]+=Xp[1+2*(l+k*COORD_DIM)];
  519. for(Integer l=0;l<COORD_DIM;l++) C[l]/=2*p1*(p1+1)+2;
  520. for(Integer l=0;l<COORD_DIM;l++) C[l]=(round(C[l]/period))*period;
  521. }
  522. for(Long i=0;i<p1+1;i++){
  523. for(Long j=0;j<2*p1;j++){
  524. for(Integer l=0;l<COORD_DIM;l++){
  525. point_coord.push_back(X[j+2*p1*(i+(p1+1)*(l+k*COORD_DIM))]-C[l]);
  526. }
  527. }
  528. }
  529. for(Integer l=0;l<COORD_DIM;l++) point_coord.push_back(Xp[0+2*(l+k*COORD_DIM)]-C[l]);
  530. for(Integer l=0;l<COORD_DIM;l++) point_coord.push_back(Xp[1+2*(l+k*COORD_DIM)]-C[l]);
  531. }
  532. if(v_ptr) {
  533. Long data__dof = V.Dim() / (2*p1*(p1+1));
  534. for(Long k=0;k<N_ves;k++){ // Set point_value
  535. for(Long i=0;i<p1+1;i++){
  536. for(Long j=0;j<2*p1;j++){
  537. for(Long l=0;l<data__dof;l++){
  538. point_value.push_back(V[j+2*p1*(i+(p1+1)*(l+k*data__dof))]);
  539. }
  540. }
  541. }
  542. for(Long l=0;l<data__dof;l++) point_value.push_back(Vp[0+2*(l+k*data__dof)]);
  543. for(Long l=0;l<data__dof;l++) point_value.push_back(Vp[1+2*(l+k*data__dof)]);
  544. }
  545. }
  546. for(Long k=0;k<N_ves;k++){
  547. for(Long j=0;j<2*p1;j++){
  548. Long i0= 0;
  549. Long i1=p1;
  550. Long j0=((j+0) );
  551. Long j1=((j+1)%(2*p1));
  552. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*(p1+1)+0);
  553. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j0);
  554. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j1);
  555. poly_offset.push_back(poly_connect.size());
  556. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*(p1+1)+1);
  557. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j0);
  558. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j1);
  559. poly_offset.push_back(poly_connect.size());
  560. }
  561. for(Long i=0;i<p1;i++){
  562. for(Long j=0;j<2*p1;j++){
  563. Long i0=((i+0) );
  564. Long i1=((i+1) );
  565. Long j0=((j+0) );
  566. Long j1=((j+1)%(2*p1));
  567. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j0);
  568. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j0);
  569. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j1);
  570. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j1);
  571. poly_offset.push_back(poly_connect.size());
  572. }
  573. }
  574. }
  575. }
  576. Integer np = comm.Size();
  577. Integer myrank = comm.Rank();
  578. std::vector<VTKReal>& coord=point_coord;
  579. std::vector<VTKReal>& value=point_value;
  580. std::vector<int32_t>& connect=poly_connect;
  581. std::vector<int32_t>& offset=poly_offset;
  582. Long pt_cnt=coord.size()/COORD_DIM;
  583. Long poly_cnt=poly_offset.size();
  584. // Open file for writing.
  585. std::stringstream vtufname;
  586. vtufname<<fname<<"_"<<std::setfill('0')<<std::setw(6)<<myrank<<".vtp";
  587. std::ofstream vtufile;
  588. vtufile.open(vtufname.str().c_str());
  589. if(vtufile.fail()) return;
  590. bool isLittleEndian;
  591. { // Set isLittleEndian
  592. uint16_t number = 0x1;
  593. uint8_t *numPtr = (uint8_t*)&number;
  594. isLittleEndian=(numPtr[0] == 1);
  595. }
  596. // Proceed to write to file.
  597. Long data_size=0;
  598. vtufile<<"<?xml version=\"1.0\"?>\n";
  599. if(isLittleEndian) vtufile<<"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">\n";
  600. else vtufile<<"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"BigEndian\">\n";
  601. //===========================================================================
  602. vtufile<<" <PolyData>\n";
  603. vtufile<<" <Piece NumberOfPoints=\""<<pt_cnt<<"\" NumberOfVerts=\"0\" NumberOfLines=\"0\" NumberOfStrips=\"0\" NumberOfPolys=\""<<poly_cnt<<"\">\n";
  604. //---------------------------------------------------------------------------
  605. vtufile<<" <Points>\n";
  606. vtufile<<" <DataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<COORD_DIM<<"\" Name=\"Position\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  607. data_size+=sizeof(uint32_t)+coord.size()*sizeof(VTKReal);
  608. vtufile<<" </Points>\n";
  609. //---------------------------------------------------------------------------
  610. if(value.size()){ // value
  611. vtufile<<" <PointData>\n";
  612. vtufile<<" <DataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<value.size()/pt_cnt<<"\" Name=\""<<"value"<<"\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  613. data_size+=sizeof(uint32_t)+value.size()*sizeof(VTKReal);
  614. vtufile<<" </PointData>\n";
  615. }
  616. //---------------------------------------------------------------------------
  617. vtufile<<" <Polys>\n";
  618. vtufile<<" <DataArray type=\"Int32\" Name=\"connectivity\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  619. data_size+=sizeof(uint32_t)+connect.size()*sizeof(int32_t);
  620. vtufile<<" <DataArray type=\"Int32\" Name=\"offsets\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  621. data_size+=sizeof(uint32_t)+offset.size() *sizeof(int32_t);
  622. vtufile<<" </Polys>\n";
  623. //---------------------------------------------------------------------------
  624. vtufile<<" </Piece>\n";
  625. vtufile<<" </PolyData>\n";
  626. //===========================================================================
  627. vtufile<<" <AppendedData encoding=\"raw\">\n";
  628. vtufile<<" _";
  629. int32_t block_size;
  630. block_size=coord.size()*sizeof(VTKReal); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&coord [0], coord.size()*sizeof(VTKReal));
  631. if(value.size()){ // value
  632. block_size=value.size()*sizeof(VTKReal); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&value [0], value.size()*sizeof(VTKReal));
  633. }
  634. block_size=connect.size()*sizeof(int32_t); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&connect[0], connect.size()*sizeof(int32_t));
  635. block_size=offset .size()*sizeof(int32_t); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&offset [0], offset .size()*sizeof(int32_t));
  636. vtufile<<"\n";
  637. vtufile<<" </AppendedData>\n";
  638. //===========================================================================
  639. vtufile<<"</VTKFile>\n";
  640. vtufile.close();
  641. if(myrank) return;
  642. std::stringstream pvtufname;
  643. pvtufname<<fname<<".pvtp";
  644. std::ofstream pvtufile;
  645. pvtufile.open(pvtufname.str().c_str());
  646. if(pvtufile.fail()) return;
  647. pvtufile<<"<?xml version=\"1.0\"?>\n";
  648. pvtufile<<"<VTKFile type=\"PPolyData\">\n";
  649. pvtufile<<" <PPolyData GhostLevel=\"0\">\n";
  650. pvtufile<<" <PPoints>\n";
  651. pvtufile<<" <PDataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<COORD_DIM<<"\" Name=\"Position\"/>\n";
  652. pvtufile<<" </PPoints>\n";
  653. if(value.size()){ // value
  654. pvtufile<<" <PPointData>\n";
  655. pvtufile<<" <PDataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<value.size()/pt_cnt<<"\" Name=\""<<"value"<<"\"/>\n";
  656. pvtufile<<" </PPointData>\n";
  657. }
  658. {
  659. // Extract filename from path.
  660. std::stringstream vtupath;
  661. vtupath<<'/'<<fname;
  662. std::string pathname = vtupath.str();
  663. auto found = pathname.find_last_of("/\\");
  664. std::string fname_ = pathname.substr(found+1);
  665. for(Integer i=0;i<np;i++) pvtufile<<" <Piece Source=\""<<fname_<<"_"<<std::setfill('0')<<std::setw(6)<<i<<".vtp\"/>\n";
  666. }
  667. pvtufile<<" </PPolyData>\n";
  668. pvtufile<<"</VTKFile>\n";
  669. pvtufile.close();
  670. }
  671. template <class Real> void SphericalHarmonics<Real>::Grid2VecSHC(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& S, SHCArrange arrange) {
  672. Long N = X.Dim() / (Np*Nt);
  673. assert(X.Dim() == N*Np*Nt);
  674. assert(N % COORD_DIM == 0);
  675. Vector<Real> B0(N*Nt*Np);
  676. { // Set B0
  677. B0 = X;
  678. const auto& Y = LegendreNodes(Nt - 1);
  679. assert(Y.Dim() == Nt);
  680. for (Long k = 0; k < N; k++) {
  681. if (k % COORD_DIM) {
  682. for (Long i = 0; i < Nt; i++) {
  683. Real s = 1/sqrt<Real>(1 - Y[i]*Y[i]);
  684. for (Long j = 0; j < Np; j++) {
  685. B0[(k*Nt+i)*Np+j] *= s;
  686. }
  687. }
  688. }
  689. }
  690. }
  691. Long M = (p+1)*(p+1);
  692. Vector<Real> B1(N*M);
  693. Grid2SHC_(B0, Nt, Np, p, B1);
  694. B1 *= 1 / sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  695. Vector<Real> B2(N*M);
  696. const Complex<Real> imag(0,1);
  697. for (Long i=0; i<N; i+=COORD_DIM) {
  698. for (Long m=0; m<=p; m++) {
  699. for (Long n=m; n<=p; n++) {
  700. auto read_coeff = [&](const Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  701. Complex<Real> c;
  702. if (0<=m && m<=n && n<=p) {
  703. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  704. Long idx_imag = idx_real + (p+1-m)*N;
  705. c.real = coeff[idx_real];
  706. if (m) c.imag = coeff[idx_imag];
  707. }
  708. return c;
  709. };
  710. auto write_coeff = [&](Complex<Real> c, Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  711. if (0<=m && m<=n && n<=p) {
  712. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  713. Long idx_imag = idx_real + (p+1-m)*N;
  714. coeff[idx_real] = c.real;
  715. if (m) coeff[idx_imag] = c.imag;
  716. }
  717. };
  718. auto gr = [&](Long n, Long m) { return read_coeff(B1, i+0, p, n, m); };
  719. auto gt = [&](Long n, Long m) { return read_coeff(B1, i+1, p, n, m); };
  720. auto gp = [&](Long n, Long m) { return read_coeff(B1, i+2, p, n, m); };
  721. Complex<Real> phiY, phiG, phiX;
  722. { // (phiG, phiX) <-- (gt, gp)
  723. auto A = [&](Long n, Long m) { return (0<=n && m<=n && n<=p ? sqrt<Real>(n*n * ((n+1)*(n+1) - m*m) / (Real)((2*n+1)*(2*n+3))) : 0); };
  724. auto B = [&](Long n, Long m) { return (0<=n && m<=n && n<=p ? sqrt<Real>((n+1)*(n+1) * (n*n - m*m) / (Real)((2*n+1)*(2*n-1))) : 0); };
  725. phiY = gr(n,m);
  726. phiG = (gt(n+1,m)*A(n,m) - gt(n-1,m)*B(n,m) - imag*m*gp(n,m)) * (1/(Real)(std::max<Long>(n,1)*(n+1)));
  727. phiX = (gp(n+1,m)*A(n,m) - gp(n-1,m)*B(n,m) + imag*m*gt(n,m)) * (1/(Real)(std::max<Long>(n,1)*(n+1)));
  728. }
  729. auto phiV = (phiG * (n + 0) - phiY) * (1/(Real)(2*n + 1));
  730. auto phiW = (phiG * (n + 1) + phiY) * (1/(Real)(2*n + 1));
  731. write_coeff(phiV, B2, i+0, p, n, m);
  732. write_coeff(phiW, B2, i+1, p, n, m);
  733. write_coeff(phiX, B2, i+2, p, n, m);
  734. }
  735. }
  736. }
  737. SHCArrange0(B2, p, S, arrange);
  738. }
  739. template <class Real> void SphericalHarmonics<Real>::VecSHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p, Long Nt, Long Np, Vector<Real>& X) {
  740. Vector<Real> B0;
  741. SHCArrange1(S, arrange, p, B0);
  742. Long M = (p+1)*(p+1);
  743. Long N = B0.Dim() / M;
  744. assert(B0.Dim() == N*M);
  745. assert(N % COORD_DIM == 0);
  746. Vector<Real> B1(N*M);
  747. const Complex<Real> imag(0,1);
  748. for (Long i=0; i<N; i+=COORD_DIM) {
  749. for (Long n=0; n<=p; n++) {
  750. for (Long m=0; m<=n; m++) {
  751. auto read_coeff = [&](const Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  752. Complex<Real> c;
  753. if (0<=m && m<=n && n<=p) {
  754. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  755. Long idx_imag = idx_real + (p+1-m)*N;
  756. c.real = coeff[idx_real];
  757. if (m) c.imag = coeff[idx_imag];
  758. }
  759. return c;
  760. };
  761. auto write_coeff = [&](Complex<Real> c, Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  762. if (0<=m && m<=n && n<=p) {
  763. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  764. Long idx_imag = idx_real + (p+1-m)*N;
  765. coeff[idx_real] = c.real;
  766. if (m) coeff[idx_imag] = c.imag;
  767. }
  768. };
  769. auto phiG = [&](Long n, Long m) {
  770. auto phiV = read_coeff(B0, i+0, p, n, m);
  771. auto phiW = read_coeff(B0, i+1, p, n, m);
  772. return phiV + phiW;
  773. };
  774. auto phiY = [&](Long n, Long m) {
  775. auto phiV = read_coeff(B0, i+0, p, n, m);
  776. auto phiW = read_coeff(B0, i+1, p, n, m);
  777. return -phiV * (n + 1) + phiW * n;
  778. };
  779. auto phiX = [&](Long n, Long m) {
  780. return read_coeff(B0, i+2, p, n, m);
  781. };
  782. Complex<Real> gr, gt, gp;
  783. { // (gt, gp) <-- (phiG, phiX)
  784. auto A = [&](Long n, Long m) { return (0<=n && m<=n && n<=p ? sqrt<Real>(n*n * ((n+1)*(n+1) - m*m) / (Real)((2*n+1)*(2*n+3))) : 0); };
  785. auto B = [&](Long n, Long m) { return (0<=n && m<=n && n<=p ? sqrt<Real>((n+1)*(n+1) * (n*n - m*m) / (Real)((2*n+1)*(2*n-1))) : 0); };
  786. gr = phiY(n,m);
  787. gt = phiG(n-1,m)*A(n-1,m) - phiG(n+1,m)*B(n+1,m) - imag*m*phiX(n,m);
  788. gp = phiX(n-1,m)*A(n-1,m) - phiX(n+1,m)*B(n+1,m) + imag*m*phiG(n,m);
  789. }
  790. write_coeff(gr, B1, i+0, p, n, m);
  791. write_coeff(gt, B1, i+1, p, n, m);
  792. write_coeff(gp, B1, i+2, p, n, m);
  793. }
  794. }
  795. }
  796. B1 *= sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  797. SHC2Grid_(B1, p, Nt, Np, &X);
  798. { // Set X
  799. const auto& Y = LegendreNodes(Nt - 1);
  800. assert(Y.Dim() == Nt);
  801. for (Long k = 0; k < N; k++) {
  802. if (k % COORD_DIM) {
  803. for (Long i = 0; i < Nt; i++) {
  804. Real s = 1/sqrt<Real>(1 - Y[i]*Y[i]);
  805. for (Long j = 0; j < Np; j++) {
  806. X[(k*Nt+i)*Np+j] *= s;
  807. }
  808. }
  809. }
  810. }
  811. }
  812. }
  813. template <class Real> void SphericalHarmonics<Real>::LegPoly(Vector<Real>& poly_val, const Vector<Real>& X, Long degree){
  814. Long N = X.Dim();
  815. Long Npoly = (degree + 1) * (degree + 2) / 2;
  816. if (poly_val.Dim() != Npoly * N) poly_val.ReInit(Npoly * N);
  817. Real fact = 1 / sqrt<Real>(4 * const_pi<Real>());
  818. Vector<Real> u(N);
  819. for (Long n = 0; n < N; n++) {
  820. u[n] = (X[n]*X[n]<1 ? sqrt<Real>(1-X[n]*X[n]) : 0);
  821. poly_val[n] = fact;
  822. }
  823. Long idx = 0;
  824. Long idx_nxt = 0;
  825. for (Long i = 1; i <= degree; i++) {
  826. idx_nxt += N*(degree-i+2);
  827. Real c = sqrt<Real>((2*i+1)/(Real)(2*i));
  828. for (Long n = 0; n < N; n++) {
  829. poly_val[idx_nxt+n] = -poly_val[idx+n] * u[n] * c;
  830. }
  831. idx = idx_nxt;
  832. }
  833. idx = 0;
  834. for (Long m = 0; m < degree; m++) {
  835. for (Long n = 0; n < N; n++) {
  836. Real pmm = 0;
  837. Real pmmp1 = poly_val[idx+n];
  838. for (Long ll = m + 1; ll <= degree; ll++) {
  839. Real a = sqrt<Real>(((2*ll-1)*(2*ll+1) ) / (Real)((ll-m)*(ll+m) ));
  840. Real b = sqrt<Real>(((2*ll+1)*(ll+m-1)*(ll-m-1)) / (Real)((ll-m)*(ll+m)*(2*ll-3)));
  841. Real pll = X[n]*a*pmmp1 - b*pmm;
  842. pmm = pmmp1;
  843. pmmp1 = pll;
  844. poly_val[idx + N*(ll-m) + n] = pll;
  845. }
  846. }
  847. idx += N * (degree - m + 1);
  848. }
  849. }
  850. template <class Real> void SphericalHarmonics<Real>::LegPolyDeriv(Vector<Real>& poly_val, const Vector<Real>& X, Long degree){
  851. Long N = X.Dim();
  852. Long Npoly = (degree + 1) * (degree + 2) / 2;
  853. if (poly_val.Dim() != N * Npoly) poly_val.ReInit(N * Npoly);
  854. Vector<Real> leg_poly(Npoly * N);
  855. LegPoly(leg_poly, X, degree);
  856. for (Long m = 0; m <= degree; m++) {
  857. for (Long n = m; n <= degree; n++) {
  858. ConstIterator<Real> Pn = leg_poly.begin() + N * ((degree * 2 - m + 1) * (m + 0) / 2 + n);
  859. ConstIterator<Real> Pn_ = leg_poly.begin() + N * ((degree * 2 - m + 0) * (m + 1) / 2 + n) * (m < n);
  860. Iterator <Real> Hn = poly_val.begin() + N * ((degree * 2 - m + 1) * (m + 0) / 2 + n);
  861. Real c2 = sqrt<Real>(m<n ? (n+m+1)*(n-m) : 0);
  862. for (Long i = 0; i < N; i++) {
  863. Real c1 = (X[i]*X[i]<1 ? m/sqrt<Real>(1-X[i]*X[i]) : 0);
  864. Hn[i] = -c1*X[i]*Pn[i] - c2*Pn_[i];
  865. }
  866. }
  867. }
  868. }
  869. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::LegendreNodes(Long p){
  870. assert(p<SCTL_SHMAXDEG);
  871. Vector<Real>& Qx=MatrixStore().Qx_[p];
  872. if(!Qx.Dim()){
  873. Vector<double> qx1(p+1);
  874. Vector<double> qw1(p+1);
  875. cgqf(p+1, 1, 0.0, 0.0, -1.0, 1.0, &qx1[0], &qw1[0]);
  876. assert(typeid(Real) == typeid(double) || typeid(Real) == typeid(float)); // TODO: works only for float and double
  877. if (Qx.Dim() != p+1) Qx.ReInit(p+1);
  878. for (Long i = 0; i < p + 1; i++) Qx[i] = -qx1[i];
  879. }
  880. return Qx;
  881. }
  882. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::LegendreWeights(Long p){
  883. assert(p<SCTL_SHMAXDEG);
  884. Vector<Real>& Qw=MatrixStore().Qw_[p];
  885. if(!Qw.Dim()){
  886. Vector<double> qx1(p+1);
  887. Vector<double> qw1(p+1);
  888. cgqf(p+1, 1, 0.0, 0.0, -1.0, 1.0, &qx1[0], &qw1[0]);
  889. assert(typeid(Real) == typeid(double) || typeid(Real) == typeid(float)); // TODO: works only for float and double
  890. if (Qw.Dim() != p+1) Qw.ReInit(p+1);
  891. for (Long i = 0; i < p + 1; i++) Qw[i] = qw1[i];
  892. }
  893. return Qw;
  894. }
  895. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::SingularWeights(Long p1){
  896. assert(p1<SCTL_SHMAXDEG);
  897. Vector<Real>& Sw=MatrixStore().Sw_[p1];
  898. if(!Sw.Dim()){
  899. const Vector<Real>& qx1 = LegendreNodes(p1);
  900. const Vector<Real>& qw1 = LegendreWeights(p1);
  901. std::vector<Real> Yf(p1+1,0);
  902. { // Set Yf
  903. Vector<Real> x0(1); x0=1.0;
  904. Vector<Real> alp0((p1+1)*(p1+2)/2);
  905. LegPoly(alp0, x0, p1);
  906. Vector<Real> alp((p1+1) * (p1+1)*(p1+2)/2);
  907. LegPoly(alp, qx1, p1);
  908. for(Long j=0;j<p1+1;j++){
  909. for(Long i=0;i<p1+1;i++){
  910. Yf[i]+=4*M_PI/(2*j+1) * alp0[j] * alp[j*(p1+1)+i];
  911. }
  912. }
  913. }
  914. Sw.ReInit(p1+1);
  915. for(Long i=0;i<p1+1;i++){
  916. Sw[i]=(qw1[i]*M_PI/p1)*Yf[i]/cos(acos(qx1[i])/2);
  917. }
  918. }
  919. return Sw;
  920. }
  921. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourier(Long p0, Long p1){
  922. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  923. Matrix<Real>& Mf =MatrixStore().Mf_ [p0*SCTL_SHMAXDEG+p1];
  924. if(!Mf.Dim(0)){
  925. const Real SQRT2PI=sqrt(2*M_PI);
  926. { // Set Mf
  927. Matrix<Real> M(2*p0,2*p1);
  928. for(Long j=0;j<2*p1;j++){
  929. M[0][j]=SQRT2PI*1.0;
  930. for(Long k=1;k<p0;k++){
  931. M[2*k-1][j]=SQRT2PI*cos(j*k*M_PI/p1);
  932. M[2*k-0][j]=SQRT2PI*sin(j*k*M_PI/p1);
  933. }
  934. M[2*p0-1][j]=SQRT2PI*cos(j*p0*M_PI/p1);
  935. }
  936. Mf=M;
  937. }
  938. }
  939. return Mf;
  940. }
  941. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourierInv(Long p0, Long p1){
  942. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  943. Matrix<Real>& Mf =MatrixStore().Mfinv_ [p0*SCTL_SHMAXDEG+p1];
  944. if(!Mf.Dim(0)){
  945. const Real INVSQRT2PI=1.0/sqrt(2*M_PI)/p0;
  946. { // Set Mf
  947. Matrix<Real> M(2*p0,2*p1);
  948. M.SetZero();
  949. if(p1>p0) p1=p0;
  950. for(Long j=0;j<2*p0;j++){
  951. M[j][0]=INVSQRT2PI*0.5;
  952. for(Long k=1;k<p1;k++){
  953. M[j][2*k-1]=INVSQRT2PI*cos(j*k*M_PI/p0);
  954. M[j][2*k-0]=INVSQRT2PI*sin(j*k*M_PI/p0);
  955. }
  956. M[j][2*p1-1]=INVSQRT2PI*cos(j*p1*M_PI/p0);
  957. }
  958. if(p1==p0) for(Long j=0;j<2*p0;j++) M[j][2*p1-1]*=0.5;
  959. Mf=M;
  960. }
  961. }
  962. return Mf;
  963. }
  964. template <class Real> const FFT<Real>& SphericalHarmonics<Real>::OpFourier(Long Np){
  965. assert(Np<SCTL_SHMAXDEG);
  966. auto& Mf =MatrixStore().Mfftinv_ [Np];
  967. #pragma omp critical (SCTL_FFT_PLAN0)
  968. if(!Mf.Dim(0)){
  969. StaticArray<Long,1> fft_dim = {Np};
  970. Mf.Setup(FFT_Type::C2R, 1, Vector<Long>(1,fft_dim,false));
  971. }
  972. return Mf;
  973. }
  974. template <class Real> const FFT<Real>& SphericalHarmonics<Real>::OpFourierInv(Long Np){
  975. assert(Np<SCTL_SHMAXDEG);
  976. auto& Mf =MatrixStore().Mfft_ [Np];
  977. #pragma omp critical (SCTL_FFT_PLAN1)
  978. if(!Mf.Dim(0)){
  979. StaticArray<Long,1> fft_dim = {Np};
  980. Mf.Setup(FFT_Type::R2C, 1, Vector<Long>(1,fft_dim,false));
  981. }
  982. return Mf;
  983. }
  984. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourierGrad(Long p0, Long p1){
  985. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  986. Matrix<Real>& Mdf=MatrixStore().Mdf_[p0*SCTL_SHMAXDEG+p1];
  987. if(!Mdf.Dim(0)){
  988. const Real SQRT2PI=sqrt(2*M_PI);
  989. { // Set Mdf_
  990. Matrix<Real> M(2*p0,2*p1);
  991. for(Long j=0;j<2*p1;j++){
  992. M[0][j]=SQRT2PI*0.0;
  993. for(Long k=1;k<p0;k++){
  994. M[2*k-1][j]=-SQRT2PI*k*sin(j*k*M_PI/p1);
  995. M[2*k-0][j]= SQRT2PI*k*cos(j*k*M_PI/p1);
  996. }
  997. M[2*p0-1][j]=-SQRT2PI*p0*sin(j*p0*M_PI/p1);
  998. }
  999. Mdf=M;
  1000. }
  1001. }
  1002. return Mdf;
  1003. }
  1004. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendre(Long p0, Long p1){
  1005. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1006. std::vector<Matrix<Real>>& Ml =MatrixStore().Ml_ [p0*SCTL_SHMAXDEG+p1];
  1007. if(!Ml.size()){
  1008. const Vector<Real>& qx1 = LegendreNodes(p1);
  1009. Vector<Real> alp(qx1.Dim()*(p0+1)*(p0+2)/2);
  1010. LegPoly(alp, qx1, p0);
  1011. Ml.resize(p0+1);
  1012. auto ptr = alp.begin();
  1013. for(Long i=0;i<=p0;i++){
  1014. Ml[i].ReInit(p0+1-i, qx1.Dim(), ptr);
  1015. ptr+=Ml[i].Dim(0)*Ml[i].Dim(1);
  1016. }
  1017. }
  1018. return Ml;
  1019. }
  1020. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendreInv(Long p0, Long p1){
  1021. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1022. std::vector<Matrix<Real>>& Ml =MatrixStore().Mlinv_ [p0*SCTL_SHMAXDEG+p1];
  1023. if(!Ml.size()){
  1024. const Vector<Real>& qx1 = LegendreNodes(p0);
  1025. const Vector<Real>& qw1 = LegendreWeights(p0);
  1026. Vector<Real> alp(qx1.Dim()*(p1+1)*(p1+2)/2);
  1027. LegPoly(alp, qx1, p1);
  1028. Ml.resize(p1+1);
  1029. auto ptr = alp.begin();
  1030. for(Long i=0;i<=p1;i++){
  1031. Ml[i].ReInit(qx1.Dim(), p1+1-i);
  1032. Matrix<Real> M(p1+1-i, qx1.Dim(), ptr, false);
  1033. for(Long j=0;j<p1+1-i;j++){ // Transpose and weights
  1034. for(Long k=0;k<qx1.Dim();k++){
  1035. Ml[i][k][j]=M[j][k]*qw1[k]*2*M_PI;
  1036. }
  1037. }
  1038. ptr+=Ml[i].Dim(0)*Ml[i].Dim(1);
  1039. }
  1040. }
  1041. return Ml;
  1042. }
  1043. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendreGrad(Long p0, Long p1){
  1044. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1045. std::vector<Matrix<Real>>& Mdl=MatrixStore().Mdl_[p0*SCTL_SHMAXDEG+p1];
  1046. if(!Mdl.size()){
  1047. const Vector<Real>& qx1 = LegendreNodes(p1);
  1048. Vector<Real> alp(qx1.Dim()*(p0+1)*(p0+2)/2);
  1049. LegPolyDeriv(alp, qx1, p0);
  1050. Mdl.resize(p0+1);
  1051. auto ptr = alp.begin();
  1052. for(Long i=0;i<=p0;i++){
  1053. Mdl[i].ReInit(p0+1-i, qx1.Dim(), ptr);
  1054. ptr+=Mdl[i].Dim(0)*Mdl[i].Dim(1);
  1055. }
  1056. }
  1057. return Mdl;
  1058. }
  1059. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatRotate(Long p0){
  1060. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1061. { // Set coeff_perm
  1062. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1063. Long itr=0;
  1064. for(Long i=0;i<2*p0;i++){
  1065. Long m=(i+1)/2;
  1066. for(Long n=m;n<=p0;n++){
  1067. coeff_perm[n][i]=itr;
  1068. itr++;
  1069. }
  1070. }
  1071. }
  1072. assert(p0<SCTL_SHMAXDEG);
  1073. std::vector<Matrix<Real>>& Mr=MatrixStore().Mr_[p0];
  1074. if(!Mr.size()){
  1075. const Real SQRT2PI=sqrt(2*M_PI);
  1076. Long Ncoef=p0*(p0+2);
  1077. Long Ngrid=2*p0*(p0+1);
  1078. Long Naleg=(p0+1)*(p0+2)/2;
  1079. Matrix<Real> Mcoord0(3,Ngrid);
  1080. const Vector<Real>& x=LegendreNodes(p0);
  1081. for(Long i=0;i<p0+1;i++){ // Set Mcoord0
  1082. for(Long j=0;j<2*p0;j++){
  1083. Mcoord0[0][i*2*p0+j]=x[i];
  1084. Mcoord0[1][i*2*p0+j]=sqrt(1-x[i]*x[i])*sin(M_PI*j/p0);
  1085. Mcoord0[2][i*2*p0+j]=sqrt(1-x[i]*x[i])*cos(M_PI*j/p0);
  1086. }
  1087. }
  1088. for(Long l=0;l<p0+1;l++){ // For each rotation angle
  1089. Matrix<Real> Mcoord1;
  1090. { // Rotate coordinates
  1091. Matrix<Real> M(COORD_DIM, COORD_DIM);
  1092. Real cos_=-x[l];
  1093. Real sin_=-sqrt(1.0-x[l]*x[l]);
  1094. M[0][0]= cos_; M[0][1]=0; M[0][2]=-sin_;
  1095. M[1][0]= 0; M[1][1]=1; M[1][2]= 0;
  1096. M[2][0]= sin_; M[2][1]=0; M[2][2]= cos_;
  1097. Mcoord1=M*Mcoord0;
  1098. }
  1099. Matrix<Real> Mleg(Naleg, Ngrid);
  1100. { // Set Mleg
  1101. const Vector<Real> Vcoord1(Mcoord1.Dim(0)*Mcoord1.Dim(1), Mcoord1.begin(), false);
  1102. Vector<Real> Vleg(Mleg.Dim(0)*Mleg.Dim(1), Mleg.begin(), false);
  1103. LegPoly(Vleg, Vcoord1, p0);
  1104. }
  1105. Vector<Real> theta(Ngrid);
  1106. for(Long i=0;i<theta.Dim();i++){ // Set theta
  1107. theta[i]=atan2(Mcoord1[1][i],Mcoord1[2][i]);
  1108. }
  1109. Matrix<Real> Mcoef2grid(Ncoef, Ngrid);
  1110. { // Build Mcoef2grid
  1111. Long offset0=0;
  1112. Long offset1=0;
  1113. for(Long i=0;i<p0+1;i++){
  1114. Long len=p0+1-i;
  1115. { // P * cos
  1116. for(Long j=0;j<len;j++){
  1117. for(Long k=0;k<Ngrid;k++){
  1118. Mcoef2grid[offset1+j][k]=SQRT2PI*Mleg[offset0+j][k]*cos(i*theta[k]);
  1119. }
  1120. }
  1121. offset1+=len;
  1122. }
  1123. if(i!=0 && i!=p0){ // P * sin
  1124. for(Long j=0;j<len;j++){
  1125. for(Long k=0;k<Ngrid;k++){
  1126. Mcoef2grid[offset1+j][k]=SQRT2PI*Mleg[offset0+j][k]*sin(i*theta[k]);
  1127. }
  1128. }
  1129. offset1+=len;
  1130. }
  1131. offset0+=len;
  1132. }
  1133. assert(offset0==Naleg);
  1134. assert(offset1==Ncoef);
  1135. }
  1136. Vector<Real> Vcoef2coef(Ncoef*Ncoef);
  1137. Vector<Real> Vcoef2grid(Ncoef*Ngrid, Mcoef2grid[0], false);
  1138. Grid2SHC(Vcoef2grid, p0+1, 2*p0, p0, Vcoef2coef, SHCArrange::COL_MAJOR_NONZERO);
  1139. Matrix<Real> Mcoef2coef(Ncoef, Ncoef, Vcoef2coef.begin(), false);
  1140. for(Long n=0;n<=p0;n++){ // Create matrices for fast rotation
  1141. Matrix<Real> M(coeff_perm[n].size(),coeff_perm[n].size());
  1142. for(Long i=0;i<(Long)coeff_perm[n].size();i++){
  1143. for(Long j=0;j<(Long)coeff_perm[n].size();j++){
  1144. M[i][j]=Mcoef2coef[coeff_perm[n][i]][coeff_perm[n][j]];
  1145. }
  1146. }
  1147. Mr.push_back(M);
  1148. }
  1149. }
  1150. }
  1151. return Mr;
  1152. }
  1153. template <class Real> void SphericalHarmonics<Real>::SHC2GridTranspose(const Vector<Real>& X, Long p0, Long p1, Vector<Real>& S){
  1154. Matrix<Real> Mf =SphericalHarmonics<Real>::MatFourier(p1,p0).Transpose();
  1155. std::vector<Matrix<Real>> Ml =SphericalHarmonics<Real>::MatLegendre(p1,p0);
  1156. for(Long i=0;i<(Long)Ml.size();i++) Ml[i]=Ml[i].Transpose();
  1157. assert(p1==(Long)Ml.size()-1);
  1158. assert(p0==Mf.Dim(0)/2);
  1159. assert(p1==Mf.Dim(1)/2);
  1160. Long N=X.Dim()/(2*p0*(p0+1));
  1161. assert(N*2*p0*(p0+1)==X.Dim());
  1162. if(S.Dim()!=N*(p1*(p1+2))) S.ReInit(N*(p1*(p1+2)));
  1163. Vector<Real> B0, B1;
  1164. B0.ReInit(N* p1*(p1+2));
  1165. B1.ReInit(N*2*p1*(p0+1));
  1166. #pragma omp parallel
  1167. { // Evaluate Fourier and transpose
  1168. Integer tid=omp_get_thread_num();
  1169. Integer omp_p=omp_get_num_threads();
  1170. Long a=(tid+0)*N*(p0+1)/omp_p;
  1171. Long b=(tid+1)*N*(p0+1)/omp_p;
  1172. const Long block_size=16;
  1173. Matrix<Real> B2(block_size,2*p1);
  1174. for(Long i0=a;i0<b;i0+=block_size){
  1175. Long i1=std::min(b,i0+block_size);
  1176. const Matrix<Real> Min (i1-i0,2*p0, (Iterator<Real>)X.begin()+i0*2*p0, false);
  1177. Matrix<Real> Mout(i1-i0,2*p1, B2.begin(), false);
  1178. Matrix<Real>::GEMM(Mout, Min, Mf);
  1179. for(Long i=i0;i<i1;i++){
  1180. for(Long j=0;j<2*p1;j++){
  1181. B1[j*N*(p0+1)+i]=B2[i-i0][j];
  1182. }
  1183. }
  1184. }
  1185. }
  1186. #pragma omp parallel
  1187. { // Evaluate Legendre polynomial
  1188. Integer tid=omp_get_thread_num();
  1189. Integer omp_p=omp_get_num_threads();
  1190. Long offset0=0;
  1191. Long offset1=0;
  1192. for(Long i=0;i<p1+1;i++){
  1193. Long N0=2*N;
  1194. if(i==0 || i==p1) N0=N;
  1195. Matrix<Real> Min (N0, p0+1 , B1.begin()+offset0, false);
  1196. Matrix<Real> Mout(N0, p1+1-i, B0.begin()+offset1, false);
  1197. { // Mout = Min * Ml[i] // split between threads
  1198. Long a=(tid+0)*N0/omp_p;
  1199. Long b=(tid+1)*N0/omp_p;
  1200. if(a<b){
  1201. Matrix<Real> Min_ (b-a, Min .Dim(1), Min [a], false);
  1202. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  1203. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  1204. }
  1205. }
  1206. offset0+=Min .Dim(0)*Min .Dim(1);
  1207. offset1+=Mout.Dim(0)*Mout.Dim(1);
  1208. }
  1209. }
  1210. #pragma omp parallel
  1211. { // S <-- Rearrange(B0)
  1212. Integer tid=omp_get_thread_num();
  1213. Integer omp_p=omp_get_num_threads();
  1214. Long a=(tid+0)*N/omp_p;
  1215. Long b=(tid+1)*N/omp_p;
  1216. for(Long i=a;i<b;i++){
  1217. Long offset=0;
  1218. for(Long j=0;j<2*p1;j++){
  1219. Long len=p1+1-(j+1)/2;
  1220. Real* B_=&B0[i*len+N*offset];
  1221. Real* S_=&S[i*p1*(p1+2)+offset];
  1222. for(Long k=0;k<len;k++) S_[k]=B_[k];
  1223. offset+=len;
  1224. }
  1225. }
  1226. }
  1227. }
  1228. template <class Real> void SphericalHarmonics<Real>::RotateAll(const Vector<Real>& S, Long p0, Long dof, Vector<Real>& S_){
  1229. const std::vector<Matrix<Real>>& Mr=MatRotate(p0);
  1230. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1231. { // Set coeff_perm
  1232. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1233. Long itr=0;
  1234. for(Long i=0;i<2*p0;i++){
  1235. Long m=(i+1)/2;
  1236. for(Long n=m;n<=p0;n++){
  1237. coeff_perm[n][i]=itr;
  1238. itr++;
  1239. }
  1240. }
  1241. }
  1242. Long Ncoef=p0*(p0+2);
  1243. Long N=S.Dim()/Ncoef/dof;
  1244. assert(N*Ncoef*dof==S.Dim());
  1245. if(S_.Dim()!=N*dof*Ncoef*p0*(p0+1)) S_.ReInit(N*dof*Ncoef*p0*(p0+1));
  1246. const Matrix<Real> S0(N*dof, Ncoef, (Iterator<Real>)S.begin(), false);
  1247. Matrix<Real> S1(N*dof*p0*(p0+1), Ncoef, S_.begin(), false);
  1248. #pragma omp parallel
  1249. { // Construct all p0*(p0+1) rotations
  1250. Integer tid=omp_get_thread_num();
  1251. Integer omp_p=omp_get_num_threads();
  1252. Matrix<Real> B0(dof*p0,Ncoef); // memory buffer
  1253. std::vector<Matrix<Real>> Bi(p0+1), Bo(p0+1); // memory buffers
  1254. for(Long i=0;i<=p0;i++){ // initialize Bi, Bo
  1255. Bi[i].ReInit(dof*p0,coeff_perm[i].size());
  1256. Bo[i].ReInit(dof*p0,coeff_perm[i].size());
  1257. }
  1258. Long a=(tid+0)*N/omp_p;
  1259. Long b=(tid+1)*N/omp_p;
  1260. for(Long i=a;i<b;i++){
  1261. for(Long d=0;d<dof;d++){
  1262. for(Long j=0;j<p0;j++){
  1263. Long offset=0;
  1264. for(Long k=0;k<p0+1;k++){
  1265. Real r[2]={cos(k*j*M_PI/p0),-sin(k*j*M_PI/p0)}; // exp(i*k*theta)
  1266. Long len=p0+1-k;
  1267. if(k!=0 && k!=p0){
  1268. for(Long l=0;l<len;l++){
  1269. Real x[2];
  1270. x[0]=S0[i*dof+d][offset+len*0+l];
  1271. x[1]=S0[i*dof+d][offset+len*1+l];
  1272. B0[j*dof+d][offset+len*0+l]=x[0]*r[0]-x[1]*r[1];
  1273. B0[j*dof+d][offset+len*1+l]=x[0]*r[1]+x[1]*r[0];
  1274. }
  1275. offset+=2*len;
  1276. }else{
  1277. for(Long l=0;l<len;l++){
  1278. B0[j*dof+d][offset+l]=S0[i*dof+d][offset+l];
  1279. }
  1280. offset+=len;
  1281. }
  1282. }
  1283. assert(offset==Ncoef);
  1284. }
  1285. }
  1286. { // Fast rotation
  1287. for(Long k=0;k<dof*p0;k++){ // forward permutation
  1288. for(Long l=0;l<=p0;l++){
  1289. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1290. Bi[l][k][j]=B0[k][coeff_perm[l][j]];
  1291. }
  1292. }
  1293. }
  1294. for(Long t=0;t<=p0;t++){
  1295. for(Long l=0;l<=p0;l++){ // mat-vec
  1296. Matrix<Real>::GEMM(Bo[l],Bi[l],Mr[t*(p0+1)+l]);
  1297. }
  1298. Matrix<Real> Mout(dof*p0,Ncoef, S1[(i*(p0+1)+t)*dof*p0], false);
  1299. for(Long k=0;k<dof*p0;k++){ // reverse permutation
  1300. for(Long l=0;l<=p0;l++){
  1301. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1302. Mout[k][coeff_perm[l][j]]=Bo[l][k][j];
  1303. }
  1304. }
  1305. }
  1306. }
  1307. }
  1308. }
  1309. }
  1310. }
  1311. template <class Real> void SphericalHarmonics<Real>::RotateTranspose(const Vector<Real>& S_, Long p0, Long dof, Vector<Real>& S){
  1312. std::vector<Matrix<Real>> Mr=MatRotate(p0);
  1313. for(Long i=0;i<(Long)Mr.size();i++) Mr[i]=Mr[i].Transpose();
  1314. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1315. { // Set coeff_perm
  1316. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1317. Long itr=0;
  1318. for(Long i=0;i<2*p0;i++){
  1319. Long m=(i+1)/2;
  1320. for(Long n=m;n<=p0;n++){
  1321. coeff_perm[n][i]=itr;
  1322. itr++;
  1323. }
  1324. }
  1325. }
  1326. Long Ncoef=p0*(p0+2);
  1327. Long N=S_.Dim()/Ncoef/dof/(p0*(p0+1));
  1328. assert(N*Ncoef*dof*(p0*(p0+1))==S_.Dim());
  1329. if(S.Dim()!=N*dof*Ncoef*p0*(p0+1)) S.ReInit(N*dof*Ncoef*p0*(p0+1));
  1330. Matrix<Real> S0(N*dof*p0*(p0+1), Ncoef, S.begin(), false);
  1331. const Matrix<Real> S1(N*dof*p0*(p0+1), Ncoef, (Iterator<Real>)S_.begin(), false);
  1332. #pragma omp parallel
  1333. { // Transpose all p0*(p0+1) rotations
  1334. Integer tid=omp_get_thread_num();
  1335. Integer omp_p=omp_get_num_threads();
  1336. Matrix<Real> B0(dof*p0,Ncoef); // memory buffer
  1337. std::vector<Matrix<Real>> Bi(p0+1), Bo(p0+1); // memory buffers
  1338. for(Long i=0;i<=p0;i++){ // initialize Bi, Bo
  1339. Bi[i].ReInit(dof*p0,coeff_perm[i].size());
  1340. Bo[i].ReInit(dof*p0,coeff_perm[i].size());
  1341. }
  1342. Long a=(tid+0)*N/omp_p;
  1343. Long b=(tid+1)*N/omp_p;
  1344. for(Long i=a;i<b;i++){
  1345. for(Long t=0;t<p0+1;t++){
  1346. Long idx0=(i*(p0+1)+t)*p0*dof;
  1347. { // Fast rotation
  1348. const Matrix<Real> Min(p0*dof, Ncoef, (Iterator<Real>)S1[idx0], false);
  1349. for(Long k=0;k<dof*p0;k++){ // forward permutation
  1350. for(Long l=0;l<=p0;l++){
  1351. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1352. Bi[l][k][j]=Min[k][coeff_perm[l][j]];
  1353. }
  1354. }
  1355. }
  1356. for(Long l=0;l<=p0;l++){ // mat-vec
  1357. Matrix<Real>::GEMM(Bo[l],Bi[l],Mr[t*(p0+1)+l]);
  1358. }
  1359. for(Long k=0;k<dof*p0;k++){ // reverse permutation
  1360. for(Long l=0;l<=p0;l++){
  1361. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1362. B0[k][coeff_perm[l][j]]=Bo[l][k][j];
  1363. }
  1364. }
  1365. }
  1366. }
  1367. for(Long j=0;j<p0;j++){
  1368. for(Long d=0;d<dof;d++){
  1369. Long idx1=idx0+j*dof+d;
  1370. Long offset=0;
  1371. for(Long k=0;k<p0+1;k++){
  1372. Real r[2]={cos(k*j*M_PI/p0),sin(k*j*M_PI/p0)}; // exp(i*k*theta)
  1373. Long len=p0+1-k;
  1374. if(k!=0 && k!=p0){
  1375. for(Long l=0;l<len;l++){
  1376. Real x[2];
  1377. x[0]=B0[j*dof+d][offset+len*0+l];
  1378. x[1]=B0[j*dof+d][offset+len*1+l];
  1379. S0[idx1][offset+len*0+l]=x[0]*r[0]-x[1]*r[1];
  1380. S0[idx1][offset+len*1+l]=x[0]*r[1]+x[1]*r[0];
  1381. }
  1382. offset+=2*len;
  1383. }else{
  1384. for(Long l=0;l<len;l++){
  1385. S0[idx1][offset+l]=B0[j*dof+d][offset+l];
  1386. }
  1387. offset+=len;
  1388. }
  1389. }
  1390. assert(offset==Ncoef);
  1391. }
  1392. }
  1393. }
  1394. }
  1395. }
  1396. }
  1397. template <class Real> void SphericalHarmonics<Real>::StokesSingularInteg(const Vector<Real>& S, Long p0, Long p1, Vector<Real>* SLMatrix, Vector<Real>* DLMatrix){
  1398. Long Ngrid=2*p0*(p0+1);
  1399. Long Ncoef= p0*(p0+2);
  1400. Long Nves=S.Dim()/(Ngrid*COORD_DIM);
  1401. if(SLMatrix) SLMatrix->ReInit(Nves*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM));
  1402. if(DLMatrix) DLMatrix->ReInit(Nves*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM));
  1403. Long BLOCK_SIZE=(Long)6e9/((3*2*p1*(p1+1))*(3*2*p0*(p0+1))*2*8); // Limit memory usage to 6GB
  1404. BLOCK_SIZE=std::min<Long>(BLOCK_SIZE,omp_get_max_threads());
  1405. BLOCK_SIZE=std::max<Long>(BLOCK_SIZE,1);
  1406. for(Long a=0;a<Nves;a+=BLOCK_SIZE){
  1407. Long b=std::min(a+BLOCK_SIZE, Nves);
  1408. Vector<Real> _SLMatrix, _DLMatrix;
  1409. if(SLMatrix) _SLMatrix.ReInit((b-a)*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), SLMatrix->begin()+a*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), false);
  1410. if(DLMatrix) _DLMatrix.ReInit((b-a)*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), DLMatrix->begin()+a*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), false);
  1411. const Vector<Real> _S ((b-a)*(Ngrid*COORD_DIM) , (Iterator<Real>)S.begin()+a*(Ngrid*COORD_DIM), false);
  1412. if(SLMatrix && DLMatrix) StokesSingularInteg_< true, true>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1413. else if(SLMatrix) StokesSingularInteg_< true, false>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1414. else if(DLMatrix) StokesSingularInteg_<false, true>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1415. }
  1416. }
  1417. template <class Real> template <bool SLayer, bool DLayer> void SphericalHarmonics<Real>::StokesSingularInteg_(const Vector<Real>& X0, Long p0, Long p1, Vector<Real>& SL, Vector<Real>& DL){
  1418. Profile::Tic("Rotate");
  1419. Vector<Real> S0, S;
  1420. SphericalHarmonics<Real>::Grid2SHC(X0, p0+1, 2*p0, p0, S0, SHCArrange::COL_MAJOR_NONZERO);
  1421. SphericalHarmonics<Real>::RotateAll(S0, p0, COORD_DIM, S);
  1422. Profile::Toc();
  1423. Profile::Tic("Upsample");
  1424. Vector<Real> X, X_theta, X_phi, trg;
  1425. SphericalHarmonics<Real>::SHC2Grid(S, SHCArrange::COL_MAJOR_NONZERO, p0, p1+1, 2*p1, &X, &X_phi, &X_theta);
  1426. SphericalHarmonics<Real>::SHC2Pole(S, SHCArrange::COL_MAJOR_NONZERO, p0, trg);
  1427. Profile::Toc();
  1428. Profile::Tic("Stokes");
  1429. Vector<Real> SL0, DL0;
  1430. { // Stokes kernel
  1431. //Long M0=2*p0*(p0+1);
  1432. Long M1=2*p1*(p1+1);
  1433. Long N=trg.Dim()/(2*COORD_DIM);
  1434. assert(X.Dim()==M1*COORD_DIM*N);
  1435. if(SLayer && SL0.Dim()!=N*2*6*M1) SL0.ReInit(2*N*6*M1);
  1436. if(DLayer && DL0.Dim()!=N*2*6*M1) DL0.ReInit(2*N*6*M1);
  1437. const Vector<Real>& qw=SphericalHarmonics<Real>::SingularWeights(p1);
  1438. const Real scal_const_dl = 3.0/(4.0*M_PI);
  1439. const Real scal_const_sl = 1.0/(8.0*M_PI);
  1440. static Real eps=-1;
  1441. if(eps<0){
  1442. eps=1;
  1443. while(eps*(Real)0.5+(Real)1.0>1.0) eps*=0.5;
  1444. }
  1445. #pragma omp parallel
  1446. {
  1447. Integer tid=omp_get_thread_num();
  1448. Integer omp_p=omp_get_num_threads();
  1449. Long a=(tid+0)*N/omp_p;
  1450. Long b=(tid+1)*N/omp_p;
  1451. for(Long i=a;i<b;i++){
  1452. for(Long t=0;t<2;t++){
  1453. Real tx, ty, tz;
  1454. { // Read target coordinates
  1455. tx=trg[i*2*COORD_DIM+0*2+t];
  1456. ty=trg[i*2*COORD_DIM+1*2+t];
  1457. tz=trg[i*2*COORD_DIM+2*2+t];
  1458. }
  1459. for(Long j0=0;j0<p1+1;j0++){
  1460. for(Long j1=0;j1<2*p1;j1++){
  1461. Long s=2*p1*j0+j1;
  1462. Real dx, dy, dz;
  1463. { // Compute dx, dy, dz
  1464. dx=tx-X[(i*COORD_DIM+0)*M1+s];
  1465. dy=ty-X[(i*COORD_DIM+1)*M1+s];
  1466. dz=tz-X[(i*COORD_DIM+2)*M1+s];
  1467. }
  1468. Real nx, ny, nz;
  1469. { // Compute source normal
  1470. Real x_theta=X_phi[(i*COORD_DIM+0)*M1+s];
  1471. Real y_theta=X_phi[(i*COORD_DIM+1)*M1+s];
  1472. Real z_theta=X_phi[(i*COORD_DIM+2)*M1+s];
  1473. Real x_phi=X_theta[(i*COORD_DIM+0)*M1+s];
  1474. Real y_phi=X_theta[(i*COORD_DIM+1)*M1+s];
  1475. Real z_phi=X_theta[(i*COORD_DIM+2)*M1+s];
  1476. nx=(y_theta*z_phi-z_theta*y_phi);
  1477. ny=(z_theta*x_phi-x_theta*z_phi);
  1478. nz=(x_theta*y_phi-y_theta*x_phi);
  1479. }
  1480. Real area_elem=1.0;
  1481. if(SLayer){ // Compute area_elem
  1482. area_elem=sqrt(nx*nx+ny*ny+nz*nz);
  1483. }
  1484. Real rinv, rinv2;
  1485. { // Compute rinv, rinv2
  1486. Real r2=dx*dx+dy*dy+dz*dz;
  1487. rinv=1.0/sqrt(r2);
  1488. if(r2<=eps) rinv=0;
  1489. rinv2=rinv*rinv;
  1490. }
  1491. if(DLayer){
  1492. Real rinv5=rinv2*rinv2*rinv;
  1493. Real r_dot_n_rinv5=scal_const_dl*qw[j0*t+(p1-j0)*(1-t)] * (nx*dx+ny*dy+nz*dz)*rinv5;
  1494. DL0[((i*2+t)*6+0)*M1+s]=dx*dx*r_dot_n_rinv5;
  1495. DL0[((i*2+t)*6+1)*M1+s]=dx*dy*r_dot_n_rinv5;
  1496. DL0[((i*2+t)*6+2)*M1+s]=dx*dz*r_dot_n_rinv5;
  1497. DL0[((i*2+t)*6+3)*M1+s]=dy*dy*r_dot_n_rinv5;
  1498. DL0[((i*2+t)*6+4)*M1+s]=dy*dz*r_dot_n_rinv5;
  1499. DL0[((i*2+t)*6+5)*M1+s]=dz*dz*r_dot_n_rinv5;
  1500. }
  1501. if(SLayer){
  1502. Real area_rinv =scal_const_sl*qw[j0*t+(p1-j0)*(1-t)] * area_elem*rinv;
  1503. Real area_rinv2=area_rinv*rinv2;
  1504. SL0[((i*2+t)*6+0)*M1+s]=area_rinv+dx*dx*area_rinv2;
  1505. SL0[((i*2+t)*6+1)*M1+s]= dx*dy*area_rinv2;
  1506. SL0[((i*2+t)*6+2)*M1+s]= dx*dz*area_rinv2;
  1507. SL0[((i*2+t)*6+3)*M1+s]=area_rinv+dy*dy*area_rinv2;
  1508. SL0[((i*2+t)*6+4)*M1+s]= dy*dz*area_rinv2;
  1509. SL0[((i*2+t)*6+5)*M1+s]=area_rinv+dz*dz*area_rinv2;
  1510. }
  1511. }
  1512. }
  1513. }
  1514. }
  1515. }
  1516. Profile::Add_FLOP(20*(2*p1)*(p1+1)*2*N);
  1517. if(SLayer) Profile::Add_FLOP((19+6)*(2*p1)*(p1+1)*2*N);
  1518. if(DLayer) Profile::Add_FLOP( 22 *(2*p1)*(p1+1)*2*N);
  1519. }
  1520. Profile::Toc();
  1521. Profile::Tic("UpsampleTranspose");
  1522. Vector<Real> SL1, DL1;
  1523. SphericalHarmonics<Real>::SHC2GridTranspose(SL0, p1, p0, SL1);
  1524. SphericalHarmonics<Real>::SHC2GridTranspose(DL0, p1, p0, DL1);
  1525. Profile::Toc();
  1526. Profile::Tic("RotateTranspose");
  1527. Vector<Real> SL2, DL2;
  1528. SphericalHarmonics<Real>::RotateTranspose(SL1, p0, 2*6, SL2);
  1529. SphericalHarmonics<Real>::RotateTranspose(DL1, p0, 2*6, DL2);
  1530. Profile::Toc();
  1531. Profile::Tic("Rearrange");
  1532. Vector<Real> SL3, DL3;
  1533. { // Transpose
  1534. Long Ncoef=p0*(p0+2);
  1535. Long Ngrid=2*p0*(p0+1);
  1536. { // Transpose SL2
  1537. Long N=SL2.Dim()/(6*Ncoef*Ngrid);
  1538. SL3.ReInit(N*COORD_DIM*Ncoef*COORD_DIM*Ngrid);
  1539. #pragma omp parallel
  1540. {
  1541. Integer tid=omp_get_thread_num();
  1542. Integer omp_p=omp_get_num_threads();
  1543. Matrix<Real> B(COORD_DIM*Ncoef,Ngrid*COORD_DIM);
  1544. Long a=(tid+0)*N/omp_p;
  1545. Long b=(tid+1)*N/omp_p;
  1546. for(Long i=a;i<b;i++){
  1547. Matrix<Real> M0(Ngrid*6, Ncoef, SL2.begin()+i*Ngrid*6*Ncoef, false);
  1548. for(Long k=0;k<Ncoef;k++){ // Transpose
  1549. for(Long j=0;j<Ngrid;j++){ // TODO: needs blocking
  1550. B[k+Ncoef*0][j*COORD_DIM+0]=M0[j*6+0][k];
  1551. B[k+Ncoef*1][j*COORD_DIM+0]=M0[j*6+1][k];
  1552. B[k+Ncoef*2][j*COORD_DIM+0]=M0[j*6+2][k];
  1553. B[k+Ncoef*0][j*COORD_DIM+1]=M0[j*6+1][k];
  1554. B[k+Ncoef*1][j*COORD_DIM+1]=M0[j*6+3][k];
  1555. B[k+Ncoef*2][j*COORD_DIM+1]=M0[j*6+4][k];
  1556. B[k+Ncoef*0][j*COORD_DIM+2]=M0[j*6+2][k];
  1557. B[k+Ncoef*1][j*COORD_DIM+2]=M0[j*6+4][k];
  1558. B[k+Ncoef*2][j*COORD_DIM+2]=M0[j*6+5][k];
  1559. }
  1560. }
  1561. Matrix<Real> M1(Ncoef*COORD_DIM, COORD_DIM*Ngrid, SL3.begin()+i*COORD_DIM*Ncoef*COORD_DIM*Ngrid, false);
  1562. for(Long k=0;k<B.Dim(0);k++){ // Rearrange
  1563. for(Long j0=0;j0<COORD_DIM;j0++){
  1564. for(Long j1=0;j1<p0+1;j1++){
  1565. for(Long j2=0;j2<p0;j2++) M1[k][((j0*(p0+1)+ j1)*2+0)*p0+j2]=B[k][((j1*p0+j2)*2+0)*COORD_DIM+j0];
  1566. for(Long j2=0;j2<p0;j2++) M1[k][((j0*(p0+1)+p0-j1)*2+1)*p0+j2]=B[k][((j1*p0+j2)*2+1)*COORD_DIM+j0];
  1567. }
  1568. }
  1569. }
  1570. }
  1571. }
  1572. }
  1573. { // Transpose DL2
  1574. Long N=DL2.Dim()/(6*Ncoef*Ngrid);
  1575. DL3.ReInit(N*COORD_DIM*Ncoef*COORD_DIM*Ngrid);
  1576. #pragma omp parallel
  1577. {
  1578. Integer tid=omp_get_thread_num();
  1579. Integer omp_p=omp_get_num_threads();
  1580. Matrix<Real> B(COORD_DIM*Ncoef,Ngrid*COORD_DIM);
  1581. Long a=(tid+0)*N/omp_p;
  1582. Long b=(tid+1)*N/omp_p;
  1583. for(Long i=a;i<b;i++){
  1584. Matrix<Real> M0(Ngrid*6, Ncoef, DL2.begin()+i*Ngrid*6*Ncoef, false);
  1585. for(Long k=0;k<Ncoef;k++){ // Transpose
  1586. for(Long j=0;j<Ngrid;j++){ // TODO: needs blocking
  1587. B[k+Ncoef*0][j*COORD_DIM+0]=M0[j*6+0][k];
  1588. B[k+Ncoef*1][j*COORD_DIM+0]=M0[j*6+1][k];
  1589. B[k+Ncoef*2][j*COORD_DIM+0]=M0[j*6+2][k];
  1590. B[k+Ncoef*0][j*COORD_DIM+1]=M0[j*6+1][k];
  1591. B[k+Ncoef*1][j*COORD_DIM+1]=M0[j*6+3][k];
  1592. B[k+Ncoef*2][j*COORD_DIM+1]=M0[j*6+4][k];
  1593. B[k+Ncoef*0][j*COORD_DIM+2]=M0[j*6+2][k];
  1594. B[k+Ncoef*1][j*COORD_DIM+2]=M0[j*6+4][k];
  1595. B[k+Ncoef*2][j*COORD_DIM+2]=M0[j*6+5][k];
  1596. }
  1597. }
  1598. Matrix<Real> M1(Ncoef*COORD_DIM, COORD_DIM*Ngrid, DL3.begin()+i*COORD_DIM*Ncoef*COORD_DIM*Ngrid, false);
  1599. for(Long k=0;k<B.Dim(0);k++){ // Rearrange
  1600. for(Long j0=0;j0<COORD_DIM;j0++){
  1601. for(Long j1=0;j1<p0+1;j1++){
  1602. for(Long j2=0;j2<p0;j2++) M1[k][((j0*(p0+1)+ j1)*2+0)*p0+j2]=B[k][((j1*p0+j2)*2+0)*COORD_DIM+j0];
  1603. for(Long j2=0;j2<p0;j2++) M1[k][((j0*(p0+1)+p0-j1)*2+1)*p0+j2]=B[k][((j1*p0+j2)*2+1)*COORD_DIM+j0];
  1604. }
  1605. }
  1606. }
  1607. }
  1608. }
  1609. }
  1610. }
  1611. Profile::Toc();
  1612. Profile::Tic("Grid2SHC");
  1613. SphericalHarmonics<Real>::Grid2SHC(SL3, p0+1, 2*p0, p0, SL, SHCArrange::COL_MAJOR_NONZERO);
  1614. SphericalHarmonics<Real>::Grid2SHC(DL3, p0+1, 2*p0, p0, DL, SHCArrange::COL_MAJOR_NONZERO);
  1615. Profile::Toc();
  1616. }
  1617. } // end namespace