sph_harm.txx 82 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379
  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. SHCArrange0(B1, p1, S, arrange);
  10. }
  11. 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_theta, Vector<Real>* X_phi){
  12. Vector<Real> B0;
  13. SHCArrange1(S, arrange, p0, B0);
  14. SHC2Grid_(B0, p0, Nt, Np, X, X_phi, X_theta);
  15. }
  16. template <class Real> void SphericalHarmonics<Real>::SHCEval(const Vector<Real>& S, SHCArrange arrange, Long p0, const Vector<Real>& cos_theta_phi, Vector<Real>& X) {
  17. Long M = (p0+1) * (p0+1);
  18. Long dof;
  19. Matrix<Real> B1;
  20. { // Set B1, dof
  21. Vector<Real> B0;
  22. SHCArrange1(S, arrange, p0, B0);
  23. dof = B0.Dim() / M;
  24. assert(B0.Dim() == dof * M);
  25. B1.ReInit(dof, M);
  26. Vector<Real> B1_(B1.Dim(0) * B1.Dim(1), B1.begin(), false);
  27. SHCArrange0(B0, p0, B1_, SHCArrange::COL_MAJOR_NONZERO);
  28. }
  29. assert(B1.Dim(0) == dof);
  30. assert(B1.Dim(1) == M);
  31. Matrix<Real> SHBasis;
  32. SHBasisEval(p0, cos_theta_phi, SHBasis);
  33. assert(SHBasis.Dim(1) == M);
  34. Long N = SHBasis.Dim(0);
  35. { // Set X
  36. if (X.Dim() != N*dof) X.ReInit(N * dof);
  37. for (Long k0 = 0; k0 < N; k0++) {
  38. for (Long k1 = 0; k1 < dof; k1++) {
  39. Real X_ = 0;
  40. for (Long i = 0; i < M; i++) X_ += B1[k1][i] * SHBasis[k0][i];
  41. X[k0 * dof + k1] = X_;
  42. }
  43. }
  44. }
  45. }
  46. template <class Real> void SphericalHarmonics<Real>::SHC2Pole(const Vector<Real>& S, SHCArrange arrange, Long p0, Vector<Real>& P){
  47. Vector<Real> QP[2];
  48. { // Set QP // TODO: store these weights
  49. Vector<Real> x(1), alp;
  50. const Real SQRT2PI = sqrt<Real>(4 * const_pi<Real>());
  51. for (Long i = 0; i < 2; i++) {
  52. x = (i ? -1 : 1);
  53. LegPoly(alp, x, p0);
  54. QP[i].ReInit(p0 + 1, alp.begin());
  55. QP[i] *= SQRT2PI;
  56. }
  57. }
  58. Long M, N;
  59. { // Set M, N
  60. M = 0;
  61. if (arrange == SHCArrange::ALL) M = 2*(p0+1)*(p0+1);
  62. if (arrange == SHCArrange::ROW_MAJOR) M = (p0+1)*(p0+2);
  63. if (arrange == SHCArrange::COL_MAJOR_NONZERO) M = (p0+1)*(p0+1);
  64. if (M == 0) return;
  65. N = S.Dim() / M;
  66. assert(S.Dim() == N * M);
  67. }
  68. if(P.Dim() != N * 2) P.ReInit(N * 2);
  69. if (arrange == SHCArrange::ALL) {
  70. #pragma omp parallel
  71. { // Compute pole
  72. Integer tid = omp_get_thread_num();
  73. Integer omp_p = omp_get_num_threads();
  74. Long a = (tid + 0) * N / omp_p;
  75. Long b = (tid + 1) * N / omp_p;
  76. for (Long i = a; i < b; i++) {
  77. Real P_[2] = {0, 0};
  78. for (Long j = 0; j < p0 + 1; j++) {
  79. P_[0] += S[i*M + j*(p0+1)*2] * QP[0][j];
  80. P_[1] += S[i*M + j*(p0+1)*2] * QP[1][j];
  81. }
  82. P[2*i+0] = P_[0];
  83. P[2*i+1] = P_[1];
  84. }
  85. }
  86. }
  87. if (arrange == SHCArrange::ROW_MAJOR) {
  88. #pragma omp parallel
  89. { // Compute pole
  90. Integer tid = omp_get_thread_num();
  91. Integer omp_p = omp_get_num_threads();
  92. Long a = (tid + 0) * N / omp_p;
  93. Long b = (tid + 1) * N / omp_p;
  94. for (Long i = a; i < b; i++) {
  95. Long idx = 0;
  96. Real P_[2] = {0, 0};
  97. for (Long j = 0; j < p0 + 1; j++) {
  98. P_[0] += S[i*M+idx] * QP[0][j];
  99. P_[1] += S[i*M+idx] * QP[1][j];
  100. idx += 2*(j+1);
  101. }
  102. P[2*i+0] = P_[0];
  103. P[2*i+1] = P_[1];
  104. }
  105. }
  106. }
  107. if (arrange == SHCArrange::COL_MAJOR_NONZERO) {
  108. #pragma omp parallel
  109. { // Compute pole
  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. Real P_[2] = {0, 0};
  116. for (Long j = 0; j < p0 + 1; j++) {
  117. P_[0] += S[i*M+j] * QP[0][j];
  118. P_[1] += S[i*M+j] * QP[1][j];
  119. }
  120. P[2*i+0] = P_[0];
  121. P[2*i+1] = P_[1];
  122. }
  123. }
  124. }
  125. }
  126. 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){
  127. typedef double VTKReal;
  128. Vector<Real> SS;
  129. if (S == nullptr) {
  130. Integer p = 2;
  131. Integer Ncoeff = (p + 1) * (p + 1);
  132. Vector<Real> SSS(COORD_DIM * Ncoeff), SSS_grid;
  133. SSS.SetZero();
  134. SSS[1+0*p+0*Ncoeff] = sqrt<Real>(2.0)/sqrt<Real>(3.0);
  135. SSS[1+1*p+1*Ncoeff] = 1/sqrt<Real>(3.0);
  136. SSS[1+2*p+2*Ncoeff] = 1/sqrt<Real>(3.0);
  137. SphericalHarmonics<Real>::SHC2Grid(SSS, SHCArrange::COL_MAJOR_NONZERO, p, p+1, 2*p+2, &SSS_grid);
  138. SphericalHarmonics<Real>::Grid2SHC(SSS_grid, p+1, 2*p+2, p0, SS, arrange);
  139. S = &SS;
  140. }
  141. Vector<Real> X, Xp, V, Vp;
  142. { // Upsample X
  143. const Vector<Real>& X0=*S;
  144. SphericalHarmonics<Real>::SHC2Grid(X0, arrange, p0, p1+1, 2*p1, &X);
  145. SphericalHarmonics<Real>::SHC2Pole(X0, arrange, p0, Xp);
  146. }
  147. if(v_ptr){ // Upsample V
  148. const Vector<Real>& X0=*v_ptr;
  149. SphericalHarmonics<Real>::SHC2Grid(X0, arrange, p0, p1+1, 2*p1, &V);
  150. SphericalHarmonics<Real>::SHC2Pole(X0, arrange, p0, Vp);
  151. }
  152. std::vector<VTKReal> point_coord;
  153. std::vector<VTKReal> point_value;
  154. std::vector<int32_t> poly_connect;
  155. std::vector<int32_t> poly_offset;
  156. { // Set point_coord, point_value, poly_connect
  157. Long N_ves = X.Dim()/(2*p1*(p1+1)*COORD_DIM); // Number of vesicles
  158. assert(Xp.Dim() == N_ves*2*COORD_DIM);
  159. for(Long k=0;k<N_ves;k++){ // Set point_coord
  160. Real C[COORD_DIM]={0,0,0};
  161. if(period>0){
  162. for(Integer l=0;l<COORD_DIM;l++) C[l]=0;
  163. for(Long i=0;i<p1+1;i++){
  164. for(Long j=0;j<2*p1;j++){
  165. for(Integer l=0;l<COORD_DIM;l++){
  166. C[l]+=X[j+2*p1*(i+(p1+1)*(l+k*COORD_DIM))];
  167. }
  168. }
  169. }
  170. for(Integer l=0;l<COORD_DIM;l++) C[l]+=Xp[0+2*(l+k*COORD_DIM)];
  171. for(Integer l=0;l<COORD_DIM;l++) C[l]+=Xp[1+2*(l+k*COORD_DIM)];
  172. for(Integer l=0;l<COORD_DIM;l++) C[l]/=2*p1*(p1+1)+2;
  173. for(Integer l=0;l<COORD_DIM;l++) C[l]=(round(C[l]/period))*period;
  174. }
  175. for(Long i=0;i<p1+1;i++){
  176. for(Long j=0;j<2*p1;j++){
  177. for(Integer l=0;l<COORD_DIM;l++){
  178. point_coord.push_back(X[j+2*p1*(i+(p1+1)*(l+k*COORD_DIM))]-C[l]);
  179. }
  180. }
  181. }
  182. for(Integer l=0;l<COORD_DIM;l++) point_coord.push_back(Xp[0+2*(l+k*COORD_DIM)]-C[l]);
  183. for(Integer l=0;l<COORD_DIM;l++) point_coord.push_back(Xp[1+2*(l+k*COORD_DIM)]-C[l]);
  184. }
  185. if(v_ptr) {
  186. Long data__dof = V.Dim() / (2*p1*(p1+1));
  187. for(Long k=0;k<N_ves;k++){ // Set point_value
  188. for(Long i=0;i<p1+1;i++){
  189. for(Long j=0;j<2*p1;j++){
  190. for(Long l=0;l<data__dof;l++){
  191. point_value.push_back(V[j+2*p1*(i+(p1+1)*(l+k*data__dof))]);
  192. }
  193. }
  194. }
  195. for(Long l=0;l<data__dof;l++) point_value.push_back(Vp[0+2*(l+k*data__dof)]);
  196. for(Long l=0;l<data__dof;l++) point_value.push_back(Vp[1+2*(l+k*data__dof)]);
  197. }
  198. }
  199. for(Long k=0;k<N_ves;k++){
  200. for(Long j=0;j<2*p1;j++){
  201. Long i0= 0;
  202. Long i1=p1;
  203. Long j0=((j+0) );
  204. Long j1=((j+1)%(2*p1));
  205. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*(p1+1)+0);
  206. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j0);
  207. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j1);
  208. poly_offset.push_back(poly_connect.size());
  209. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*(p1+1)+1);
  210. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j0);
  211. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j1);
  212. poly_offset.push_back(poly_connect.size());
  213. }
  214. for(Long i=0;i<p1;i++){
  215. for(Long j=0;j<2*p1;j++){
  216. Long i0=((i+0) );
  217. Long i1=((i+1) );
  218. Long j0=((j+0) );
  219. Long j1=((j+1)%(2*p1));
  220. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j0);
  221. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j0);
  222. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i1+j1);
  223. poly_connect.push_back((2*p1*(p1+1)+2)*k + 2*p1*i0+j1);
  224. poly_offset.push_back(poly_connect.size());
  225. }
  226. }
  227. }
  228. }
  229. Integer np = comm.Size();
  230. Integer myrank = comm.Rank();
  231. std::vector<VTKReal>& coord=point_coord;
  232. std::vector<VTKReal>& value=point_value;
  233. std::vector<int32_t>& connect=poly_connect;
  234. std::vector<int32_t>& offset=poly_offset;
  235. Long pt_cnt=coord.size()/COORD_DIM;
  236. Long poly_cnt=poly_offset.size();
  237. // Open file for writing.
  238. std::stringstream vtufname;
  239. vtufname<<fname<<"_"<<std::setfill('0')<<std::setw(6)<<myrank<<".vtp";
  240. std::ofstream vtufile;
  241. vtufile.open(vtufname.str().c_str());
  242. if(vtufile.fail()) return;
  243. bool isLittleEndian;
  244. { // Set isLittleEndian
  245. uint16_t number = 0x1;
  246. uint8_t *numPtr = (uint8_t*)&number;
  247. isLittleEndian=(numPtr[0] == 1);
  248. }
  249. // Proceed to write to file.
  250. Long data_size=0;
  251. vtufile<<"<?xml version=\"1.0\"?>\n";
  252. if(isLittleEndian) vtufile<<"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">\n";
  253. else vtufile<<"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"BigEndian\">\n";
  254. //===========================================================================
  255. vtufile<<" <PolyData>\n";
  256. vtufile<<" <Piece NumberOfPoints=\""<<pt_cnt<<"\" NumberOfVerts=\"0\" NumberOfLines=\"0\" NumberOfStrips=\"0\" NumberOfPolys=\""<<poly_cnt<<"\">\n";
  257. //---------------------------------------------------------------------------
  258. vtufile<<" <Points>\n";
  259. vtufile<<" <DataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<COORD_DIM<<"\" Name=\"Position\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  260. data_size+=sizeof(uint32_t)+coord.size()*sizeof(VTKReal);
  261. vtufile<<" </Points>\n";
  262. //---------------------------------------------------------------------------
  263. if(value.size()){ // value
  264. vtufile<<" <PointData>\n";
  265. vtufile<<" <DataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<value.size()/pt_cnt<<"\" Name=\""<<"value"<<"\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  266. data_size+=sizeof(uint32_t)+value.size()*sizeof(VTKReal);
  267. vtufile<<" </PointData>\n";
  268. }
  269. //---------------------------------------------------------------------------
  270. vtufile<<" <Polys>\n";
  271. vtufile<<" <DataArray type=\"Int32\" Name=\"connectivity\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  272. data_size+=sizeof(uint32_t)+connect.size()*sizeof(int32_t);
  273. vtufile<<" <DataArray type=\"Int32\" Name=\"offsets\" format=\"appended\" offset=\""<<data_size<<"\" />\n";
  274. data_size+=sizeof(uint32_t)+offset.size() *sizeof(int32_t);
  275. vtufile<<" </Polys>\n";
  276. //---------------------------------------------------------------------------
  277. vtufile<<" </Piece>\n";
  278. vtufile<<" </PolyData>\n";
  279. //===========================================================================
  280. vtufile<<" <AppendedData encoding=\"raw\">\n";
  281. vtufile<<" _";
  282. int32_t block_size;
  283. block_size=coord.size()*sizeof(VTKReal); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&coord [0], coord.size()*sizeof(VTKReal));
  284. if(value.size()){ // value
  285. block_size=value.size()*sizeof(VTKReal); vtufile.write((char*)&block_size, sizeof(int32_t)); vtufile.write((char*)&value [0], value.size()*sizeof(VTKReal));
  286. }
  287. 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));
  288. 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));
  289. vtufile<<"\n";
  290. vtufile<<" </AppendedData>\n";
  291. //===========================================================================
  292. vtufile<<"</VTKFile>\n";
  293. vtufile.close();
  294. if(myrank) return;
  295. std::stringstream pvtufname;
  296. pvtufname<<fname<<".pvtp";
  297. std::ofstream pvtufile;
  298. pvtufile.open(pvtufname.str().c_str());
  299. if(pvtufile.fail()) return;
  300. pvtufile<<"<?xml version=\"1.0\"?>\n";
  301. pvtufile<<"<VTKFile type=\"PPolyData\">\n";
  302. pvtufile<<" <PPolyData GhostLevel=\"0\">\n";
  303. pvtufile<<" <PPoints>\n";
  304. pvtufile<<" <PDataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<COORD_DIM<<"\" Name=\"Position\"/>\n";
  305. pvtufile<<" </PPoints>\n";
  306. if(value.size()){ // value
  307. pvtufile<<" <PPointData>\n";
  308. pvtufile<<" <PDataArray type=\"Float"<<sizeof(VTKReal)*8<<"\" NumberOfComponents=\""<<value.size()/pt_cnt<<"\" Name=\""<<"value"<<"\"/>\n";
  309. pvtufile<<" </PPointData>\n";
  310. }
  311. {
  312. // Extract filename from path.
  313. std::stringstream vtupath;
  314. vtupath<<'/'<<fname;
  315. std::string pathname = vtupath.str();
  316. auto found = pathname.find_last_of("/\\");
  317. std::string fname_ = pathname.substr(found+1);
  318. for(Integer i=0;i<np;i++) pvtufile<<" <Piece Source=\""<<fname_<<"_"<<std::setfill('0')<<std::setw(6)<<i<<".vtp\"/>\n";
  319. }
  320. pvtufile<<" </PPolyData>\n";
  321. pvtufile<<"</VTKFile>\n";
  322. pvtufile.close();
  323. }
  324. template <class Real> void SphericalHarmonics<Real>::Grid2VecSHC(const Vector<Real>& X, Long Nt, Long Np, Long p0, Vector<Real>& S, SHCArrange arrange) {
  325. Long N = X.Dim() / (Np*Nt);
  326. assert(X.Dim() == N*Np*Nt);
  327. assert(N % COORD_DIM == 0);
  328. Vector<Real> B0(N*Nt*Np);
  329. { // Set B0
  330. Vector<Real> sin_phi(Np), cos_phi(Np);
  331. for (Long i = 0; i < Np; i++) {
  332. sin_phi[i] = sin(2 * const_pi<Real>() * i / Np);
  333. cos_phi[i] = cos(2 * const_pi<Real>() * i / Np);
  334. }
  335. const auto& Y = LegendreNodes(Nt - 1);
  336. assert(Y.Dim() == Nt);
  337. Long Ngrid = Nt * Np;
  338. for (Long k = 0; k < N; k+=COORD_DIM) {
  339. for (Long i = 0; i < Nt; i++) {
  340. Real sin_theta = sqrt<Real>(1 - Y[i]*Y[i]);
  341. Real cos_theta = Y[i];
  342. Real s = 1 / sin_theta;
  343. const auto X_ = X.begin() + (k*Nt+i)*Np;
  344. auto B0_ = B0.begin() + (k*Nt+i)*Np;
  345. for (Long j = 0; j < Np; j++) {
  346. StaticArray<Real,3> in;
  347. in[0] = X_[0*Ngrid+j];
  348. in[1] = X_[1*Ngrid+j];
  349. in[2] = X_[2*Ngrid+j];
  350. StaticArray<Real,9> M;
  351. M[0] = sin_theta*cos_phi[j]; M[1] = sin_theta*sin_phi[j]; M[2] = cos_theta;
  352. M[3] = cos_theta*cos_phi[j]; M[4] = cos_theta*sin_phi[j]; M[5] =-sin_theta;
  353. M[6] = -sin_phi[j]; M[7] = cos_phi[j]; M[8] = 0;
  354. B0_[0*Ngrid+j] = ( M[0] * in[0] + M[1] * in[1] + M[2] * in[2] );
  355. B0_[1*Ngrid+j] = ( M[3] * in[0] + M[4] * in[1] + M[5] * in[2] ) * s;
  356. B0_[2*Ngrid+j] = ( M[6] * in[0] + M[7] * in[1] + M[8] * in[2] ) * s;
  357. }
  358. }
  359. }
  360. }
  361. Long p_ = p0 + 1;
  362. Long M0 = (p0+1)*(p0+1);
  363. Long M_ = (p_+1)*(p_+1);
  364. Vector<Real> B1(N*M_);
  365. Grid2SHC_(B0, Nt, Np, p_, B1);
  366. Vector<Real> B2(N*M0);
  367. const Complex<Real> imag(0,1);
  368. for (Long i=0; i<N; i+=COORD_DIM) {
  369. for (Long m=0; m<=p0; m++) {
  370. for (Long n=m; n<=p0; n++) {
  371. auto read_coeff = [&](const Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  372. Complex<Real> c;
  373. if (0<=m && m<=n && n<=p) {
  374. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  375. Long idx_imag = idx_real + (p+1-m)*N;
  376. c.real = coeff[idx_real];
  377. if (m) c.imag = coeff[idx_imag];
  378. }
  379. return c;
  380. };
  381. auto write_coeff = [&](Complex<Real> c, Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  382. if (0<=m && m<=n && n<=p) {
  383. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  384. Long idx_imag = idx_real + (p+1-m)*N;
  385. coeff[idx_real] = c.real;
  386. if (m) coeff[idx_imag] = c.imag;
  387. }
  388. };
  389. auto gr = [&](Long n, Long m) { return read_coeff(B1, i+0, p_, n, m); };
  390. auto gt = [&](Long n, Long m) { return read_coeff(B1, i+1, p_, n, m); };
  391. auto gp = [&](Long n, Long m) { return read_coeff(B1, i+2, p_, n, m); };
  392. Complex<Real> phiY, phiG, phiX;
  393. { // (phiG, phiX) <-- (gt, gp)
  394. 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); };
  395. 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); };
  396. phiY = gr(n,m);
  397. 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)));
  398. 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)));
  399. }
  400. auto phiV = (phiG * (n + 0) - phiY) * (1/(Real)(2*n + 1));
  401. auto phiW = (phiG * (n + 1) + phiY) * (1/(Real)(2*n + 1));
  402. if (n==0) {
  403. phiW = 0;
  404. phiX = 0;
  405. }
  406. write_coeff(phiV, B2, i+0, p0, n, m);
  407. write_coeff(phiW, B2, i+1, p0, n, m);
  408. write_coeff(phiX, B2, i+2, p0, n, m);
  409. }
  410. }
  411. }
  412. SHCArrange0(B2, p0, S, arrange);
  413. }
  414. template <class Real> void SphericalHarmonics<Real>::VecSHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p0, Long Nt, Long Np, Vector<Real>& X) {
  415. Vector<Real> B0;
  416. SHCArrange1(S, arrange, p0, B0);
  417. Long p_ = p0 + 1;
  418. Long M0 = (p0+1)*(p0+1);
  419. Long M_ = (p_+1)*(p_+1);
  420. Long N = B0.Dim() / M0;
  421. assert(B0.Dim() == N*M0);
  422. assert(N % COORD_DIM == 0);
  423. Vector<Real> B1(N*M_);
  424. const Complex<Real> imag(0,1);
  425. for (Long i=0; i<N; i+=COORD_DIM) {
  426. for (Long m=0; m<=p_; m++) {
  427. for (Long n=m; n<=p_; n++) {
  428. auto read_coeff = [&](const Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  429. Complex<Real> c;
  430. if (0<=m && m<=n && n<=p) {
  431. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  432. Long idx_imag = idx_real + (p+1-m)*N;
  433. c.real = coeff[idx_real];
  434. if (m) c.imag = coeff[idx_imag];
  435. }
  436. return c;
  437. };
  438. auto write_coeff = [&](Complex<Real> c, Vector<Real>& coeff, Long i, Long p, Long n, Long m) {
  439. if (0<=m && m<=n && n<=p) {
  440. Long idx_real = ((2*p-m+3)*m - (m?p+1:0))*N + (p+1-m)*i - m + n;
  441. Long idx_imag = idx_real + (p+1-m)*N;
  442. coeff[idx_real] = c.real;
  443. if (m) coeff[idx_imag] = c.imag;
  444. }
  445. };
  446. auto phiG = [&](Long n, Long m) {
  447. auto phiV = read_coeff(B0, i+0, p0, n, m);
  448. auto phiW = read_coeff(B0, i+1, p0, n, m);
  449. return phiV + phiW;
  450. };
  451. auto phiY = [&](Long n, Long m) {
  452. auto phiV = read_coeff(B0, i+0, p0, n, m);
  453. auto phiW = read_coeff(B0, i+1, p0, n, m);
  454. return phiW * n - phiV * (n + 1);
  455. };
  456. auto phiX = [&](Long n, Long m) {
  457. return read_coeff(B0, i+2, p0, n, m);
  458. };
  459. Complex<Real> gr, gt, gp;
  460. { // (gt, gp) <-- (phiG, phiX)
  461. 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); };
  462. 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); };
  463. gr = phiY(n,m);
  464. gt = phiG(n-1,m)*A(n-1,m) - phiG(n+1,m)*B(n+1,m) - imag*m*phiX(n,m);
  465. gp = phiX(n-1,m)*A(n-1,m) - phiX(n+1,m)*B(n+1,m) + imag*m*phiG(n,m);
  466. }
  467. write_coeff(gr, B1, i+0, p_, n, m);
  468. write_coeff(gt, B1, i+1, p_, n, m);
  469. write_coeff(gp, B1, i+2, p_, n, m);
  470. }
  471. }
  472. }
  473. { // Set X
  474. SHC2Grid_(B1, p_, Nt, Np, &X);
  475. Vector<Real> sin_phi(Np), cos_phi(Np);
  476. for (Long i = 0; i < Np; i++) {
  477. sin_phi[i] = sin(2 * const_pi<Real>() * i / Np);
  478. cos_phi[i] = cos(2 * const_pi<Real>() * i / Np);
  479. }
  480. const auto& Y = LegendreNodes(Nt - 1);
  481. assert(Y.Dim() == Nt);
  482. Long Ngrid = Nt * Np;
  483. for (Long k = 0; k < N; k+=COORD_DIM) {
  484. for (Long i = 0; i < Nt; i++) {
  485. Real sin_theta = sqrt<Real>(1 - Y[i]*Y[i]);
  486. Real cos_theta = Y[i];
  487. Real s = 1 / sin_theta;
  488. auto X_ = X.begin() + (k*Nt+i)*Np;
  489. for (Long j = 0; j < Np; j++) {
  490. StaticArray<Real,3> in;
  491. in[0] = X_[0*Ngrid+j];
  492. in[1] = X_[1*Ngrid+j] * s;
  493. in[2] = X_[2*Ngrid+j] * s;
  494. StaticArray<Real,9> M;
  495. M[0] = sin_theta*cos_phi[j]; M[1] = sin_theta*sin_phi[j]; M[2] = cos_theta;
  496. M[3] = cos_theta*cos_phi[j]; M[4] = cos_theta*sin_phi[j]; M[5] =-sin_theta;
  497. M[6] = -sin_phi[j]; M[7] = cos_phi[j]; M[8] = 0;
  498. X_[0*Ngrid+j] = ( M[0] * in[0] + M[3] * in[1] + M[6] * in[2] );
  499. X_[1*Ngrid+j] = ( M[1] * in[0] + M[4] * in[1] + M[7] * in[2] );
  500. X_[2*Ngrid+j] = ( M[2] * in[0] + M[5] * in[1] + M[8] * in[2] );
  501. }
  502. }
  503. }
  504. }
  505. }
  506. template <class Real> void SphericalHarmonics<Real>::VecSHCEval(const Vector<Real>& S, SHCArrange arrange, Long p0, const Vector<Real>& cos_theta_phi, Vector<Real>& X) {
  507. Long M = (p0+1) * (p0+1);
  508. Long dof;
  509. Matrix<Real> B1;
  510. { // Set B1, dof
  511. Vector<Real> B0;
  512. SHCArrange1(S, arrange, p0, B0);
  513. dof = B0.Dim() / M / COORD_DIM;
  514. assert(B0.Dim() == dof * COORD_DIM * M);
  515. B1.ReInit(dof, COORD_DIM * M);
  516. Vector<Real> B1_(B1.Dim(0) * B1.Dim(1), B1.begin(), false);
  517. SHCArrange0(B0, p0, B1_, SHCArrange::COL_MAJOR_NONZERO);
  518. }
  519. assert(B1.Dim(1) == COORD_DIM * M);
  520. assert(B1.Dim(0) == dof);
  521. Matrix<Real> SHBasis;
  522. VecSHBasisEval(p0, cos_theta_phi, SHBasis);
  523. assert(SHBasis.Dim(1) == COORD_DIM * M);
  524. Long N = SHBasis.Dim(0) / COORD_DIM;
  525. { // Set X
  526. if (X.Dim() != N * dof * COORD_DIM) X.ReInit(N * dof * COORD_DIM);
  527. for (Long k0 = 0; k0 < N; k0++) {
  528. for (Long k1 = 0; k1 < dof; k1++) {
  529. StaticArray<Real,COORD_DIM> in;
  530. for (Long j = 0; j < COORD_DIM; j++) {
  531. in[j] = 0;
  532. for (Long i = 0; i < COORD_DIM * M; i++) {
  533. in[j] += B1[k1][i] * SHBasis[k0 * COORD_DIM + j][i];
  534. }
  535. }
  536. StaticArray<Real,9> M;
  537. Real cos_theta = cos_theta_phi[k0 * 2 + 0];
  538. Real sin_theta = sqrt<Real>(1 - cos_theta * cos_theta);
  539. Real cos_phi = cos(cos_theta_phi[k0 * 2 + 1]);
  540. Real sin_phi = sin(cos_theta_phi[k0 * 2 + 1]);
  541. M[0] = sin_theta*cos_phi; M[1] = sin_theta*sin_phi; M[2] = cos_theta;
  542. M[3] = cos_theta*cos_phi; M[4] = cos_theta*sin_phi; M[5] =-sin_theta;
  543. M[6] = -sin_phi; M[7] = cos_phi; M[8] = 0;
  544. X[(k0 * dof + k1) * COORD_DIM + 0] = M[0] * in[0] + M[3] * in[1] + M[6] * in[2];
  545. X[(k0 * dof + k1) * COORD_DIM + 1] = M[1] * in[0] + M[4] * in[1] + M[7] * in[2];
  546. X[(k0 * dof + k1) * COORD_DIM + 2] = M[2] * in[0] + M[5] * in[1] + M[8] * in[2];
  547. }
  548. }
  549. }
  550. }
  551. template <class Real> void SphericalHarmonics<Real>::StokesEvalSL(const Vector<Real>& S, SHCArrange arrange, Long p0, const Vector<Real>& coord, bool interior, Vector<Real>& X) {
  552. Long M = (p0+1) * (p0+1);
  553. Long dof;
  554. Matrix<Real> B1;
  555. { // Set B1, dof
  556. Vector<Real> B0;
  557. SHCArrange1(S, arrange, p0, B0);
  558. dof = B0.Dim() / M / COORD_DIM;
  559. assert(B0.Dim() == dof * COORD_DIM * M);
  560. B1.ReInit(dof, COORD_DIM * M);
  561. Vector<Real> B1_(B1.Dim(0) * B1.Dim(1), B1.begin(), false);
  562. SHCArrange0(B0, p0, B1_, SHCArrange::COL_MAJOR_NONZERO);
  563. }
  564. assert(B1.Dim(1) == COORD_DIM * M);
  565. assert(B1.Dim(0) == dof);
  566. Long N;
  567. Matrix<Real> SHBasis;
  568. Vector<Real> R, cos_theta_phi;
  569. { // Set N, R, SHBasis
  570. N = coord.Dim() / COORD_DIM;
  571. assert(coord.Dim() == N * COORD_DIM);
  572. R.ReInit(N);
  573. cos_theta_phi.ReInit(2 * N);
  574. for (Long i = 0; i < N; i++) { // Set R, cos_theta_phi
  575. ConstIterator<Real> x = coord.begin() + i * COORD_DIM;
  576. R[i] = sqrt<Real>(x[0]*x[0] + x[1]*x[1] + x[2]*x[2]);
  577. cos_theta_phi[i * 2 + 0] = x[2] / R[i];
  578. cos_theta_phi[i * 2 + 1] = atan2(x[1], x[0]); // TODO: works only for float and double
  579. }
  580. VecSHBasisEval(p0, cos_theta_phi, SHBasis);
  581. assert(SHBasis.Dim(1) == COORD_DIM * M);
  582. assert(SHBasis.Dim(0) == N * COORD_DIM);
  583. }
  584. Matrix<Real> StokesOp(SHBasis.Dim(0), SHBasis.Dim(1));
  585. for (Long i = 0; i < N; i++) { // Set StokesOp
  586. for (Long m = 0; m <= p0; m++) {
  587. for (Long n = m; n <= p0; n++) {
  588. auto read_coeff = [&](Long n, Long m, Long k0, Long k1) {
  589. Complex<Real> c;
  590. if (0 <= m && m <= n && n <= p0 && 0 <= k0 && k0 < COORD_DIM && 0 <= k1 && k1 < COORD_DIM) {
  591. Long idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  592. c.real = SHBasis[i * COORD_DIM + k1][k0 * M + idx];
  593. if (m) {
  594. idx += (p0+1-m);
  595. c.imag = SHBasis[i * COORD_DIM + k1][k0 * M + idx];
  596. }
  597. }
  598. return c;
  599. };
  600. auto write_coeff = [&](Complex<Real> c, Long n, Long m, Long k0, Long k1) {
  601. if (0 <= m && m <= n && n <= p0 && 0 <= k0 && k0 < COORD_DIM && 0 <= k1 && k1 < COORD_DIM) {
  602. Long idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  603. StokesOp[i * COORD_DIM + k1][k0 * M + idx] = c.real;
  604. if (m) {
  605. idx += (p0+1-m);
  606. StokesOp[i * COORD_DIM + k1][k0 * M + idx] = c.imag;
  607. }
  608. }
  609. };
  610. auto Vr = read_coeff(n, m, 0, 0);
  611. auto Vt = read_coeff(n, m, 0, 1);
  612. auto Vp = read_coeff(n, m, 0, 2);
  613. auto Wr = read_coeff(n, m, 1, 0);
  614. auto Wt = read_coeff(n, m, 1, 1);
  615. auto Wp = read_coeff(n, m, 1, 2);
  616. auto Xr = read_coeff(n, m, 2, 0);
  617. auto Xt = read_coeff(n, m, 2, 1);
  618. auto Xp = read_coeff(n, m, 2, 2);
  619. Complex<Real> SVr, SVt, SVp;
  620. Complex<Real> SWr, SWt, SWp;
  621. Complex<Real> SXr, SXt, SXp;
  622. if (interior) {
  623. Real a,b;
  624. a = n / (Real)((2*n+1) * (2*n+3)) * pow<Real>(R[i], n+1);
  625. b = -(n+1) / (Real)(4*n+2) * (pow<Real>(R[i], n-1) - pow<Real>(R[i], n+1));
  626. SVr = a * Vr + b * Wr;
  627. SVt = a * Vt + b * Wt;
  628. SVp = a * Vp + b * Wp;
  629. a = (n+1) / (Real)((2*n+1) * (2*n-1)) * pow<Real>(R[i], n-1);
  630. SWr = a * Wr;
  631. SWt = a * Wt;
  632. SWp = a * Wp;
  633. a = 1 / (Real)(2*n+1) * pow<Real>(R[i], n);
  634. SXr = a * Xr;
  635. SXt = a * Xt;
  636. SXp = a * Xp;
  637. } else {
  638. Real a,b;
  639. a = n / (Real)((2*n+1) * (2*n+3)) * pow<Real>(R[i], -n-2);
  640. SVr = a * Vr;
  641. SVt = a * Vt;
  642. SVp = a * Vp;
  643. a = (n+1) / (Real)((2*n+1) * (2*n-1)) * pow<Real>(R[i], -n);
  644. b = n / (Real)(4*n+2) * (pow<Real>(R[i], -n-2) - pow<Real>(R[i], -n));
  645. SWr = a * Wr + b * Vr;
  646. SWt = a * Wt + b * Vt;
  647. SWp = a * Wp + b * Vp;
  648. a = 1 / (Real)(2*n+1) * pow<Real>(R[i], -n-1);
  649. SXr = a * Xr;
  650. SXt = a * Xt;
  651. SXp = a * Xp;
  652. }
  653. write_coeff(SVr, n, m, 0, 0);
  654. write_coeff(SVt, n, m, 0, 1);
  655. write_coeff(SVp, n, m, 0, 2);
  656. write_coeff(SWr, n, m, 1, 0);
  657. write_coeff(SWt, n, m, 1, 1);
  658. write_coeff(SWp, n, m, 1, 2);
  659. write_coeff(SXr, n, m, 2, 0);
  660. write_coeff(SXt, n, m, 2, 1);
  661. write_coeff(SXp, n, m, 2, 2);
  662. }
  663. }
  664. }
  665. { // Set X
  666. if (X.Dim() != N * dof * COORD_DIM) X.ReInit(N * dof * COORD_DIM);
  667. for (Long k0 = 0; k0 < N; k0++) {
  668. for (Long k1 = 0; k1 < dof; k1++) {
  669. StaticArray<Real,COORD_DIM> in;
  670. for (Long j = 0; j < COORD_DIM; j++) {
  671. in[j] = 0;
  672. for (Long i = 0; i < COORD_DIM * M; i++) {
  673. in[j] += B1[k1][i] * StokesOp[k0 * COORD_DIM + j][i];
  674. }
  675. }
  676. StaticArray<Real,9> M;
  677. Real cos_theta = cos_theta_phi[k0 * 2 + 0];
  678. Real sin_theta = sqrt<Real>(1 - cos_theta * cos_theta);
  679. Real cos_phi = cos(cos_theta_phi[k0 * 2 + 1]);
  680. Real sin_phi = sin(cos_theta_phi[k0 * 2 + 1]);
  681. M[0] = sin_theta*cos_phi; M[1] = sin_theta*sin_phi; M[2] = cos_theta;
  682. M[3] = cos_theta*cos_phi; M[4] = cos_theta*sin_phi; M[5] =-sin_theta;
  683. M[6] = -sin_phi; M[7] = cos_phi; M[8] = 0;
  684. X[(k0 * dof + k1) * COORD_DIM + 0] = M[0] * in[0] + M[3] * in[1] + M[6] * in[2];
  685. X[(k0 * dof + k1) * COORD_DIM + 1] = M[1] * in[0] + M[4] * in[1] + M[7] * in[2];
  686. X[(k0 * dof + k1) * COORD_DIM + 2] = M[2] * in[0] + M[5] * in[1] + M[8] * in[2];
  687. }
  688. }
  689. }
  690. }
  691. template <class Real> void SphericalHarmonics<Real>::StokesEvalDL(const Vector<Real>& S, SHCArrange arrange, Long p0, const Vector<Real>& coord, bool interior, Vector<Real>& X) {
  692. Long M = (p0+1) * (p0+1);
  693. Long dof;
  694. Matrix<Real> B1;
  695. { // Set B1, dof
  696. Vector<Real> B0;
  697. SHCArrange1(S, arrange, p0, B0);
  698. dof = B0.Dim() / M / COORD_DIM;
  699. assert(B0.Dim() == dof * COORD_DIM * M);
  700. B1.ReInit(dof, COORD_DIM * M);
  701. Vector<Real> B1_(B1.Dim(0) * B1.Dim(1), B1.begin(), false);
  702. SHCArrange0(B0, p0, B1_, SHCArrange::COL_MAJOR_NONZERO);
  703. }
  704. assert(B1.Dim(1) == COORD_DIM * M);
  705. assert(B1.Dim(0) == dof);
  706. Long N;
  707. Matrix<Real> SHBasis;
  708. Vector<Real> R, cos_theta_phi;
  709. { // Set N, R, SHBasis
  710. N = coord.Dim() / COORD_DIM;
  711. assert(coord.Dim() == N * COORD_DIM);
  712. R.ReInit(N);
  713. cos_theta_phi.ReInit(2 * N);
  714. for (Long i = 0; i < N; i++) { // Set R, cos_theta_phi
  715. ConstIterator<Real> x = coord.begin() + i * COORD_DIM;
  716. R[i] = sqrt<Real>(x[0]*x[0] + x[1]*x[1] + x[2]*x[2]);
  717. cos_theta_phi[i * 2 + 0] = x[2] / R[i];
  718. cos_theta_phi[i * 2 + 1] = atan2(x[1], x[0]); // TODO: works only for float and double
  719. }
  720. VecSHBasisEval(p0, cos_theta_phi, SHBasis);
  721. assert(SHBasis.Dim(1) == COORD_DIM * M);
  722. assert(SHBasis.Dim(0) == N * COORD_DIM);
  723. }
  724. Matrix<Real> StokesOp(SHBasis.Dim(0), SHBasis.Dim(1));
  725. for (Long i = 0; i < N; i++) { // Set StokesOp
  726. for (Long m = 0; m <= p0; m++) {
  727. for (Long n = m; n <= p0; n++) {
  728. auto read_coeff = [&](Long n, Long m, Long k0, Long k1) {
  729. Complex<Real> c;
  730. if (0 <= m && m <= n && n <= p0 && 0 <= k0 && k0 < COORD_DIM && 0 <= k1 && k1 < COORD_DIM) {
  731. Long idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  732. c.real = SHBasis[i * COORD_DIM + k1][k0 * M + idx];
  733. if (m) {
  734. idx += (p0+1-m);
  735. c.imag = SHBasis[i * COORD_DIM + k1][k0 * M + idx];
  736. }
  737. }
  738. return c;
  739. };
  740. auto write_coeff = [&](Complex<Real> c, Long n, Long m, Long k0, Long k1) {
  741. if (0 <= m && m <= n && n <= p0 && 0 <= k0 && k0 < COORD_DIM && 0 <= k1 && k1 < COORD_DIM) {
  742. Long idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  743. StokesOp[i * COORD_DIM + k1][k0 * M + idx] = c.real;
  744. if (m) {
  745. idx += (p0+1-m);
  746. StokesOp[i * COORD_DIM + k1][k0 * M + idx] = c.imag;
  747. }
  748. }
  749. };
  750. auto Vr = read_coeff(n, m, 0, 0);
  751. auto Vt = read_coeff(n, m, 0, 1);
  752. auto Vp = read_coeff(n, m, 0, 2);
  753. auto Wr = read_coeff(n, m, 1, 0);
  754. auto Wt = read_coeff(n, m, 1, 1);
  755. auto Wp = read_coeff(n, m, 1, 2);
  756. auto Xr = read_coeff(n, m, 2, 0);
  757. auto Xt = read_coeff(n, m, 2, 1);
  758. auto Xp = read_coeff(n, m, 2, 2);
  759. Complex<Real> SVr, SVt, SVp;
  760. Complex<Real> SWr, SWt, SWp;
  761. Complex<Real> SXr, SXt, SXp;
  762. if (interior) {
  763. Real a,b;
  764. a = -2*n*(n+2) / (Real)((2*n+1) * (2*n+3)) * pow<Real>(R[i], n+1);
  765. b = -(n+1)*(n+2) / (Real)(2*n+1) * (pow<Real>(R[i], n+1) - pow<Real>(R[i], n-1));
  766. SVr = a * Vr + b * Wr;
  767. SVt = a * Vt + b * Wt;
  768. SVp = a * Vp + b * Wp;
  769. a = -(2*n*n+1) / (Real)((2*n+1) * (2*n-1)) * pow<Real>(R[i], n-1);
  770. SWr = a * Wr;
  771. SWt = a * Wt;
  772. SWp = a * Wp;
  773. a = -(n+2) / (Real)(2*n+1) * pow<Real>(R[i], n);
  774. SXr = a * Xr;
  775. SXt = a * Xt;
  776. SXp = a * Xp;
  777. } else {
  778. Real a,b;
  779. a = (2*n*n+4*n+3) / (Real)((2*n+1) * (2*n+3)) * pow<Real>(R[i], -n-2);
  780. SVr = a * Vr;
  781. SVt = a * Vt;
  782. SVp = a * Vp;
  783. a = 2*(n+1)*(n-1) / (Real)((2*n+1) * (2*n-1)) * pow<Real>(R[i], -n);
  784. b = 2*n*(n-1) / (Real)(4*n+2) * (pow<Real>(R[i], -n-2) - pow<Real>(R[i], -n));
  785. SWr = a * Wr + b * Vr;
  786. SWt = a * Wt + b * Vt;
  787. SWp = a * Wp + b * Vp;
  788. a = (n-1) / (Real)(2*n+1) * pow<Real>(R[i], -n-1);
  789. SXr = a * Xr;
  790. SXt = a * Xt;
  791. SXp = a * Xp;
  792. }
  793. write_coeff(SVr, n, m, 0, 0);
  794. write_coeff(SVt, n, m, 0, 1);
  795. write_coeff(SVp, n, m, 0, 2);
  796. write_coeff(SWr, n, m, 1, 0);
  797. write_coeff(SWt, n, m, 1, 1);
  798. write_coeff(SWp, n, m, 1, 2);
  799. write_coeff(SXr, n, m, 2, 0);
  800. write_coeff(SXt, n, m, 2, 1);
  801. write_coeff(SXp, n, m, 2, 2);
  802. }
  803. }
  804. }
  805. { // Set X
  806. if (X.Dim() != N * dof * COORD_DIM) X.ReInit(N * dof * COORD_DIM);
  807. for (Long k0 = 0; k0 < N; k0++) {
  808. for (Long k1 = 0; k1 < dof; k1++) {
  809. StaticArray<Real,COORD_DIM> in;
  810. for (Long j = 0; j < COORD_DIM; j++) {
  811. in[j] = 0;
  812. for (Long i = 0; i < COORD_DIM * M; i++) {
  813. in[j] += B1[k1][i] * StokesOp[k0 * COORD_DIM + j][i];
  814. }
  815. }
  816. StaticArray<Real,9> M;
  817. Real cos_theta = cos_theta_phi[k0 * 2 + 0];
  818. Real sin_theta = sqrt<Real>(1 - cos_theta * cos_theta);
  819. Real cos_phi = cos(cos_theta_phi[k0 * 2 + 1]);
  820. Real sin_phi = sin(cos_theta_phi[k0 * 2 + 1]);
  821. M[0] = sin_theta*cos_phi; M[1] = sin_theta*sin_phi; M[2] = cos_theta;
  822. M[3] = cos_theta*cos_phi; M[4] = cos_theta*sin_phi; M[5] =-sin_theta;
  823. M[6] = -sin_phi; M[7] = cos_phi; M[8] = 0;
  824. X[(k0 * dof + k1) * COORD_DIM + 0] = M[0] * in[0] + M[3] * in[1] + M[6] * in[2];
  825. X[(k0 * dof + k1) * COORD_DIM + 1] = M[1] * in[0] + M[4] * in[1] + M[7] * in[2];
  826. X[(k0 * dof + k1) * COORD_DIM + 2] = M[2] * in[0] + M[5] * in[1] + M[8] * in[2];
  827. }
  828. }
  829. }
  830. }
  831. template <class Real> void SphericalHarmonics<Real>::Grid2SHC_(const Vector<Real>& X, Long Nt, Long Np, Long p1, Vector<Real>& B1){
  832. const auto& Mf = OpFourierInv(Np);
  833. assert(Mf.Dim(0) == Np);
  834. const std::vector<Matrix<Real>>& Ml = SphericalHarmonics<Real>::MatLegendreInv(Nt-1,p1);
  835. assert((Long)Ml.size() == p1+1);
  836. Long N = X.Dim() / (Np*Nt);
  837. assert(X.Dim() == N*Np*Nt);
  838. Vector<Real> B0((2*p1+1) * N*Nt);
  839. #pragma omp parallel
  840. { // B0 <-- Transpose(FFT(X))
  841. Integer tid=omp_get_thread_num();
  842. Integer omp_p=omp_get_num_threads();
  843. Long a=(tid+0)*N*Nt/omp_p;
  844. Long b=(tid+1)*N*Nt/omp_p;
  845. Vector<Real> buff(Mf.Dim(1));
  846. Long fft_coeff_len = std::min(buff.Dim(), 2*p1+2);
  847. Matrix<Real> B0_(2*p1+1, N*Nt, B0.begin(), false);
  848. const Matrix<Real> MX(N * Nt, Np, (Iterator<Real>)X.begin(), false);
  849. for (Long i = a; i < b; i++) {
  850. { // buff <-- FFT(Xi)
  851. const Vector<Real> Xi(Np, (Iterator<Real>)X.begin() + Np * i, false);
  852. Mf.Execute(Xi, buff);
  853. }
  854. { // B0 <-- Transpose(buff)
  855. B0_[0][i] = buff[0]; // skipping buff[1] == 0
  856. for (Long j = 2; j < fft_coeff_len; j++) B0_[j-1][i] = buff[j];
  857. for (Long j = fft_coeff_len; j < 2*p1+2; j++) B0_[j-1][i] = 0;
  858. }
  859. }
  860. }
  861. if (B1.Dim() != N*(p1+1)*(p1+1)) B1.ReInit(N*(p1+1)*(p1+1));
  862. #pragma omp parallel
  863. { // Evaluate Legendre polynomial
  864. Integer tid=omp_get_thread_num();
  865. Integer omp_p=omp_get_num_threads();
  866. Long offset0=0;
  867. Long offset1=0;
  868. for (Long i = 0; i < p1+1; i++) {
  869. Long N_ = (i==0 ? N : 2*N);
  870. Matrix<Real> Min (N_, Nt , B0.begin()+offset0, false);
  871. Matrix<Real> Mout(N_, p1+1-i, B1.begin()+offset1, false);
  872. { // Mout = Min * Ml[i] // split between threads
  873. Long a=(tid+0)*N_/omp_p;
  874. Long b=(tid+1)*N_/omp_p;
  875. if (a < b) {
  876. Matrix<Real> Min_ (b-a, Min .Dim(1), Min [a], false);
  877. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  878. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  879. }
  880. }
  881. offset0+=Min .Dim(0)*Min .Dim(1);
  882. offset1+=Mout.Dim(0)*Mout.Dim(1);
  883. }
  884. assert(offset0 == B0.Dim());
  885. assert(offset1 == B1.Dim());
  886. }
  887. B1 *= 1 / sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  888. }
  889. template <class Real> void SphericalHarmonics<Real>::SHCArrange0(const Vector<Real>& B1, Long p1, Vector<Real>& S, SHCArrange arrange){
  890. Long M = (p1+1)*(p1+1);
  891. Long N = B1.Dim() / M;
  892. assert(B1.Dim() == N*M);
  893. if (arrange == SHCArrange::ALL) { // S <-- Rearrange(B1)
  894. Long M = 2*(p1+1)*(p1+1);
  895. if(S.Dim() != N * M) S.ReInit(N * M);
  896. #pragma omp parallel
  897. { // S <-- Rearrange(B1)
  898. Integer tid=omp_get_thread_num();
  899. Integer omp_p=omp_get_num_threads();
  900. Long a=(tid+0)*N/omp_p;
  901. Long b=(tid+1)*N/omp_p;
  902. for (Long i = a; i < b; i++) {
  903. Long offset = 0;
  904. for (Long j = 0; j < p1+1; j++) {
  905. Long len = p1+1 - j;
  906. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  907. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  908. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 0;
  909. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = B_[k];
  910. offset += len;
  911. }
  912. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  913. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  914. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 1;
  915. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = B_[k];
  916. offset += len;
  917. } else {
  918. Iterator<Real> S_ = S .begin() + i*M + j*(p1+1)*2 + j*2 + 1;
  919. for (Long k = 0; k < len; k++) S_[k * (p1+1)*2] = 0;
  920. }
  921. }
  922. }
  923. }
  924. }
  925. if (arrange == SHCArrange::ROW_MAJOR) { // S <-- Rearrange(B1)
  926. Long M = (p1+1)*(p1+2);
  927. if(S.Dim() != N * M) S.ReInit(N * M);
  928. #pragma omp parallel
  929. { // S <-- Rearrange(B1)
  930. Integer tid=omp_get_thread_num();
  931. Integer omp_p=omp_get_num_threads();
  932. Long a=(tid+0)*N/omp_p;
  933. Long b=(tid+1)*N/omp_p;
  934. for (Long i = a; i < b; i++) {
  935. Long offset = 0;
  936. for (Long j = 0; j < p1+1; j++) {
  937. Long len = p1+1 - j;
  938. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  939. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  940. Iterator<Real> S_ = S .begin() + i*M + 0;
  941. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = B_[k];
  942. offset += len;
  943. }
  944. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  945. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  946. Iterator<Real> S_ = S .begin() + i*M + 1;
  947. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = B_[k];
  948. offset += len;
  949. } else {
  950. Iterator<Real> S_ = S .begin() + i*M + 1;
  951. for (Long k=0;k<len;k++) S_[(j+k)*(j+k+1) + 2*j] = 0;
  952. }
  953. }
  954. }
  955. }
  956. }
  957. if (arrange == SHCArrange::COL_MAJOR_NONZERO) { // S <-- Rearrange(B1)
  958. Long M = (p1+1)*(p1+1);
  959. if(S.Dim() != N * M) S.ReInit(N * M);
  960. #pragma omp parallel
  961. { // S <-- Rearrange(B1)
  962. Integer tid=omp_get_thread_num();
  963. Integer omp_p=omp_get_num_threads();
  964. Long a=(tid+0)*N/omp_p;
  965. Long b=(tid+1)*N/omp_p;
  966. for (Long i = a; i < b; i++) {
  967. Long offset = 0;
  968. for (Long j = 0; j < p1+1; j++) {
  969. Long len = p1+1 - j;
  970. if (1) { // Set Real(S_n^m) for m=j and n=j..p
  971. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  972. Iterator<Real> S_ = S .begin() + i*M + offset;
  973. for (Long k = 0; k < len; k++) S_[k] = B_[k];
  974. offset += len;
  975. }
  976. if (j) { // Set Imag(S_n^m) for m=j and n=j..p
  977. ConstIterator<Real> B_ = B1.begin() + i*len + N*offset;
  978. Iterator<Real> S_ = S .begin() + i*M + offset;
  979. for (Long k = 0; k < len; k++) S_[k] = B_[k];
  980. offset += len;
  981. }
  982. }
  983. }
  984. }
  985. }
  986. }
  987. template <class Real> void SphericalHarmonics<Real>::SHCArrange1(const Vector<Real>& S, SHCArrange arrange, Long p0, Vector<Real>& B0){
  988. Long M, N;
  989. { // Set M, N
  990. M = 0;
  991. if (arrange == SHCArrange::ALL) M = 2*(p0+1)*(p0+1);
  992. if (arrange == SHCArrange::ROW_MAJOR) M = (p0+1)*(p0+2);
  993. if (arrange == SHCArrange::COL_MAJOR_NONZERO) M = (p0+1)*(p0+1);
  994. if (M == 0) return;
  995. N = S.Dim() / M;
  996. assert(S.Dim() == N * M);
  997. }
  998. if (B0.Dim() != N*(p0+1)*(p0+1)) B0.ReInit(N*(p0+1)*(p0+1));
  999. if (arrange == SHCArrange::ALL) { // B0 <-- Rearrange(S)
  1000. #pragma omp parallel
  1001. { // B0 <-- Rearrange(S)
  1002. Integer tid=omp_get_thread_num();
  1003. Integer omp_p=omp_get_num_threads();
  1004. Long a=(tid+0)*N/omp_p;
  1005. Long b=(tid+1)*N/omp_p;
  1006. for (Long i = a; i < b; i++) {
  1007. Long offset = 0;
  1008. for (Long j = 0; j < p0+1; j++) {
  1009. Long len = p0+1 - j;
  1010. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  1011. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1012. ConstIterator<Real> S_ = S .begin() + i*M + j*(p0+1)*2 + j*2 + 0;
  1013. for (Long k = 0; k < len; k++) B_[k] = S_[k * (p0+1)*2];
  1014. offset += len;
  1015. }
  1016. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  1017. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1018. ConstIterator<Real> S_ = S .begin() + i*M + j*(p0+1)*2 + j*2 + 1;
  1019. for (Long k = 0; k < len; k++) B_[k] = S_[k * (p0+1)*2];
  1020. offset += len;
  1021. }
  1022. }
  1023. }
  1024. }
  1025. }
  1026. if (arrange == SHCArrange::ROW_MAJOR) { // B0 <-- Rearrange(S)
  1027. #pragma omp parallel
  1028. { // B0 <-- Rearrange(S)
  1029. Integer tid=omp_get_thread_num();
  1030. Integer omp_p=omp_get_num_threads();
  1031. Long a=(tid+0)*N/omp_p;
  1032. Long b=(tid+1)*N/omp_p;
  1033. for (Long i = a; i < b; i++) {
  1034. Long offset = 0;
  1035. for (Long j = 0; j < p0+1; j++) {
  1036. Long len = p0+1 - j;
  1037. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  1038. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1039. ConstIterator<Real> S_ = S .begin() + i*M + 0;
  1040. for (Long k=0;k<len;k++) B_[k] = S_[(j+k)*(j+k+1) + 2*j];
  1041. offset += len;
  1042. }
  1043. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  1044. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1045. ConstIterator<Real> S_ = S .begin() + i*M + 1;
  1046. for (Long k=0;k<len;k++) B_[k] = S_[(j+k)*(j+k+1) + 2*j];
  1047. offset += len;
  1048. }
  1049. }
  1050. }
  1051. }
  1052. }
  1053. if (arrange == SHCArrange::COL_MAJOR_NONZERO) { // B0 <-- Rearrange(S)
  1054. #pragma omp parallel
  1055. { // B0 <-- Rearrange(S)
  1056. Integer tid=omp_get_thread_num();
  1057. Integer omp_p=omp_get_num_threads();
  1058. Long a=(tid+0)*N/omp_p;
  1059. Long b=(tid+1)*N/omp_p;
  1060. for (Long i = a; i < b; i++) {
  1061. Long offset = 0;
  1062. for (Long j = 0; j < p0+1; j++) {
  1063. Long len = p0+1 - j;
  1064. if (1) { // Get Real(S_n^m) for m=j and n=j..p
  1065. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1066. ConstIterator<Real> S_ = S .begin() + i*M + offset;
  1067. for (Long k = 0; k < len; k++) B_[k] = S_[k];
  1068. offset += len;
  1069. }
  1070. if (j) { // Get Imag(S_n^m) for m=j and n=j..p
  1071. Iterator<Real> B_ = B0.begin() + i*len + N*offset;
  1072. ConstIterator<Real> S_ = S .begin() + i*M + offset;
  1073. for (Long k = 0; k < len; k++) B_[k] = S_[k];
  1074. offset += len;
  1075. }
  1076. }
  1077. }
  1078. }
  1079. }
  1080. }
  1081. 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){
  1082. const auto& Mf = OpFourier(Np);
  1083. assert(Mf.Dim(1) == Np);
  1084. const std::vector<Matrix<Real>>& Ml =SphericalHarmonics<Real>::MatLegendre (p0,Nt-1);
  1085. const std::vector<Matrix<Real>>& Mdl=SphericalHarmonics<Real>::MatLegendreGrad(p0,Nt-1);
  1086. assert((Long)Ml .size() == p0+1);
  1087. assert((Long)Mdl.size() == p0+1);
  1088. Long N = B0.Dim() / ((p0+1)*(p0+1));
  1089. assert(B0.Dim() == N*(p0+1)*(p0+1));
  1090. if(X && X ->Dim()!=N*Np*Nt) X ->ReInit(N*Np*Nt);
  1091. if(X_theta && X_theta->Dim()!=N*Np*Nt) X_theta->ReInit(N*Np*Nt);
  1092. if(X_phi && X_phi ->Dim()!=N*Np*Nt) X_phi ->ReInit(N*Np*Nt);
  1093. Vector<Real> B1(N*(2*p0+1)*Nt);
  1094. if(X || X_phi){
  1095. #pragma omp parallel
  1096. { // Evaluate Legendre polynomial
  1097. Integer tid=omp_get_thread_num();
  1098. Integer omp_p=omp_get_num_threads();
  1099. Long offset0=0;
  1100. Long offset1=0;
  1101. for(Long i=0;i<p0+1;i++){
  1102. Long N_ = (i==0 ? N : 2*N);
  1103. const Matrix<Real> Min (N_, p0+1-i, (Iterator<Real>)B0.begin()+offset0, false);
  1104. Matrix<Real> Mout(N_, Nt , B1.begin()+offset1, false);
  1105. { // Mout = Min * Ml[i] // split between threads
  1106. Long a=(tid+0)*N_/omp_p;
  1107. Long b=(tid+1)*N_/omp_p;
  1108. if(a<b){
  1109. const Matrix<Real> Min_ (b-a, Min .Dim(1), (Iterator<Real>)Min [a], false);
  1110. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  1111. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  1112. }
  1113. }
  1114. offset0+=Min .Dim(0)*Min .Dim(1);
  1115. offset1+=Mout.Dim(0)*Mout.Dim(1);
  1116. }
  1117. }
  1118. B1 *= sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  1119. #pragma omp parallel
  1120. { // Transpose and evaluate Fourier
  1121. Integer tid=omp_get_thread_num();
  1122. Integer omp_p=omp_get_num_threads();
  1123. Long a=(tid+0)*N*Nt/omp_p;
  1124. Long b=(tid+1)*N*Nt/omp_p;
  1125. Vector<Real> buff(Mf.Dim(0)); buff = 0;
  1126. Long fft_coeff_len = std::min(buff.Dim(), 2*p0+2);
  1127. Matrix<Real> B1_(2*p0+1, N*Nt, B1.begin(), false);
  1128. for (Long i = a; i < b; i++) {
  1129. { // buff <-- Transpose(B1)
  1130. buff[0] = B1_[0][i];
  1131. buff[1] = 0;
  1132. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  1133. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  1134. }
  1135. { // X <-- FFT(buff)
  1136. Vector<Real> Xi(Np, X->begin() + Np * i, false);
  1137. Mf.Execute(buff, Xi);
  1138. }
  1139. if(X_phi){ // Evaluate Fourier gradient
  1140. { // buff <-- Transpose(B1)
  1141. buff[0] = 0;
  1142. buff[1] = 0;
  1143. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  1144. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  1145. for (Long j = 1; j < buff.Dim()/2; j++) {
  1146. Real x = buff[2*j+0];
  1147. Real y = buff[2*j+1];
  1148. buff[2*j+0] = -j*y;
  1149. buff[2*j+1] = j*x;
  1150. }
  1151. }
  1152. { // X_phi <-- FFT(buff)
  1153. Vector<Real> Xi(Np, X_phi->begin() + Np * i, false);
  1154. Mf.Execute(buff, Xi);
  1155. }
  1156. }
  1157. }
  1158. }
  1159. }
  1160. if(X_theta){
  1161. #pragma omp parallel
  1162. { // Evaluate Legendre gradient
  1163. Integer tid=omp_get_thread_num();
  1164. Integer omp_p=omp_get_num_threads();
  1165. Long offset0=0;
  1166. Long offset1=0;
  1167. for(Long i=0;i<p0+1;i++){
  1168. Long N_ = (i==0 ? N : 2*N);
  1169. const Matrix<Real> Min (N_, p0+1-i, (Iterator<Real>)B0.begin()+offset0, false);
  1170. Matrix<Real> Mout(N_, Nt , B1.begin()+offset1, false);
  1171. { // Mout = Min * Mdl[i] // split between threads
  1172. Long a=(tid+0)*N_/omp_p;
  1173. Long b=(tid+1)*N_/omp_p;
  1174. if(a<b){
  1175. const Matrix<Real> Min_ (b-a, Min .Dim(1), (Iterator<Real>)Min [a], false);
  1176. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  1177. Matrix<Real>::GEMM(Mout_,Min_,Mdl[i]);
  1178. }
  1179. }
  1180. offset0+=Min .Dim(0)*Min .Dim(1);
  1181. offset1+=Mout.Dim(0)*Mout.Dim(1);
  1182. }
  1183. }
  1184. B1 *= sqrt<Real>(4 * const_pi<Real>() * Np); // Scaling to match Zydrunas Fortran code.
  1185. #pragma omp parallel
  1186. { // Transpose and evaluate Fourier
  1187. Integer tid=omp_get_thread_num();
  1188. Integer omp_p=omp_get_num_threads();
  1189. Long a=(tid+0)*N*Nt/omp_p;
  1190. Long b=(tid+1)*N*Nt/omp_p;
  1191. Vector<Real> buff(Mf.Dim(0)); buff = 0;
  1192. Long fft_coeff_len = std::min(buff.Dim(), 2*p0+2);
  1193. Matrix<Real> B1_(2*p0+1, N*Nt, B1.begin(), false);
  1194. for (Long i = a; i < b; i++) {
  1195. { // buff <-- Transpose(B1)
  1196. buff[0] = B1_[0][i];
  1197. buff[1] = 0;
  1198. for (Long j = 2; j < fft_coeff_len; j++) buff[j] = B1_[j-1][i];
  1199. for (Long j = fft_coeff_len; j < buff.Dim(); j++) buff[j] = 0;
  1200. }
  1201. { // Xi <-- FFT(buff)
  1202. Vector<Real> Xi(Np, X_theta->begin() + Np * i, false);
  1203. Mf.Execute(buff, Xi);
  1204. }
  1205. }
  1206. }
  1207. }
  1208. }
  1209. template <class Real> void SphericalHarmonics<Real>::LegPoly(Vector<Real>& poly_val, const Vector<Real>& X, Long degree){
  1210. Long N = X.Dim();
  1211. Long Npoly = (degree + 1) * (degree + 2) / 2;
  1212. if (poly_val.Dim() != Npoly * N) poly_val.ReInit(Npoly * N);
  1213. Real fact = 1 / sqrt<Real>(4 * const_pi<Real>());
  1214. Vector<Real> u(N);
  1215. for (Long n = 0; n < N; n++) {
  1216. u[n] = (X[n]*X[n]<1 ? sqrt<Real>(1-X[n]*X[n]) : 0);
  1217. poly_val[n] = fact;
  1218. }
  1219. Long idx = 0;
  1220. Long idx_nxt = 0;
  1221. for (Long i = 1; i <= degree; i++) {
  1222. idx_nxt += N*(degree-i+2);
  1223. Real c = sqrt<Real>((2*i+1)/(Real)(2*i));
  1224. for (Long n = 0; n < N; n++) {
  1225. poly_val[idx_nxt+n] = -poly_val[idx+n] * u[n] * c;
  1226. }
  1227. idx = idx_nxt;
  1228. }
  1229. idx = 0;
  1230. for (Long m = 0; m < degree; m++) {
  1231. for (Long n = 0; n < N; n++) {
  1232. Real pmm = 0;
  1233. Real pmmp1 = poly_val[idx+n];
  1234. for (Long ll = m + 1; ll <= degree; ll++) {
  1235. Real a = sqrt<Real>(((2*ll-1)*(2*ll+1) ) / (Real)((ll-m)*(ll+m) ));
  1236. Real b = sqrt<Real>(((2*ll+1)*(ll+m-1)*(ll-m-1)) / (Real)((ll-m)*(ll+m)*(2*ll-3)));
  1237. Real pll = X[n]*a*pmmp1 - b*pmm;
  1238. pmm = pmmp1;
  1239. pmmp1 = pll;
  1240. poly_val[idx + N*(ll-m) + n] = pll;
  1241. }
  1242. }
  1243. idx += N * (degree - m + 1);
  1244. }
  1245. }
  1246. template <class Real> void SphericalHarmonics<Real>::LegPolyDeriv(Vector<Real>& poly_val, const Vector<Real>& X, Long degree){
  1247. Long N = X.Dim();
  1248. Long Npoly = (degree + 1) * (degree + 2) / 2;
  1249. if (poly_val.Dim() != N * Npoly) poly_val.ReInit(N * Npoly);
  1250. Vector<Real> leg_poly(Npoly * N);
  1251. LegPoly(leg_poly, X, degree);
  1252. for (Long m = 0; m <= degree; m++) {
  1253. for (Long n = m; n <= degree; n++) {
  1254. ConstIterator<Real> Pn = leg_poly.begin() + N * ((degree * 2 - m + 1) * (m + 0) / 2 + n);
  1255. ConstIterator<Real> Pn_ = leg_poly.begin() + N * ((degree * 2 - m + 0) * (m + 1) / 2 + n) * (m < n);
  1256. Iterator <Real> Hn = poly_val.begin() + N * ((degree * 2 - m + 1) * (m + 0) / 2 + n);
  1257. Real c2 = sqrt<Real>(m<n ? (n+m+1)*(n-m) : 0);
  1258. for (Long i = 0; i < N; i++) {
  1259. Real c1 = (X[i]*X[i]<1 ? m/sqrt<Real>(1-X[i]*X[i]) : 0);
  1260. Hn[i] = c1*X[i]*Pn[i] + c2*Pn_[i];
  1261. }
  1262. }
  1263. }
  1264. }
  1265. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::LegendreNodes(Long p){
  1266. assert(p<SCTL_SHMAXDEG);
  1267. Vector<Real>& Qx=MatrixStore().Qx_[p];
  1268. if(!Qx.Dim()){
  1269. Vector<double> qx1(p+1);
  1270. Vector<double> qw1(p+1);
  1271. cgqf(p+1, 1, 0.0, 0.0, -1.0, 1.0, &qx1[0], &qw1[0]);
  1272. assert(typeid(Real) == typeid(double) || typeid(Real) == typeid(float)); // TODO: works only for float and double
  1273. if (Qx.Dim() != p+1) Qx.ReInit(p+1);
  1274. for (Long i = 0; i < p + 1; i++) Qx[i] = -qx1[i];
  1275. }
  1276. return Qx;
  1277. }
  1278. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::LegendreWeights(Long p){
  1279. assert(p<SCTL_SHMAXDEG);
  1280. Vector<Real>& Qw=MatrixStore().Qw_[p];
  1281. if(!Qw.Dim()){
  1282. Vector<double> qx1(p+1);
  1283. Vector<double> qw1(p+1);
  1284. cgqf(p+1, 1, 0.0, 0.0, -1.0, 1.0, &qx1[0], &qw1[0]);
  1285. assert(typeid(Real) == typeid(double) || typeid(Real) == typeid(float)); // TODO: works only for float and double
  1286. if (Qw.Dim() != p+1) Qw.ReInit(p+1);
  1287. for (Long i = 0; i < p + 1; i++) Qw[i] = qw1[i];
  1288. }
  1289. return Qw;
  1290. }
  1291. template <class Real> const Vector<Real>& SphericalHarmonics<Real>::SingularWeights(Long p1){
  1292. assert(p1<SCTL_SHMAXDEG);
  1293. Vector<Real>& Sw=MatrixStore().Sw_[p1];
  1294. if(!Sw.Dim()){
  1295. const Vector<Real>& qx1 = LegendreNodes(p1);
  1296. const Vector<Real>& qw1 = LegendreWeights(p1);
  1297. std::vector<Real> Yf(p1+1,0);
  1298. { // Set Yf
  1299. Vector<Real> x0(1); x0=1.0;
  1300. Vector<Real> alp0((p1+1)*(p1+2)/2);
  1301. LegPoly(alp0, x0, p1);
  1302. Vector<Real> alp((p1+1) * (p1+1)*(p1+2)/2);
  1303. LegPoly(alp, qx1, p1);
  1304. for(Long j=0;j<p1+1;j++){
  1305. for(Long i=0;i<p1+1;i++){
  1306. Yf[i]+=4*M_PI/(2*j+1) * alp0[j] * alp[j*(p1+1)+i];
  1307. }
  1308. }
  1309. }
  1310. Sw.ReInit(p1+1);
  1311. for(Long i=0;i<p1+1;i++){
  1312. Sw[i]=(qw1[i]*M_PI/p1)*Yf[i]/cos(acos(qx1[i])/2);
  1313. }
  1314. }
  1315. return Sw;
  1316. }
  1317. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourier(Long p0, Long p1){
  1318. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1319. Matrix<Real>& Mf =MatrixStore().Mf_ [p0*SCTL_SHMAXDEG+p1];
  1320. if(!Mf.Dim(0)){
  1321. const Real SQRT2PI=sqrt(2*M_PI);
  1322. { // Set Mf
  1323. Matrix<Real> M(2*p0,2*p1);
  1324. for(Long j=0;j<2*p1;j++){
  1325. M[0][j]=SQRT2PI*1.0;
  1326. for(Long k=1;k<p0;k++){
  1327. M[2*k-1][j]=SQRT2PI*cos(j*k*M_PI/p1);
  1328. M[2*k-0][j]=SQRT2PI*sin(j*k*M_PI/p1);
  1329. }
  1330. M[2*p0-1][j]=SQRT2PI*cos(j*p0*M_PI/p1);
  1331. }
  1332. Mf=M;
  1333. }
  1334. }
  1335. return Mf;
  1336. }
  1337. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourierInv(Long p0, Long p1){
  1338. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1339. Matrix<Real>& Mf =MatrixStore().Mfinv_ [p0*SCTL_SHMAXDEG+p1];
  1340. if(!Mf.Dim(0)){
  1341. const Real INVSQRT2PI=1.0/sqrt(2*M_PI)/p0;
  1342. { // Set Mf
  1343. Matrix<Real> M(2*p0,2*p1);
  1344. M.SetZero();
  1345. if(p1>p0) p1=p0;
  1346. for(Long j=0;j<2*p0;j++){
  1347. M[j][0]=INVSQRT2PI*0.5;
  1348. for(Long k=1;k<p1;k++){
  1349. M[j][2*k-1]=INVSQRT2PI*cos(j*k*M_PI/p0);
  1350. M[j][2*k-0]=INVSQRT2PI*sin(j*k*M_PI/p0);
  1351. }
  1352. M[j][2*p1-1]=INVSQRT2PI*cos(j*p1*M_PI/p0);
  1353. }
  1354. if(p1==p0) for(Long j=0;j<2*p0;j++) M[j][2*p1-1]*=0.5;
  1355. Mf=M;
  1356. }
  1357. }
  1358. return Mf;
  1359. }
  1360. template <class Real> const Matrix<Real>& SphericalHarmonics<Real>::MatFourierGrad(Long p0, Long p1){
  1361. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1362. Matrix<Real>& Mdf=MatrixStore().Mdf_[p0*SCTL_SHMAXDEG+p1];
  1363. if(!Mdf.Dim(0)){
  1364. const Real SQRT2PI=sqrt(2*M_PI);
  1365. { // Set Mdf_
  1366. Matrix<Real> M(2*p0,2*p1);
  1367. for(Long j=0;j<2*p1;j++){
  1368. M[0][j]=SQRT2PI*0.0;
  1369. for(Long k=1;k<p0;k++){
  1370. M[2*k-1][j]=-SQRT2PI*k*sin(j*k*M_PI/p1);
  1371. M[2*k-0][j]= SQRT2PI*k*cos(j*k*M_PI/p1);
  1372. }
  1373. M[2*p0-1][j]=-SQRT2PI*p0*sin(j*p0*M_PI/p1);
  1374. }
  1375. Mdf=M;
  1376. }
  1377. }
  1378. return Mdf;
  1379. }
  1380. template <class Real> const FFT<Real>& SphericalHarmonics<Real>::OpFourier(Long Np){
  1381. assert(Np<SCTL_SHMAXDEG);
  1382. auto& Mf =MatrixStore().Mfftinv_ [Np];
  1383. #pragma omp critical (SCTL_FFT_PLAN0)
  1384. if(!Mf.Dim(0)){
  1385. StaticArray<Long,1> fft_dim = {Np};
  1386. Mf.Setup(FFT_Type::C2R, 1, Vector<Long>(1,fft_dim,false));
  1387. }
  1388. return Mf;
  1389. }
  1390. template <class Real> const FFT<Real>& SphericalHarmonics<Real>::OpFourierInv(Long Np){
  1391. assert(Np<SCTL_SHMAXDEG);
  1392. auto& Mf =MatrixStore().Mfft_ [Np];
  1393. #pragma omp critical (SCTL_FFT_PLAN1)
  1394. if(!Mf.Dim(0)){
  1395. StaticArray<Long,1> fft_dim = {Np};
  1396. Mf.Setup(FFT_Type::R2C, 1, Vector<Long>(1,fft_dim,false));
  1397. }
  1398. return Mf;
  1399. }
  1400. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendre(Long p0, Long p1){
  1401. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1402. std::vector<Matrix<Real>>& Ml =MatrixStore().Ml_ [p0*SCTL_SHMAXDEG+p1];
  1403. if(!Ml.size()){
  1404. const Vector<Real>& qx1 = LegendreNodes(p1);
  1405. Vector<Real> alp(qx1.Dim()*(p0+1)*(p0+2)/2);
  1406. LegPoly(alp, qx1, p0);
  1407. Ml.resize(p0+1);
  1408. auto ptr = alp.begin();
  1409. for(Long i=0;i<=p0;i++){
  1410. Ml[i].ReInit(p0+1-i, qx1.Dim(), ptr);
  1411. ptr+=Ml[i].Dim(0)*Ml[i].Dim(1);
  1412. }
  1413. }
  1414. return Ml;
  1415. }
  1416. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendreInv(Long p0, Long p1){
  1417. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1418. std::vector<Matrix<Real>>& Ml =MatrixStore().Mlinv_ [p0*SCTL_SHMAXDEG+p1];
  1419. if(!Ml.size()){
  1420. const Vector<Real>& qx1 = LegendreNodes(p0);
  1421. const Vector<Real>& qw1 = LegendreWeights(p0);
  1422. Vector<Real> alp(qx1.Dim()*(p1+1)*(p1+2)/2);
  1423. LegPoly(alp, qx1, p1);
  1424. Ml.resize(p1+1);
  1425. auto ptr = alp.begin();
  1426. for(Long i=0;i<=p1;i++){
  1427. Ml[i].ReInit(qx1.Dim(), p1+1-i);
  1428. Matrix<Real> M(p1+1-i, qx1.Dim(), ptr, false);
  1429. for(Long j=0;j<p1+1-i;j++){ // Transpose and weights
  1430. for(Long k=0;k<qx1.Dim();k++){
  1431. Ml[i][k][j]=M[j][k]*qw1[k]*2*M_PI;
  1432. }
  1433. }
  1434. ptr+=Ml[i].Dim(0)*Ml[i].Dim(1);
  1435. }
  1436. }
  1437. return Ml;
  1438. }
  1439. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatLegendreGrad(Long p0, Long p1){
  1440. assert(p0<SCTL_SHMAXDEG && p1<SCTL_SHMAXDEG);
  1441. std::vector<Matrix<Real>>& Mdl=MatrixStore().Mdl_[p0*SCTL_SHMAXDEG+p1];
  1442. if(!Mdl.size()){
  1443. const Vector<Real>& qx1 = LegendreNodes(p1);
  1444. Vector<Real> alp(qx1.Dim()*(p0+1)*(p0+2)/2);
  1445. LegPolyDeriv(alp, qx1, p0);
  1446. Mdl.resize(p0+1);
  1447. auto ptr = alp.begin();
  1448. for(Long i=0;i<=p0;i++){
  1449. Mdl[i].ReInit(p0+1-i, qx1.Dim(), ptr);
  1450. ptr+=Mdl[i].Dim(0)*Mdl[i].Dim(1);
  1451. }
  1452. }
  1453. return Mdl;
  1454. }
  1455. template <class Real> void SphericalHarmonics<Real>::SHBasisEval(Long p0, const Vector<Real>& cos_theta_phi, Matrix<Real>& SHBasis) {
  1456. Long M = (p0+1) * (p0+1);
  1457. Long N = cos_theta_phi.Dim() / 2;
  1458. assert(cos_theta_phi.Dim() == N * 2);
  1459. Vector<Complex<Real>> exp_phi(N);
  1460. Matrix<Real> LegP((p0+1)*(p0+2)/2, N);
  1461. { // Set exp_phi, LegP
  1462. Vector<Real> cos_theta(N);
  1463. for (Long i = 0; i < N; i++) { // Set cos_theta, exp_phi
  1464. cos_theta[i] = cos_theta_phi[i*2+0];
  1465. exp_phi[i].real = cos(cos_theta_phi[i*2+1]);
  1466. exp_phi[i].imag = sin(cos_theta_phi[i*2+1]);
  1467. }
  1468. Vector<Real> alp(LegP.Dim(0) * LegP.Dim(1), LegP.begin(), false);
  1469. LegPoly(alp, cos_theta, p0);
  1470. }
  1471. { // Set SHBasis
  1472. SHBasis.ReInit(N, M);
  1473. Real s = 4 * sqrt<Real>(const_pi<Real>());
  1474. for (Long k0 = 0; k0 < N; k0++) {
  1475. Complex<Real> exp_phi_ = 1;
  1476. Complex<Real> exp_phi1 = exp_phi[k0];
  1477. for (Long m = 0; m <= p0; m++) {
  1478. for (Long n = m; n <= p0; n++) {
  1479. Long poly_idx = (2 * p0 - m + 1) * m / 2 + n;
  1480. Long basis_idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  1481. SHBasis[k0][basis_idx] = LegP[poly_idx][k0] * exp_phi_.real * s;
  1482. if (m) { // imaginary part
  1483. basis_idx += (p0+1-m);
  1484. SHBasis[k0][basis_idx] = -LegP[poly_idx][k0] * exp_phi_.imag * s;
  1485. } else {
  1486. SHBasis[k0][basis_idx] = SHBasis[k0][basis_idx] * 0.5;
  1487. }
  1488. }
  1489. exp_phi_ = exp_phi_ * exp_phi1;
  1490. }
  1491. }
  1492. }
  1493. assert(SHBasis.Dim(0) == N);
  1494. assert(SHBasis.Dim(1) == M);
  1495. }
  1496. template <class Real> void SphericalHarmonics<Real>::VecSHBasisEval(Long p0, const Vector<Real>& cos_theta_phi, Matrix<Real>& SHBasis) {
  1497. Long M = (p0+1) * (p0+1);
  1498. Long N = cos_theta_phi.Dim() / 2;
  1499. assert(cos_theta_phi.Dim() == N * 2);
  1500. Long p_ = p0 + 1;
  1501. Long M_ = (p_+1) * (p_+1);
  1502. Matrix<Real> Ynm(N, M_);
  1503. SHBasisEval(p_, cos_theta_phi, Ynm);
  1504. Vector<Real> cos_theta(N);
  1505. for (Long i = 0; i < N; i++) { // Set cos_theta
  1506. cos_theta[i] = cos_theta_phi[i*2+0];
  1507. }
  1508. { // Set SHBasis
  1509. SHBasis.ReInit(N * COORD_DIM, COORD_DIM * M);
  1510. SHBasis = 0;
  1511. const Complex<Real> imag(0,1);
  1512. for (Long i = 0; i < N; i++) {
  1513. Real s = 1 / sqrt<Real>(1 - cos_theta[i] * cos_theta[i]);
  1514. for (Long m = 0; m <= p0; m++) {
  1515. for (Long n = m; n <= p0; n++) {
  1516. auto Y = [&](Long n, Long m) {
  1517. Complex<Real> c;
  1518. if (0 <= m && m <= n && n <= p_) {
  1519. Long idx = (2 * p_ - m + 2) * m - (m ? p_+1 : 0) + n;
  1520. c.real = Ynm[i][idx];
  1521. if (m) {
  1522. idx += (p_+1-m);
  1523. c.imag = Ynm[i][idx];
  1524. }
  1525. }
  1526. return c;
  1527. };
  1528. auto write_coeff = [&](Complex<Real> c, Long n, Long m, Long k0, Long k1) {
  1529. if (0 <= m && m <= n && n <= p0 && 0 <= k0 && k0 < COORD_DIM && 0 <= k1 && k1 < COORD_DIM) {
  1530. Long idx = (2 * p0 - m + 2) * m - (m ? p0+1 : 0) + n;
  1531. SHBasis[i * COORD_DIM + k1][k0 * M + idx] = c.real;
  1532. if (m) {
  1533. idx += (p0+1-m);
  1534. SHBasis[i * COORD_DIM + k1][k0 * M + idx] = c.imag;
  1535. }
  1536. }
  1537. };
  1538. 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); };
  1539. 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); };
  1540. Complex<Real> AYBY = A(n,m) * Y(n+1,m) - B(n,m) * Y(n-1,m);
  1541. Complex<Real> Fv2r = Y(n,m) * (-n-1);
  1542. Complex<Real> Fw2r = Y(n,m) * n;
  1543. Complex<Real> Fx2r = 0;
  1544. Complex<Real> Fv2t = AYBY * s;
  1545. Complex<Real> Fw2t = AYBY * s;
  1546. Complex<Real> Fx2t = imag * m * Y(n,m) * s;
  1547. Complex<Real> Fv2p = -imag * m * Y(n,m) * s;
  1548. Complex<Real> Fw2p = -imag * m * Y(n,m) * s;
  1549. Complex<Real> Fx2p = AYBY * s;
  1550. write_coeff(Fv2r, n, m, 0, 0);
  1551. write_coeff(Fw2r, n, m, 1, 0);
  1552. write_coeff(Fx2r, n, m, 2, 0);
  1553. write_coeff(Fv2t, n, m, 0, 1);
  1554. write_coeff(Fw2t, n, m, 1, 1);
  1555. write_coeff(Fx2t, n, m, 2, 1);
  1556. write_coeff(Fv2p, n, m, 0, 2);
  1557. write_coeff(Fw2p, n, m, 1, 2);
  1558. write_coeff(Fx2p, n, m, 2, 2);
  1559. }
  1560. }
  1561. }
  1562. }
  1563. assert(SHBasis.Dim(0) == N * COORD_DIM);
  1564. assert(SHBasis.Dim(1) == COORD_DIM * M);
  1565. }
  1566. template <class Real> const std::vector<Matrix<Real>>& SphericalHarmonics<Real>::MatRotate(Long p0){
  1567. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1568. { // Set coeff_perm
  1569. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1570. Long itr=0;
  1571. for(Long i=0;i<2*p0;i++){
  1572. Long m=(i+1)/2;
  1573. for(Long n=m;n<=p0;n++){
  1574. coeff_perm[n][i]=itr;
  1575. itr++;
  1576. }
  1577. }
  1578. }
  1579. assert(p0<SCTL_SHMAXDEG);
  1580. std::vector<Matrix<Real>>& Mr=MatrixStore().Mr_[p0];
  1581. if(!Mr.size()){
  1582. const Real SQRT2PI=sqrt(2*M_PI);
  1583. Long Ncoef=p0*(p0+2);
  1584. Long Ngrid=2*p0*(p0+1);
  1585. Long Naleg=(p0+1)*(p0+2)/2;
  1586. Matrix<Real> Mcoord0(3,Ngrid);
  1587. const Vector<Real>& x=LegendreNodes(p0);
  1588. for(Long i=0;i<p0+1;i++){ // Set Mcoord0
  1589. for(Long j=0;j<2*p0;j++){
  1590. Mcoord0[0][i*2*p0+j]=x[i];
  1591. Mcoord0[1][i*2*p0+j]=sqrt(1-x[i]*x[i])*sin(M_PI*j/p0);
  1592. Mcoord0[2][i*2*p0+j]=sqrt(1-x[i]*x[i])*cos(M_PI*j/p0);
  1593. }
  1594. }
  1595. for(Long l=0;l<p0+1;l++){ // For each rotation angle
  1596. Matrix<Real> Mcoord1;
  1597. { // Rotate coordinates
  1598. Matrix<Real> M(COORD_DIM, COORD_DIM);
  1599. Real cos_=-x[l];
  1600. Real sin_=-sqrt(1.0-x[l]*x[l]);
  1601. M[0][0]= cos_; M[0][1]=0; M[0][2]=-sin_;
  1602. M[1][0]= 0; M[1][1]=1; M[1][2]= 0;
  1603. M[2][0]= sin_; M[2][1]=0; M[2][2]= cos_;
  1604. Mcoord1=M*Mcoord0;
  1605. }
  1606. Matrix<Real> Mleg(Naleg, Ngrid);
  1607. { // Set Mleg
  1608. const Vector<Real> Vcoord1(Mcoord1.Dim(0)*Mcoord1.Dim(1), Mcoord1.begin(), false);
  1609. Vector<Real> Vleg(Mleg.Dim(0)*Mleg.Dim(1), Mleg.begin(), false);
  1610. LegPoly(Vleg, Vcoord1, p0);
  1611. }
  1612. Vector<Real> theta(Ngrid);
  1613. for(Long i=0;i<theta.Dim();i++){ // Set theta
  1614. theta[i]=atan2(Mcoord1[1][i],Mcoord1[2][i]); // TODO: works only for float and double
  1615. }
  1616. Matrix<Real> Mcoef2grid(Ncoef, Ngrid);
  1617. { // Build Mcoef2grid
  1618. Long offset0=0;
  1619. Long offset1=0;
  1620. for(Long i=0;i<p0+1;i++){
  1621. Long len=p0+1-i;
  1622. { // P * cos
  1623. for(Long j=0;j<len;j++){
  1624. for(Long k=0;k<Ngrid;k++){
  1625. Mcoef2grid[offset1+j][k]=SQRT2PI*Mleg[offset0+j][k]*cos(i*theta[k]);
  1626. }
  1627. }
  1628. offset1+=len;
  1629. }
  1630. if(i!=0 && i!=p0){ // P * sin
  1631. for(Long j=0;j<len;j++){
  1632. for(Long k=0;k<Ngrid;k++){
  1633. Mcoef2grid[offset1+j][k]=SQRT2PI*Mleg[offset0+j][k]*sin(i*theta[k]);
  1634. }
  1635. }
  1636. offset1+=len;
  1637. }
  1638. offset0+=len;
  1639. }
  1640. assert(offset0==Naleg);
  1641. assert(offset1==Ncoef);
  1642. }
  1643. Vector<Real> Vcoef2coef(Ncoef*Ncoef);
  1644. Vector<Real> Vcoef2grid(Ncoef*Ngrid, Mcoef2grid[0], false);
  1645. Grid2SHC(Vcoef2grid, p0+1, 2*p0, p0, Vcoef2coef, SHCArrange::COL_MAJOR_NONZERO);
  1646. Matrix<Real> Mcoef2coef(Ncoef, Ncoef, Vcoef2coef.begin(), false);
  1647. for(Long n=0;n<=p0;n++){ // Create matrices for fast rotation
  1648. Matrix<Real> M(coeff_perm[n].size(),coeff_perm[n].size());
  1649. for(Long i=0;i<(Long)coeff_perm[n].size();i++){
  1650. for(Long j=0;j<(Long)coeff_perm[n].size();j++){
  1651. M[i][j]=Mcoef2coef[coeff_perm[n][i]][coeff_perm[n][j]];
  1652. }
  1653. }
  1654. Mr.push_back(M);
  1655. }
  1656. }
  1657. }
  1658. return Mr;
  1659. }
  1660. template <class Real> void SphericalHarmonics<Real>::SHC2GridTranspose(const Vector<Real>& X, Long p0, Long p1, Vector<Real>& S){
  1661. Matrix<Real> Mf =SphericalHarmonics<Real>::MatFourier(p1,p0).Transpose();
  1662. std::vector<Matrix<Real>> Ml =SphericalHarmonics<Real>::MatLegendre(p1,p0);
  1663. for(Long i=0;i<(Long)Ml.size();i++) Ml[i]=Ml[i].Transpose();
  1664. assert(p1==(Long)Ml.size()-1);
  1665. assert(p0==Mf.Dim(0)/2);
  1666. assert(p1==Mf.Dim(1)/2);
  1667. Long N=X.Dim()/(2*p0*(p0+1));
  1668. assert(N*2*p0*(p0+1)==X.Dim());
  1669. if(S.Dim()!=N*(p1*(p1+2))) S.ReInit(N*(p1*(p1+2)));
  1670. Vector<Real> B0, B1;
  1671. B0.ReInit(N* p1*(p1+2));
  1672. B1.ReInit(N*2*p1*(p0+1));
  1673. #pragma omp parallel
  1674. { // Evaluate Fourier and transpose
  1675. Integer tid=omp_get_thread_num();
  1676. Integer omp_p=omp_get_num_threads();
  1677. Long a=(tid+0)*N*(p0+1)/omp_p;
  1678. Long b=(tid+1)*N*(p0+1)/omp_p;
  1679. const Long block_size=16;
  1680. Matrix<Real> B2(block_size,2*p1);
  1681. for(Long i0=a;i0<b;i0+=block_size){
  1682. Long i1=std::min(b,i0+block_size);
  1683. const Matrix<Real> Min (i1-i0,2*p0, (Iterator<Real>)X.begin()+i0*2*p0, false);
  1684. Matrix<Real> Mout(i1-i0,2*p1, B2.begin(), false);
  1685. Matrix<Real>::GEMM(Mout, Min, Mf);
  1686. for(Long i=i0;i<i1;i++){
  1687. for(Long j=0;j<2*p1;j++){
  1688. B1[j*N*(p0+1)+i]=B2[i-i0][j];
  1689. }
  1690. }
  1691. }
  1692. }
  1693. #pragma omp parallel
  1694. { // Evaluate Legendre polynomial
  1695. Integer tid=omp_get_thread_num();
  1696. Integer omp_p=omp_get_num_threads();
  1697. Long offset0=0;
  1698. Long offset1=0;
  1699. for(Long i=0;i<p1+1;i++){
  1700. Long N0=2*N;
  1701. if(i==0 || i==p1) N0=N;
  1702. Matrix<Real> Min (N0, p0+1 , B1.begin()+offset0, false);
  1703. Matrix<Real> Mout(N0, p1+1-i, B0.begin()+offset1, false);
  1704. { // Mout = Min * Ml[i] // split between threads
  1705. Long a=(tid+0)*N0/omp_p;
  1706. Long b=(tid+1)*N0/omp_p;
  1707. if(a<b){
  1708. Matrix<Real> Min_ (b-a, Min .Dim(1), Min [a], false);
  1709. Matrix<Real> Mout_(b-a, Mout.Dim(1), Mout[a], false);
  1710. Matrix<Real>::GEMM(Mout_,Min_,Ml[i]);
  1711. }
  1712. }
  1713. offset0+=Min .Dim(0)*Min .Dim(1);
  1714. offset1+=Mout.Dim(0)*Mout.Dim(1);
  1715. }
  1716. }
  1717. #pragma omp parallel
  1718. { // S <-- Rearrange(B0)
  1719. Integer tid=omp_get_thread_num();
  1720. Integer omp_p=omp_get_num_threads();
  1721. Long a=(tid+0)*N/omp_p;
  1722. Long b=(tid+1)*N/omp_p;
  1723. for(Long i=a;i<b;i++){
  1724. Long offset=0;
  1725. for(Long j=0;j<2*p1;j++){
  1726. Long len=p1+1-(j+1)/2;
  1727. Real* B_=&B0[i*len+N*offset];
  1728. Real* S_=&S[i*p1*(p1+2)+offset];
  1729. for(Long k=0;k<len;k++) S_[k]=B_[k];
  1730. offset+=len;
  1731. }
  1732. }
  1733. }
  1734. }
  1735. template <class Real> void SphericalHarmonics<Real>::RotateAll(const Vector<Real>& S, Long p0, Long dof, Vector<Real>& S_){
  1736. const std::vector<Matrix<Real>>& Mr=MatRotate(p0);
  1737. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1738. { // Set coeff_perm
  1739. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1740. Long itr=0;
  1741. for(Long i=0;i<2*p0;i++){
  1742. Long m=(i+1)/2;
  1743. for(Long n=m;n<=p0;n++){
  1744. coeff_perm[n][i]=itr;
  1745. itr++;
  1746. }
  1747. }
  1748. }
  1749. Long Ncoef=p0*(p0+2);
  1750. Long N=S.Dim()/Ncoef/dof;
  1751. assert(N*Ncoef*dof==S.Dim());
  1752. if(S_.Dim()!=N*dof*Ncoef*p0*(p0+1)) S_.ReInit(N*dof*Ncoef*p0*(p0+1));
  1753. const Matrix<Real> S0(N*dof, Ncoef, (Iterator<Real>)S.begin(), false);
  1754. Matrix<Real> S1(N*dof*p0*(p0+1), Ncoef, S_.begin(), false);
  1755. #pragma omp parallel
  1756. { // Construct all p0*(p0+1) rotations
  1757. Integer tid=omp_get_thread_num();
  1758. Integer omp_p=omp_get_num_threads();
  1759. Matrix<Real> B0(dof*p0,Ncoef); // memory buffer
  1760. std::vector<Matrix<Real>> Bi(p0+1), Bo(p0+1); // memory buffers
  1761. for(Long i=0;i<=p0;i++){ // initialize Bi, Bo
  1762. Bi[i].ReInit(dof*p0,coeff_perm[i].size());
  1763. Bo[i].ReInit(dof*p0,coeff_perm[i].size());
  1764. }
  1765. Long a=(tid+0)*N/omp_p;
  1766. Long b=(tid+1)*N/omp_p;
  1767. for(Long i=a;i<b;i++){
  1768. for(Long d=0;d<dof;d++){
  1769. for(Long j=0;j<p0;j++){
  1770. Long offset=0;
  1771. for(Long k=0;k<p0+1;k++){
  1772. Real r[2]={cos(k*j*M_PI/p0),-sin(k*j*M_PI/p0)}; // exp(i*k*theta)
  1773. Long len=p0+1-k;
  1774. if(k!=0 && k!=p0){
  1775. for(Long l=0;l<len;l++){
  1776. Real x[2];
  1777. x[0]=S0[i*dof+d][offset+len*0+l];
  1778. x[1]=S0[i*dof+d][offset+len*1+l];
  1779. B0[j*dof+d][offset+len*0+l]=x[0]*r[0]-x[1]*r[1];
  1780. B0[j*dof+d][offset+len*1+l]=x[0]*r[1]+x[1]*r[0];
  1781. }
  1782. offset+=2*len;
  1783. }else{
  1784. for(Long l=0;l<len;l++){
  1785. B0[j*dof+d][offset+l]=S0[i*dof+d][offset+l];
  1786. }
  1787. offset+=len;
  1788. }
  1789. }
  1790. assert(offset==Ncoef);
  1791. }
  1792. }
  1793. { // Fast rotation
  1794. for(Long k=0;k<dof*p0;k++){ // forward permutation
  1795. for(Long l=0;l<=p0;l++){
  1796. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1797. Bi[l][k][j]=B0[k][coeff_perm[l][j]];
  1798. }
  1799. }
  1800. }
  1801. for(Long t=0;t<=p0;t++){
  1802. for(Long l=0;l<=p0;l++){ // mat-vec
  1803. Matrix<Real>::GEMM(Bo[l],Bi[l],Mr[t*(p0+1)+l]);
  1804. }
  1805. Matrix<Real> Mout(dof*p0,Ncoef, S1[(i*(p0+1)+t)*dof*p0], false);
  1806. for(Long k=0;k<dof*p0;k++){ // reverse permutation
  1807. for(Long l=0;l<=p0;l++){
  1808. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1809. Mout[k][coeff_perm[l][j]]=Bo[l][k][j];
  1810. }
  1811. }
  1812. }
  1813. }
  1814. }
  1815. }
  1816. }
  1817. }
  1818. template <class Real> void SphericalHarmonics<Real>::RotateTranspose(const Vector<Real>& S_, Long p0, Long dof, Vector<Real>& S){
  1819. std::vector<Matrix<Real>> Mr=MatRotate(p0);
  1820. for(Long i=0;i<(Long)Mr.size();i++) Mr[i]=Mr[i].Transpose();
  1821. std::vector<std::vector<Long>> coeff_perm(p0+1);
  1822. { // Set coeff_perm
  1823. for(Long n=0;n<=p0;n++) coeff_perm[n].resize(std::min(2*n+1,2*p0));
  1824. Long itr=0;
  1825. for(Long i=0;i<2*p0;i++){
  1826. Long m=(i+1)/2;
  1827. for(Long n=m;n<=p0;n++){
  1828. coeff_perm[n][i]=itr;
  1829. itr++;
  1830. }
  1831. }
  1832. }
  1833. Long Ncoef=p0*(p0+2);
  1834. Long N=S_.Dim()/Ncoef/dof/(p0*(p0+1));
  1835. assert(N*Ncoef*dof*(p0*(p0+1))==S_.Dim());
  1836. if(S.Dim()!=N*dof*Ncoef*p0*(p0+1)) S.ReInit(N*dof*Ncoef*p0*(p0+1));
  1837. Matrix<Real> S0(N*dof*p0*(p0+1), Ncoef, S.begin(), false);
  1838. const Matrix<Real> S1(N*dof*p0*(p0+1), Ncoef, (Iterator<Real>)S_.begin(), false);
  1839. #pragma omp parallel
  1840. { // Transpose all p0*(p0+1) rotations
  1841. Integer tid=omp_get_thread_num();
  1842. Integer omp_p=omp_get_num_threads();
  1843. Matrix<Real> B0(dof*p0,Ncoef); // memory buffer
  1844. std::vector<Matrix<Real>> Bi(p0+1), Bo(p0+1); // memory buffers
  1845. for(Long i=0;i<=p0;i++){ // initialize Bi, Bo
  1846. Bi[i].ReInit(dof*p0,coeff_perm[i].size());
  1847. Bo[i].ReInit(dof*p0,coeff_perm[i].size());
  1848. }
  1849. Long a=(tid+0)*N/omp_p;
  1850. Long b=(tid+1)*N/omp_p;
  1851. for(Long i=a;i<b;i++){
  1852. for(Long t=0;t<p0+1;t++){
  1853. Long idx0=(i*(p0+1)+t)*p0*dof;
  1854. { // Fast rotation
  1855. const Matrix<Real> Min(p0*dof, Ncoef, (Iterator<Real>)S1[idx0], false);
  1856. for(Long k=0;k<dof*p0;k++){ // forward permutation
  1857. for(Long l=0;l<=p0;l++){
  1858. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1859. Bi[l][k][j]=Min[k][coeff_perm[l][j]];
  1860. }
  1861. }
  1862. }
  1863. for(Long l=0;l<=p0;l++){ // mat-vec
  1864. Matrix<Real>::GEMM(Bo[l],Bi[l],Mr[t*(p0+1)+l]);
  1865. }
  1866. for(Long k=0;k<dof*p0;k++){ // reverse permutation
  1867. for(Long l=0;l<=p0;l++){
  1868. for(Long j=0;j<(Long)coeff_perm[l].size();j++){
  1869. B0[k][coeff_perm[l][j]]=Bo[l][k][j];
  1870. }
  1871. }
  1872. }
  1873. }
  1874. for(Long j=0;j<p0;j++){
  1875. for(Long d=0;d<dof;d++){
  1876. Long idx1=idx0+j*dof+d;
  1877. Long offset=0;
  1878. for(Long k=0;k<p0+1;k++){
  1879. Real r[2]={cos(k*j*M_PI/p0),sin(k*j*M_PI/p0)}; // exp(i*k*theta)
  1880. Long len=p0+1-k;
  1881. if(k!=0 && k!=p0){
  1882. for(Long l=0;l<len;l++){
  1883. Real x[2];
  1884. x[0]=B0[j*dof+d][offset+len*0+l];
  1885. x[1]=B0[j*dof+d][offset+len*1+l];
  1886. S0[idx1][offset+len*0+l]=x[0]*r[0]-x[1]*r[1];
  1887. S0[idx1][offset+len*1+l]=x[0]*r[1]+x[1]*r[0];
  1888. }
  1889. offset+=2*len;
  1890. }else{
  1891. for(Long l=0;l<len;l++){
  1892. S0[idx1][offset+l]=B0[j*dof+d][offset+l];
  1893. }
  1894. offset+=len;
  1895. }
  1896. }
  1897. assert(offset==Ncoef);
  1898. }
  1899. }
  1900. }
  1901. }
  1902. }
  1903. }
  1904. template <class Real> void SphericalHarmonics<Real>::StokesSingularInteg(const Vector<Real>& S, Long p0, Long p1, Vector<Real>* SLMatrix, Vector<Real>* DLMatrix){
  1905. Long Ngrid=2*p0*(p0+1);
  1906. Long Ncoef= p0*(p0+2);
  1907. Long Nves=S.Dim()/(Ngrid*COORD_DIM);
  1908. if(SLMatrix) SLMatrix->ReInit(Nves*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM));
  1909. if(DLMatrix) DLMatrix->ReInit(Nves*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM));
  1910. Long BLOCK_SIZE=(Long)6e9/((3*2*p1*(p1+1))*(3*2*p0*(p0+1))*2*8); // Limit memory usage to 6GB
  1911. BLOCK_SIZE=std::min<Long>(BLOCK_SIZE,omp_get_max_threads());
  1912. BLOCK_SIZE=std::max<Long>(BLOCK_SIZE,1);
  1913. for(Long a=0;a<Nves;a+=BLOCK_SIZE){
  1914. Long b=std::min(a+BLOCK_SIZE, Nves);
  1915. Vector<Real> _SLMatrix, _DLMatrix;
  1916. if(SLMatrix) _SLMatrix.ReInit((b-a)*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), SLMatrix->begin()+a*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), false);
  1917. if(DLMatrix) _DLMatrix.ReInit((b-a)*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), DLMatrix->begin()+a*(Ncoef*COORD_DIM)*(Ncoef*COORD_DIM), false);
  1918. const Vector<Real> _S ((b-a)*(Ngrid*COORD_DIM) , (Iterator<Real>)S.begin()+a*(Ngrid*COORD_DIM), false);
  1919. if(SLMatrix && DLMatrix) StokesSingularInteg_< true, true>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1920. else if(SLMatrix) StokesSingularInteg_< true, false>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1921. else if(DLMatrix) StokesSingularInteg_<false, true>(_S, p0, p1, _SLMatrix, _DLMatrix);
  1922. }
  1923. }
  1924. 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){
  1925. Profile::Tic("Rotate");
  1926. Vector<Real> S0, S;
  1927. SphericalHarmonics<Real>::Grid2SHC(X0, p0+1, 2*p0, p0, S0, SHCArrange::COL_MAJOR_NONZERO);
  1928. SphericalHarmonics<Real>::RotateAll(S0, p0, COORD_DIM, S);
  1929. Profile::Toc();
  1930. Profile::Tic("Upsample");
  1931. Vector<Real> X, X_theta, X_phi, trg;
  1932. SphericalHarmonics<Real>::SHC2Grid(S, SHCArrange::COL_MAJOR_NONZERO, p0, p1+1, 2*p1, &X, &X_theta, &X_phi);
  1933. SphericalHarmonics<Real>::SHC2Pole(S, SHCArrange::COL_MAJOR_NONZERO, p0, trg);
  1934. Profile::Toc();
  1935. Profile::Tic("Stokes");
  1936. Vector<Real> SL0, DL0;
  1937. { // Stokes kernel
  1938. //Long M0=2*p0*(p0+1);
  1939. Long M1=2*p1*(p1+1);
  1940. Long N=trg.Dim()/(2*COORD_DIM);
  1941. assert(X.Dim()==M1*COORD_DIM*N);
  1942. if(SLayer && SL0.Dim()!=N*2*6*M1) SL0.ReInit(2*N*6*M1);
  1943. if(DLayer && DL0.Dim()!=N*2*6*M1) DL0.ReInit(2*N*6*M1);
  1944. const Vector<Real>& qw=SphericalHarmonics<Real>::SingularWeights(p1);
  1945. const Real scal_const_dl = 3.0/(4.0*M_PI);
  1946. const Real scal_const_sl = 1.0/(8.0*M_PI);
  1947. static Real eps=-1;
  1948. if(eps<0){
  1949. eps=1;
  1950. while(eps*(Real)0.5+(Real)1.0>1.0) eps*=0.5;
  1951. }
  1952. #pragma omp parallel
  1953. {
  1954. Integer tid=omp_get_thread_num();
  1955. Integer omp_p=omp_get_num_threads();
  1956. Long a=(tid+0)*N/omp_p;
  1957. Long b=(tid+1)*N/omp_p;
  1958. for(Long i=a;i<b;i++){
  1959. for(Long t=0;t<2;t++){
  1960. Real tx, ty, tz;
  1961. { // Read target coordinates
  1962. tx=trg[i*2*COORD_DIM+0*2+t];
  1963. ty=trg[i*2*COORD_DIM+1*2+t];
  1964. tz=trg[i*2*COORD_DIM+2*2+t];
  1965. }
  1966. for(Long j0=0;j0<p1+1;j0++){
  1967. for(Long j1=0;j1<2*p1;j1++){
  1968. Long s=2*p1*j0+j1;
  1969. Real dx, dy, dz;
  1970. { // Compute dx, dy, dz
  1971. dx=tx-X[(i*COORD_DIM+0)*M1+s];
  1972. dy=ty-X[(i*COORD_DIM+1)*M1+s];
  1973. dz=tz-X[(i*COORD_DIM+2)*M1+s];
  1974. }
  1975. Real nx, ny, nz;
  1976. { // Compute source normal
  1977. Real x_theta=X_phi[(i*COORD_DIM+0)*M1+s];
  1978. Real y_theta=X_phi[(i*COORD_DIM+1)*M1+s];
  1979. Real z_theta=X_phi[(i*COORD_DIM+2)*M1+s];
  1980. Real x_phi=X_theta[(i*COORD_DIM+0)*M1+s];
  1981. Real y_phi=X_theta[(i*COORD_DIM+1)*M1+s];
  1982. Real z_phi=X_theta[(i*COORD_DIM+2)*M1+s];
  1983. nx=(y_theta*z_phi-z_theta*y_phi);
  1984. ny=(z_theta*x_phi-x_theta*z_phi);
  1985. nz=(x_theta*y_phi-y_theta*x_phi);
  1986. }
  1987. Real area_elem=1.0;
  1988. if(SLayer){ // Compute area_elem
  1989. area_elem=sqrt(nx*nx+ny*ny+nz*nz);
  1990. }
  1991. Real rinv, rinv2;
  1992. { // Compute rinv, rinv2
  1993. Real r2=dx*dx+dy*dy+dz*dz;
  1994. rinv=1.0/sqrt(r2);
  1995. if(r2<=eps) rinv=0;
  1996. rinv2=rinv*rinv;
  1997. }
  1998. if(DLayer){
  1999. Real rinv5=rinv2*rinv2*rinv;
  2000. Real r_dot_n_rinv5=scal_const_dl*qw[j0*t+(p1-j0)*(1-t)] * (nx*dx+ny*dy+nz*dz)*rinv5;
  2001. DL0[((i*2+t)*6+0)*M1+s]=dx*dx*r_dot_n_rinv5;
  2002. DL0[((i*2+t)*6+1)*M1+s]=dx*dy*r_dot_n_rinv5;
  2003. DL0[((i*2+t)*6+2)*M1+s]=dx*dz*r_dot_n_rinv5;
  2004. DL0[((i*2+t)*6+3)*M1+s]=dy*dy*r_dot_n_rinv5;
  2005. DL0[((i*2+t)*6+4)*M1+s]=dy*dz*r_dot_n_rinv5;
  2006. DL0[((i*2+t)*6+5)*M1+s]=dz*dz*r_dot_n_rinv5;
  2007. }
  2008. if(SLayer){
  2009. Real area_rinv =scal_const_sl*qw[j0*t+(p1-j0)*(1-t)] * area_elem*rinv;
  2010. Real area_rinv2=area_rinv*rinv2;
  2011. SL0[((i*2+t)*6+0)*M1+s]=area_rinv+dx*dx*area_rinv2;
  2012. SL0[((i*2+t)*6+1)*M1+s]= dx*dy*area_rinv2;
  2013. SL0[((i*2+t)*6+2)*M1+s]= dx*dz*area_rinv2;
  2014. SL0[((i*2+t)*6+3)*M1+s]=area_rinv+dy*dy*area_rinv2;
  2015. SL0[((i*2+t)*6+4)*M1+s]= dy*dz*area_rinv2;
  2016. SL0[((i*2+t)*6+5)*M1+s]=area_rinv+dz*dz*area_rinv2;
  2017. }
  2018. }
  2019. }
  2020. }
  2021. }
  2022. }
  2023. Profile::Add_FLOP(20*(2*p1)*(p1+1)*2*N);
  2024. if(SLayer) Profile::Add_FLOP((19+6)*(2*p1)*(p1+1)*2*N);
  2025. if(DLayer) Profile::Add_FLOP( 22 *(2*p1)*(p1+1)*2*N);
  2026. }
  2027. Profile::Toc();
  2028. Profile::Tic("UpsampleTranspose");
  2029. Vector<Real> SL1, DL1;
  2030. SphericalHarmonics<Real>::SHC2GridTranspose(SL0, p1, p0, SL1);
  2031. SphericalHarmonics<Real>::SHC2GridTranspose(DL0, p1, p0, DL1);
  2032. Profile::Toc();
  2033. Profile::Tic("RotateTranspose");
  2034. Vector<Real> SL2, DL2;
  2035. SphericalHarmonics<Real>::RotateTranspose(SL1, p0, 2*6, SL2);
  2036. SphericalHarmonics<Real>::RotateTranspose(DL1, p0, 2*6, DL2);
  2037. Profile::Toc();
  2038. Profile::Tic("Rearrange");
  2039. Vector<Real> SL3, DL3;
  2040. { // Transpose
  2041. Long Ncoef=p0*(p0+2);
  2042. Long Ngrid=2*p0*(p0+1);
  2043. { // Transpose SL2
  2044. Long N=SL2.Dim()/(6*Ncoef*Ngrid);
  2045. SL3.ReInit(N*COORD_DIM*Ncoef*COORD_DIM*Ngrid);
  2046. #pragma omp parallel
  2047. {
  2048. Integer tid=omp_get_thread_num();
  2049. Integer omp_p=omp_get_num_threads();
  2050. Matrix<Real> B(COORD_DIM*Ncoef,Ngrid*COORD_DIM);
  2051. Long a=(tid+0)*N/omp_p;
  2052. Long b=(tid+1)*N/omp_p;
  2053. for(Long i=a;i<b;i++){
  2054. Matrix<Real> M0(Ngrid*6, Ncoef, SL2.begin()+i*Ngrid*6*Ncoef, false);
  2055. for(Long k=0;k<Ncoef;k++){ // Transpose
  2056. for(Long j=0;j<Ngrid;j++){ // TODO: needs blocking
  2057. B[k+Ncoef*0][j*COORD_DIM+0]=M0[j*6+0][k];
  2058. B[k+Ncoef*1][j*COORD_DIM+0]=M0[j*6+1][k];
  2059. B[k+Ncoef*2][j*COORD_DIM+0]=M0[j*6+2][k];
  2060. B[k+Ncoef*0][j*COORD_DIM+1]=M0[j*6+1][k];
  2061. B[k+Ncoef*1][j*COORD_DIM+1]=M0[j*6+3][k];
  2062. B[k+Ncoef*2][j*COORD_DIM+1]=M0[j*6+4][k];
  2063. B[k+Ncoef*0][j*COORD_DIM+2]=M0[j*6+2][k];
  2064. B[k+Ncoef*1][j*COORD_DIM+2]=M0[j*6+4][k];
  2065. B[k+Ncoef*2][j*COORD_DIM+2]=M0[j*6+5][k];
  2066. }
  2067. }
  2068. Matrix<Real> M1(Ncoef*COORD_DIM, COORD_DIM*Ngrid, SL3.begin()+i*COORD_DIM*Ncoef*COORD_DIM*Ngrid, false);
  2069. for(Long k=0;k<B.Dim(0);k++){ // Rearrange
  2070. for(Long j0=0;j0<COORD_DIM;j0++){
  2071. for(Long j1=0;j1<p0+1;j1++){
  2072. 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];
  2073. 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];
  2074. }
  2075. }
  2076. }
  2077. }
  2078. }
  2079. }
  2080. { // Transpose DL2
  2081. Long N=DL2.Dim()/(6*Ncoef*Ngrid);
  2082. DL3.ReInit(N*COORD_DIM*Ncoef*COORD_DIM*Ngrid);
  2083. #pragma omp parallel
  2084. {
  2085. Integer tid=omp_get_thread_num();
  2086. Integer omp_p=omp_get_num_threads();
  2087. Matrix<Real> B(COORD_DIM*Ncoef,Ngrid*COORD_DIM);
  2088. Long a=(tid+0)*N/omp_p;
  2089. Long b=(tid+1)*N/omp_p;
  2090. for(Long i=a;i<b;i++){
  2091. Matrix<Real> M0(Ngrid*6, Ncoef, DL2.begin()+i*Ngrid*6*Ncoef, false);
  2092. for(Long k=0;k<Ncoef;k++){ // Transpose
  2093. for(Long j=0;j<Ngrid;j++){ // TODO: needs blocking
  2094. B[k+Ncoef*0][j*COORD_DIM+0]=M0[j*6+0][k];
  2095. B[k+Ncoef*1][j*COORD_DIM+0]=M0[j*6+1][k];
  2096. B[k+Ncoef*2][j*COORD_DIM+0]=M0[j*6+2][k];
  2097. B[k+Ncoef*0][j*COORD_DIM+1]=M0[j*6+1][k];
  2098. B[k+Ncoef*1][j*COORD_DIM+1]=M0[j*6+3][k];
  2099. B[k+Ncoef*2][j*COORD_DIM+1]=M0[j*6+4][k];
  2100. B[k+Ncoef*0][j*COORD_DIM+2]=M0[j*6+2][k];
  2101. B[k+Ncoef*1][j*COORD_DIM+2]=M0[j*6+4][k];
  2102. B[k+Ncoef*2][j*COORD_DIM+2]=M0[j*6+5][k];
  2103. }
  2104. }
  2105. Matrix<Real> M1(Ncoef*COORD_DIM, COORD_DIM*Ngrid, DL3.begin()+i*COORD_DIM*Ncoef*COORD_DIM*Ngrid, false);
  2106. for(Long k=0;k<B.Dim(0);k++){ // Rearrange
  2107. for(Long j0=0;j0<COORD_DIM;j0++){
  2108. for(Long j1=0;j1<p0+1;j1++){
  2109. 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];
  2110. 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];
  2111. }
  2112. }
  2113. }
  2114. }
  2115. }
  2116. }
  2117. }
  2118. Profile::Toc();
  2119. Profile::Tic("Grid2SHC");
  2120. SphericalHarmonics<Real>::Grid2SHC(SL3, p0+1, 2*p0, p0, SL, SHCArrange::COL_MAJOR_NONZERO);
  2121. SphericalHarmonics<Real>::Grid2SHC(DL3, p0+1, 2*p0, p0, DL, SHCArrange::COL_MAJOR_NONZERO);
  2122. Profile::Toc();
  2123. }
  2124. } // end namespace