kernel.txx 84 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379238023812382238323842385238623872388238923902391239223932394239523962397239823992400240124022403240424052406240724082409241024112412241324142415241624172418241924202421242224232424242524262427242824292430243124322433243424352436243724382439244024412442244324442445244624472448244924502451245224532454245524562457245824592460246124622463246424652466246724682469247024712472247324742475247624772478247924802481248224832484248524862487248824892490249124922493249424952496249724982499250025012502250325042505250625072508250925102511251225132514251525162517251825192520252125222523252425252526252725282529253025312532253325342535253625372538253925402541254225432544254525462547254825492550255125522553255425552556255725582559256025612562256325642565256625672568256925702571257225732574257525762577257825792580258125822583258425852586258725882589259025912592259325942595259625972598259926002601260226032604260526062607260826092610261126122613261426152616
  1. /**
  2. * \file kernel.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 12-20-2011
  5. * \brief This file contains the implementation of the struct Kernel and also the
  6. * implementation of various kernels for FMM.
  7. */
  8. #include <cmath>
  9. #include <cstdlib>
  10. #include <vector>
  11. #include <mem_mgr.hpp>
  12. #include <profile.hpp>
  13. #include <vector.hpp>
  14. #include <matrix.hpp>
  15. #include <precomp_mat.hpp>
  16. #ifdef __SSE__
  17. #include <xmmintrin.h>
  18. #endif
  19. #ifdef __SSE2__
  20. #include <emmintrin.h>
  21. #endif
  22. #ifdef __SSE3__
  23. #include <pmmintrin.h>
  24. #endif
  25. #ifdef __AVX__
  26. #include <immintrin.h>
  27. #endif
  28. #if defined(__MIC__)
  29. #include <immintrin.h>
  30. #endif
  31. namespace pvfmm{
  32. /**
  33. * \brief Constructor.
  34. */
  35. template <class T>
  36. Kernel<T>::Kernel(Ker_t poten, Ker_t dbl_poten, const char* name, int dim_, std::pair<int,int> k_dim,
  37. size_t dev_poten, size_t dev_dbl_poten){
  38. dim=dim_;
  39. ker_dim[0]=k_dim.first;
  40. ker_dim[1]=k_dim.second;
  41. ker_poten=poten;
  42. dbl_layer_poten=dbl_poten;
  43. ker_name=std::string(name);
  44. dev_ker_poten=dev_poten;
  45. dev_dbl_layer_poten=dev_dbl_poten;
  46. k_s2m=NULL;
  47. k_s2l=NULL;
  48. k_s2t=NULL;
  49. k_m2m=NULL;
  50. k_m2l=NULL;
  51. k_m2t=NULL;
  52. k_l2l=NULL;
  53. k_l2t=NULL;
  54. homogen=true;
  55. src_scal.Resize(ker_dim[0]); src_scal.SetZero();
  56. trg_scal.Resize(ker_dim[1]); trg_scal.SetZero();
  57. perm_vec.Resize(Perm_Count);
  58. for(size_t p_type=0;p_type<C_Perm;p_type++){
  59. perm_vec[p_type ]=Permutation<T>(ker_dim[0]);
  60. perm_vec[p_type+C_Perm]=Permutation<T>(ker_dim[1]);
  61. }
  62. init=false;
  63. }
  64. /**
  65. * \brief Initialize the kernel.
  66. */
  67. template <class T>
  68. void Kernel<T>::Initialize(bool verbose) const{
  69. if(init) return;
  70. init=true;
  71. T eps=1.0;
  72. while(eps+(T)1.0>1.0) eps*=0.5;
  73. if(ker_dim[0]*ker_dim[1]>0){ // Determine scal
  74. Matrix<T> M_scal(ker_dim[0],ker_dim[1]);
  75. size_t N=1024;
  76. T eps_=N*eps;
  77. T src_coord[3]={0,0,0};
  78. std::vector<T> trg_coord1(N*COORD_DIM);
  79. std::vector<T> trg_coord2(N*COORD_DIM);
  80. for(size_t i=0;i<N;i++){
  81. T x,y,z,r;
  82. do{
  83. x=(drand48()-0.5);
  84. y=(drand48()-0.5);
  85. z=(drand48()-0.5);
  86. r=sqrt(x*x+y*y+z*z);
  87. }while(r<0.25);
  88. trg_coord1[i*COORD_DIM+0]=x;
  89. trg_coord1[i*COORD_DIM+1]=y;
  90. trg_coord1[i*COORD_DIM+2]=z;
  91. }
  92. for(size_t i=0;i<N*COORD_DIM;i++){
  93. trg_coord2[i]=trg_coord1[i]*0.5;
  94. }
  95. T max_val=0;
  96. Matrix<T> M1(N,ker_dim[0]*ker_dim[1]);
  97. Matrix<T> M2(N,ker_dim[0]*ker_dim[1]);
  98. for(size_t i=0;i<N;i++){
  99. BuildMatrix(&src_coord [ 0], 1,
  100. &trg_coord1[i*COORD_DIM], 1, &(M1[i][0]));
  101. BuildMatrix(&src_coord [ 0], 1,
  102. &trg_coord2[i*COORD_DIM], 1, &(M2[i][0]));
  103. for(size_t j=0;j<ker_dim[0]*ker_dim[1];j++){
  104. max_val=std::max<T>(max_val,M1[i][j]);
  105. max_val=std::max<T>(max_val,M2[i][j]);
  106. }
  107. }
  108. for(size_t i=0;i<ker_dim[0]*ker_dim[1];i++){
  109. T dot11=0, dot12=0, dot22=0;
  110. for(size_t j=0;j<N;j++){
  111. dot11+=M1[j][i]*M1[j][i];
  112. dot12+=M1[j][i]*M2[j][i];
  113. dot22+=M2[j][i]*M2[j][i];
  114. }
  115. if(dot11>max_val*max_val*eps_ &&
  116. dot22>max_val*max_val*eps_ ){
  117. T s=dot12/dot11;
  118. M_scal[0][i]=log(s)/log(2.0);
  119. T err=sqrt(0.5*(dot22/dot11)/(s*s)-0.5);
  120. if(err>eps_){
  121. homogen=false;
  122. M_scal[0][i]=0.0;
  123. }
  124. assert(M_scal[0][i]>=0.0); // Kernel function must decay
  125. }else M_scal[0][i]=-1;
  126. }
  127. src_scal.Resize(ker_dim[0]); src_scal.SetZero();
  128. trg_scal.Resize(ker_dim[1]); trg_scal.SetZero();
  129. if(homogen){
  130. Matrix<T> b(ker_dim[0]*ker_dim[1]+1,1); b.SetZero();
  131. mem::memcopy(&b[0][0],&M_scal[0][0],ker_dim[0]*ker_dim[1]*sizeof(T));
  132. Matrix<T> M(ker_dim[0]*ker_dim[1]+1,ker_dim[0]+ker_dim[1]); M.SetZero();
  133. M[ker_dim[0]*ker_dim[1]][0]=1;
  134. for(size_t i0=0;i0<ker_dim[0];i0++)
  135. for(size_t i1=0;i1<ker_dim[1];i1++){
  136. size_t j=i0*ker_dim[1]+i1;
  137. if(b[j][0]>=0){
  138. M[j][ 0+ i0]=1;
  139. M[j][i1+ker_dim[0]]=1;
  140. }
  141. }
  142. Matrix<T> x=M.pinv()*b;
  143. for(size_t i=0;i<ker_dim[0];i++){
  144. src_scal[i]=x[i][0];
  145. }
  146. for(size_t i=0;i<ker_dim[1];i++){
  147. trg_scal[i]=x[ker_dim[0]+i][0];
  148. }
  149. for(size_t i0=0;i0<ker_dim[0];i0++)
  150. for(size_t i1=0;i1<ker_dim[1];i1++){
  151. if(M_scal[i0][i1]>=0){
  152. if(fabs(src_scal[i0]+trg_scal[i1]-M_scal[i0][i1])>eps_){
  153. homogen=false;
  154. }
  155. }
  156. }
  157. }
  158. if(!homogen){
  159. src_scal.SetZero();
  160. trg_scal.SetZero();
  161. //std::cout<<ker_name<<" not-scale-invariant\n";
  162. }
  163. }
  164. if(ker_dim[0]*ker_dim[1]>0){ // Determine symmetry
  165. size_t N=1024;
  166. T eps_=N*eps;
  167. T src_coord[3]={0,0,0};
  168. std::vector<T> trg_coord1(N*COORD_DIM);
  169. std::vector<T> trg_coord2(N*COORD_DIM);
  170. for(size_t i=0;i<N;i++){
  171. T x,y,z,r;
  172. do{
  173. x=(drand48()-0.5);
  174. y=(drand48()-0.5);
  175. z=(drand48()-0.5);
  176. r=sqrt(x*x+y*y+z*z);
  177. }while(r<0.25);
  178. trg_coord1[i*COORD_DIM+0]=x;
  179. trg_coord1[i*COORD_DIM+1]=y;
  180. trg_coord1[i*COORD_DIM+2]=z;
  181. }
  182. for(size_t p_type=0;p_type<C_Perm;p_type++){ // For each symmetry transform
  183. switch(p_type){ // Set trg_coord2
  184. case ReflecX:
  185. for(size_t i=0;i<N;i++){
  186. trg_coord2[i*COORD_DIM+0]=-trg_coord1[i*COORD_DIM+0];
  187. trg_coord2[i*COORD_DIM+1]= trg_coord1[i*COORD_DIM+1];
  188. trg_coord2[i*COORD_DIM+2]= trg_coord1[i*COORD_DIM+2];
  189. }
  190. break;
  191. case ReflecY:
  192. for(size_t i=0;i<N;i++){
  193. trg_coord2[i*COORD_DIM+0]= trg_coord1[i*COORD_DIM+0];
  194. trg_coord2[i*COORD_DIM+1]=-trg_coord1[i*COORD_DIM+1];
  195. trg_coord2[i*COORD_DIM+2]= trg_coord1[i*COORD_DIM+2];
  196. }
  197. break;
  198. case ReflecZ:
  199. for(size_t i=0;i<N;i++){
  200. trg_coord2[i*COORD_DIM+0]= trg_coord1[i*COORD_DIM+0];
  201. trg_coord2[i*COORD_DIM+1]= trg_coord1[i*COORD_DIM+1];
  202. trg_coord2[i*COORD_DIM+2]=-trg_coord1[i*COORD_DIM+2];
  203. }
  204. break;
  205. case SwapXY:
  206. for(size_t i=0;i<N;i++){
  207. trg_coord2[i*COORD_DIM+0]= trg_coord1[i*COORD_DIM+1];
  208. trg_coord2[i*COORD_DIM+1]= trg_coord1[i*COORD_DIM+0];
  209. trg_coord2[i*COORD_DIM+2]= trg_coord1[i*COORD_DIM+2];
  210. }
  211. break;
  212. case SwapXZ:
  213. for(size_t i=0;i<N;i++){
  214. trg_coord2[i*COORD_DIM+0]= trg_coord1[i*COORD_DIM+2];
  215. trg_coord2[i*COORD_DIM+1]= trg_coord1[i*COORD_DIM+1];
  216. trg_coord2[i*COORD_DIM+2]= trg_coord1[i*COORD_DIM+0];
  217. }
  218. break;
  219. default:
  220. for(size_t i=0;i<N;i++){
  221. trg_coord2[i*COORD_DIM+0]= trg_coord1[i*COORD_DIM+0];
  222. trg_coord2[i*COORD_DIM+1]= trg_coord1[i*COORD_DIM+1];
  223. trg_coord2[i*COORD_DIM+2]= trg_coord1[i*COORD_DIM+2];
  224. }
  225. }
  226. Matrix<long long> M11, M22;
  227. {
  228. Matrix<T> M1(N,ker_dim[0]*ker_dim[1]); M1.SetZero();
  229. Matrix<T> M2(N,ker_dim[0]*ker_dim[1]); M2.SetZero();
  230. for(size_t i=0;i<N;i++){
  231. BuildMatrix(&src_coord [ 0], 1,
  232. &trg_coord1[i*COORD_DIM], 1, &(M1[i][0]));
  233. BuildMatrix(&src_coord [ 0], 1,
  234. &trg_coord2[i*COORD_DIM], 1, &(M2[i][0]));
  235. }
  236. Matrix<T> dot11(ker_dim[0]*ker_dim[1],ker_dim[0]*ker_dim[1]);dot11.SetZero();
  237. Matrix<T> dot12(ker_dim[0]*ker_dim[1],ker_dim[0]*ker_dim[1]);dot12.SetZero();
  238. Matrix<T> dot22(ker_dim[0]*ker_dim[1],ker_dim[0]*ker_dim[1]);dot22.SetZero();
  239. std::vector<T> norm1(ker_dim[0]*ker_dim[1]);
  240. std::vector<T> norm2(ker_dim[0]*ker_dim[1]);
  241. {
  242. for(size_t k=0;k<N;k++)
  243. for(size_t i=0;i<ker_dim[0]*ker_dim[1];i++)
  244. for(size_t j=0;j<ker_dim[0]*ker_dim[1];j++){
  245. dot11[i][j]+=M1[k][i]*M1[k][j];
  246. dot12[i][j]+=M1[k][i]*M2[k][j];
  247. dot22[i][j]+=M2[k][i]*M2[k][j];
  248. }
  249. for(size_t i=0;i<ker_dim[0]*ker_dim[1];i++){
  250. norm1[i]=sqrt(dot11[i][i]);
  251. norm2[i]=sqrt(dot22[i][i]);
  252. }
  253. for(size_t i=0;i<ker_dim[0]*ker_dim[1];i++)
  254. for(size_t j=0;j<ker_dim[0]*ker_dim[1];j++){
  255. dot11[i][j]/=(norm1[i]*norm1[j]);
  256. dot12[i][j]/=(norm1[i]*norm2[j]);
  257. dot22[i][j]/=(norm2[i]*norm2[j]);
  258. }
  259. }
  260. long long flag=1;
  261. M11.Resize(ker_dim[0],ker_dim[1]); M11.SetZero();
  262. M22.Resize(ker_dim[0],ker_dim[1]); M22.SetZero();
  263. for(size_t i=0;i<ker_dim[0]*ker_dim[1];i++){
  264. if(norm1[i]>eps_ && M11[0][i]==0){
  265. for(size_t j=0;j<ker_dim[0]*ker_dim[1];j++){
  266. if(fabs(norm1[i]-norm1[j])<eps_ && fabs(fabs(dot11[i][j])-1.0)<eps_){
  267. M11[0][j]=(dot11[i][j]>0?flag:-flag);
  268. }
  269. if(fabs(norm1[i]-norm2[j])<eps_ && fabs(fabs(dot12[i][j])-1.0)<eps_){
  270. M22[0][j]=(dot12[i][j]>0?flag:-flag);
  271. }
  272. }
  273. flag++;
  274. }
  275. }
  276. }
  277. Matrix<long long> P1, P2;
  278. { // P1
  279. Matrix<long long>& P=P1;
  280. Matrix<long long> M1=M11;
  281. Matrix<long long> M2=M22;
  282. for(size_t i=0;i<M1.Dim(0);i++){
  283. for(size_t j=0;j<M1.Dim(1);j++){
  284. if(M1[i][j]<0) M1[i][j]=-M1[i][j];
  285. if(M2[i][j]<0) M2[i][j]=-M2[i][j];
  286. }
  287. std::sort(&M1[i][0],&M1[i][M1.Dim(1)]);
  288. std::sort(&M2[i][0],&M2[i][M2.Dim(1)]);
  289. }
  290. P.Resize(M1.Dim(0),M1.Dim(0));
  291. for(size_t i=0;i<M1.Dim(0);i++)
  292. for(size_t j=0;j<M1.Dim(0);j++){
  293. P[i][j]=1;
  294. for(size_t k=0;k<M1.Dim(1);k++)
  295. if(M1[i][k]!=M2[j][k]){
  296. P[i][j]=0;
  297. break;
  298. }
  299. }
  300. }
  301. { // P2
  302. Matrix<long long>& P=P2;
  303. Matrix<long long> M1=M11.Transpose();
  304. Matrix<long long> M2=M22.Transpose();
  305. for(size_t i=0;i<M1.Dim(0);i++){
  306. for(size_t j=0;j<M1.Dim(1);j++){
  307. if(M1[i][j]<0) M1[i][j]=-M1[i][j];
  308. if(M2[i][j]<0) M2[i][j]=-M2[i][j];
  309. }
  310. std::sort(&M1[i][0],&M1[i][M1.Dim(1)]);
  311. std::sort(&M2[i][0],&M2[i][M2.Dim(1)]);
  312. }
  313. P.Resize(M1.Dim(0),M1.Dim(0));
  314. for(size_t i=0;i<M1.Dim(0);i++)
  315. for(size_t j=0;j<M1.Dim(0);j++){
  316. P[i][j]=1;
  317. for(size_t k=0;k<M1.Dim(1);k++)
  318. if(M1[i][k]!=M2[j][k]){
  319. P[i][j]=0;
  320. break;
  321. }
  322. }
  323. }
  324. std::vector<Permutation<long long> > P1vec, P2vec;
  325. int dbg_cnt=0;
  326. { // P1vec
  327. Matrix<long long>& Pmat=P1;
  328. std::vector<Permutation<long long> >& Pvec=P1vec;
  329. Permutation<long long> P(Pmat.Dim(0));
  330. Vector<PERM_INT_T>& perm=P.perm;
  331. perm.SetZero();
  332. // First permutation
  333. for(size_t i=0;i<P.Dim();i++)
  334. for(size_t j=0;j<P.Dim();j++){
  335. if(Pmat[i][j]){
  336. perm[i]=j;
  337. break;
  338. }
  339. }
  340. Vector<PERM_INT_T> perm_tmp;
  341. while(true){ // Next permutation
  342. perm_tmp=perm;
  343. std::sort(&perm_tmp[0],&perm_tmp[0]+perm_tmp.Dim());
  344. for(size_t i=0;i<perm_tmp.Dim();i++){
  345. if(perm_tmp[i]!=i) break;
  346. if(i==perm_tmp.Dim()-1){
  347. Pvec.push_back(P);
  348. }
  349. }
  350. bool last=false;
  351. for(size_t i=0;i<P.Dim();i++){
  352. PERM_INT_T tmp=perm[i];
  353. for(size_t j=perm[i]+1;j<P.Dim();j++){
  354. if(Pmat[i][j]){
  355. perm[i]=j;
  356. break;
  357. }
  358. }
  359. if(perm[i]>tmp) break;
  360. for(size_t j=0;j<P.Dim();j++){
  361. if(Pmat[i][j]){
  362. perm[i]=j;
  363. break;
  364. }
  365. }
  366. if(i==P.Dim()-1) last=true;
  367. }
  368. if(last) break;
  369. }
  370. }
  371. { // P2vec
  372. Matrix<long long>& Pmat=P2;
  373. std::vector<Permutation<long long> >& Pvec=P2vec;
  374. Permutation<long long> P(Pmat.Dim(0));
  375. Vector<PERM_INT_T>& perm=P.perm;
  376. perm.SetZero();
  377. // First permutation
  378. for(size_t i=0;i<P.Dim();i++)
  379. for(size_t j=0;j<P.Dim();j++){
  380. if(Pmat[i][j]){
  381. perm[i]=j;
  382. break;
  383. }
  384. }
  385. Vector<PERM_INT_T> perm_tmp;
  386. while(true){ // Next permutation
  387. perm_tmp=perm;
  388. std::sort(&perm_tmp[0],&perm_tmp[0]+perm_tmp.Dim());
  389. for(size_t i=0;i<perm_tmp.Dim();i++){
  390. if(perm_tmp[i]!=i) break;
  391. if(i==perm_tmp.Dim()-1){
  392. Pvec.push_back(P);
  393. }
  394. }
  395. bool last=false;
  396. for(size_t i=0;i<P.Dim();i++){
  397. PERM_INT_T tmp=perm[i];
  398. for(size_t j=perm[i]+1;j<P.Dim();j++){
  399. if(Pmat[i][j]){
  400. perm[i]=j;
  401. break;
  402. }
  403. }
  404. if(perm[i]>tmp) break;
  405. for(size_t j=0;j<P.Dim();j++){
  406. if(Pmat[i][j]){
  407. perm[i]=j;
  408. break;
  409. }
  410. }
  411. if(i==P.Dim()-1) last=true;
  412. }
  413. if(last) break;
  414. }
  415. }
  416. { // Find pairs which acutally work (neglect scaling)
  417. std::vector<Permutation<long long> > P1vec_, P2vec_;
  418. Matrix<long long> M1=M11;
  419. Matrix<long long> M2=M22;
  420. for(size_t i=0;i<M1.Dim(0);i++){
  421. for(size_t j=0;j<M1.Dim(1);j++){
  422. if(M1[i][j]<0) M1[i][j]=-M1[i][j];
  423. if(M2[i][j]<0) M2[i][j]=-M2[i][j];
  424. }
  425. }
  426. Matrix<long long> M;
  427. for(size_t i=0;i<P1vec.size();i++)
  428. for(size_t j=0;j<P2vec.size();j++){
  429. M=P1vec[i]*M2*P2vec[j];
  430. for(size_t k=0;k<M.Dim(0)*M.Dim(1);k++){
  431. if(M[0][k]!=M1[0][k]) break;
  432. if(k==M.Dim(0)*M.Dim(1)-1){
  433. P1vec_.push_back(P1vec[i]);
  434. P2vec_.push_back(P2vec[j]);
  435. }
  436. }
  437. }
  438. P1vec=P1vec_;
  439. P2vec=P2vec_;
  440. }
  441. Permutation<T> P1_, P2_;
  442. { // Find pairs which acutally work
  443. for(size_t k=0;k<P1vec.size();k++){
  444. Permutation<long long> P1=P1vec[k];
  445. Permutation<long long> P2=P2vec[k];
  446. Matrix<long long> M1= M11 ;
  447. Matrix<long long> M2=P1*M22*P2;
  448. Matrix<T> M(M1.Dim(0)*M1.Dim(1)+1,M1.Dim(0)+M1.Dim(1));
  449. M.SetZero(); M[M1.Dim(0)*M1.Dim(1)][0]=1.0;
  450. for(size_t i=0;i<M1.Dim(0);i++)
  451. for(size_t j=0;j<M1.Dim(1);j++){
  452. size_t k=i*M1.Dim(1)+j;
  453. M[k][ i]= M1[i][j];
  454. M[k][M1.Dim(0)+j]=-M2[i][j];
  455. }
  456. M=M.pinv();
  457. { // Construct new permutation
  458. Permutation<long long> P1_(M1.Dim(0));
  459. Permutation<long long> P2_(M1.Dim(1));
  460. for(size_t i=0;i<M1.Dim(0);i++){
  461. P1_.scal[i]=(M[i][M1.Dim(0)*M1.Dim(1)]>0?1:-1);
  462. }
  463. for(size_t i=0;i<M1.Dim(1);i++){
  464. P2_.scal[i]=(M[M1.Dim(0)+i][M1.Dim(0)*M1.Dim(1)]>0?1:-1);
  465. }
  466. P1=P1_*P1 ;
  467. P2=P2 *P2_;
  468. }
  469. bool done=true;
  470. Matrix<long long> Merr=P1*M22*P2-M11;
  471. for(size_t i=0;i<Merr.Dim(0)*Merr.Dim(1);i++){
  472. if(Merr[0][i]){
  473. done=false;
  474. break;
  475. }
  476. }
  477. if(done){
  478. P1_=Permutation<T>(P1.Dim());
  479. P2_=Permutation<T>(P2.Dim());
  480. for(size_t i=0;i<P1.Dim();i++){
  481. P1_.perm[i]=P1.perm[i];
  482. P1_.scal[i]=P1.scal[i];
  483. }
  484. for(size_t i=0;i<P2.Dim();i++){
  485. P2_.perm[i]=P2.perm[i];
  486. P2_.scal[i]=P2.scal[i];
  487. }
  488. break;
  489. }
  490. }
  491. }
  492. //std::cout<<P1_<<'\n';
  493. //std::cout<<P2_<<'\n';
  494. perm_vec[p_type ]=P1_.Transpose();
  495. perm_vec[p_type+C_Perm]=P2_;
  496. }
  497. for(size_t i=0;i<2*C_Perm;i++){
  498. if(perm_vec[i].Dim()==0){
  499. perm_vec.Resize(0);
  500. std::cout<<"no-symmetry for: "<<ker_name<<'\n';
  501. break;
  502. }
  503. }
  504. }
  505. if(verbose){ // Display kernel information
  506. std::cout<<"\n";
  507. std::cout<<"Kernel Name : "<<ker_name<<'\n';
  508. std::cout<<"Precision : "<<(double)eps<<'\n';
  509. std::cout<<"Symmetry : "<<(perm_vec.Dim()>0?"yes":"no")<<'\n';
  510. std::cout<<"Scale Invariant: "<<(homogen?"yes":"no")<<'\n';
  511. if(homogen && ker_dim[0]*ker_dim[1]>0){
  512. std::cout<<"Scaling Matrix :\n";
  513. Matrix<T> Src(ker_dim[0],1);
  514. Matrix<T> Trg(1,ker_dim[1]);
  515. for(size_t i=0;i<ker_dim[0];i++) Src[i][0]=pow(2.0,src_scal[i]);
  516. for(size_t i=0;i<ker_dim[1];i++) Trg[0][i]=pow(2.0,trg_scal[i]);
  517. std::cout<<Src*Trg;
  518. }
  519. if(ker_dim[0]*ker_dim[1]>0){ // Accuracy of multipole expansion
  520. std::cout<<"Error : ";
  521. for(T rad=1.0; rad>1.0e-2; rad*=0.5){
  522. int m=8; // multipole order
  523. std::vector<T> equiv_surf;
  524. std::vector<T> check_surf;
  525. for(int i0=0;i0<m;i0++){
  526. for(int i1=0;i1<m;i1++){
  527. for(int i2=0;i2<m;i2++){
  528. if(i0== 0 || i1== 0 || i2== 0 ||
  529. i0==m-1 || i1==m-1 || i2==m-1){
  530. // Range: [-1/3,1/3]^3
  531. T x=((T)2*i0-(m-1))/(m-1)/3;
  532. T y=((T)2*i1-(m-1))/(m-1)/3;
  533. T z=((T)2*i2-(m-1))/(m-1)/3;
  534. equiv_surf.push_back(x*RAD0*rad);
  535. equiv_surf.push_back(y*RAD0*rad);
  536. equiv_surf.push_back(z*RAD0*rad);
  537. check_surf.push_back(x*RAD1*rad);
  538. check_surf.push_back(y*RAD1*rad);
  539. check_surf.push_back(z*RAD1*rad);
  540. }
  541. }
  542. }
  543. }
  544. size_t n_equiv=equiv_surf.size()/COORD_DIM;
  545. size_t n_check=equiv_surf.size()/COORD_DIM;
  546. size_t n_src=m;
  547. size_t n_trg=m;
  548. std::vector<T> src_coord;
  549. std::vector<T> trg_coord;
  550. for(size_t i=0;i<n_src*COORD_DIM;i++){
  551. src_coord.push_back((2*drand48()-1)/3*rad);
  552. }
  553. for(size_t i=0;i<n_trg;i++){
  554. T x,y,z,r;
  555. do{
  556. x=(drand48()-0.5);
  557. y=(drand48()-0.5);
  558. z=(drand48()-0.5);
  559. r=sqrt(x*x+y*y+z*z);
  560. }while(r==0.0);
  561. trg_coord.push_back(x/r*sqrt((T)COORD_DIM)*rad);
  562. trg_coord.push_back(y/r*sqrt((T)COORD_DIM)*rad);
  563. trg_coord.push_back(z/r*sqrt((T)COORD_DIM)*rad);
  564. }
  565. Matrix<T> M_s2c(n_src*ker_dim[0],n_check*ker_dim[1]);
  566. BuildMatrix( &src_coord[0], n_src,
  567. &check_surf[0], n_check, &(M_s2c[0][0]));
  568. Matrix<T> M_e2c(n_equiv*ker_dim[0],n_check*ker_dim[1]);
  569. BuildMatrix(&equiv_surf[0], n_equiv,
  570. &check_surf[0], n_check, &(M_e2c[0][0]));
  571. Matrix<T> M_c2e=M_e2c.pinv();
  572. Matrix<T> M_e2t(n_equiv*ker_dim[0],n_trg*ker_dim[1]);
  573. BuildMatrix(&equiv_surf[0], n_equiv,
  574. &trg_coord[0], n_trg , &(M_e2t[0][0]));
  575. Matrix<T> M_s2t(n_src*ker_dim[0],n_trg*ker_dim[1]);
  576. BuildMatrix( &src_coord[0], n_src,
  577. &trg_coord[0], n_trg , &(M_s2t[0][0]));
  578. Matrix<T> M=M_s2c*M_c2e*M_e2t-M_s2t;
  579. T max_error=0, max_value=0;
  580. for(size_t i=0;i<M.Dim(0);i++)
  581. for(size_t j=0;j<M.Dim(1);j++){
  582. max_error=std::max<T>(max_error,fabs(M [i][j]));
  583. max_value=std::max<T>(max_value,fabs(M_s2t[i][j]));
  584. }
  585. std::cout<<(double)(max_error/max_value)<<' ';
  586. if(homogen) break;
  587. }
  588. std::cout<<"\n";
  589. }
  590. std::cout<<"\n";
  591. }
  592. { // Initialize auxiliary FMM kernels
  593. if(!k_s2m) k_s2m=this;
  594. if(!k_s2l) k_s2l=this;
  595. if(!k_s2t) k_s2t=this;
  596. if(!k_m2m) k_m2m=this;
  597. if(!k_m2l) k_m2l=this;
  598. if(!k_m2t) k_m2t=this;
  599. if(!k_l2l) k_l2l=this;
  600. if(!k_l2t) k_l2t=this;
  601. assert(k_s2t->ker_dim[0]==ker_dim[0]);
  602. assert(k_s2m->ker_dim[0]==k_s2l->ker_dim[0]);
  603. assert(k_s2m->ker_dim[0]==k_s2t->ker_dim[0]);
  604. assert(k_m2m->ker_dim[0]==k_m2l->ker_dim[0]);
  605. assert(k_m2m->ker_dim[0]==k_m2t->ker_dim[0]);
  606. assert(k_l2l->ker_dim[0]==k_l2t->ker_dim[0]);
  607. assert(k_s2t->ker_dim[1]==ker_dim[1]);
  608. assert(k_s2m->ker_dim[1]==k_m2m->ker_dim[1]);
  609. assert(k_s2l->ker_dim[1]==k_l2l->ker_dim[1]);
  610. assert(k_m2l->ker_dim[1]==k_l2l->ker_dim[1]);
  611. assert(k_s2t->ker_dim[1]==k_m2t->ker_dim[1]);
  612. assert(k_s2t->ker_dim[1]==k_l2t->ker_dim[1]);
  613. k_s2m->Initialize(verbose);
  614. k_s2l->Initialize(verbose);
  615. k_s2t->Initialize(verbose);
  616. k_m2m->Initialize(verbose);
  617. k_m2l->Initialize(verbose);
  618. k_m2t->Initialize(verbose);
  619. k_l2l->Initialize(verbose);
  620. k_l2t->Initialize(verbose);
  621. }
  622. }
  623. /**
  624. * \brief Compute the transformation matrix (on the source strength vector)
  625. * to get potential at target coordinates due to sources at the given
  626. * coordinates.
  627. * \param[in] r_src Coordinates of source points.
  628. * \param[in] src_cnt Number of source points.
  629. * \param[in] r_trg Coordinates of target points.
  630. * \param[in] trg_cnt Number of target points.
  631. * \param[out] k_out Output array with potential values.
  632. */
  633. template <class T>
  634. void Kernel<T>::BuildMatrix(T* r_src, int src_cnt,
  635. T* r_trg, int trg_cnt, T* k_out) const{
  636. int dim=3; //Only supporting 3D
  637. memset(k_out, 0, src_cnt*ker_dim[0]*trg_cnt*ker_dim[1]*sizeof(T));
  638. for(int i=0;i<src_cnt;i++) //TODO Optimize this.
  639. for(int j=0;j<ker_dim[0];j++){
  640. std::vector<T> v_src(ker_dim[0],0);
  641. v_src[j]=1.0;
  642. ker_poten(&r_src[i*dim], 1, &v_src[0], 1, r_trg, trg_cnt,
  643. &k_out[(i*ker_dim[0]+j)*trg_cnt*ker_dim[1]], NULL);
  644. }
  645. }
  646. ////////////////////////////////////////////////////////////////////////////////
  647. //////// LAPLACE KERNEL ////////
  648. ////////////////////////////////////////////////////////////////////////////////
  649. /**
  650. * \brief Green's function for the Poisson's equation. Kernel tensor
  651. * dimension = 1x1.
  652. */
  653. template <class T>
  654. void laplace_poten(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  655. #ifndef __MIC__
  656. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(12*dof));
  657. #endif
  658. const T OOFP = 1.0/(4.0*const_pi<T>());
  659. for(int t=0;t<trg_cnt;t++){
  660. for(int i=0;i<dof;i++){
  661. T p=0;
  662. for(int s=0;s<src_cnt;s++){
  663. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  664. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  665. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  666. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  667. if (invR!=0) invR = 1.0/sqrt(invR);
  668. p += v_src[s*dof+i]*invR;
  669. }
  670. k_out[t*dof+i] += p*OOFP;
  671. }
  672. }
  673. }
  674. template <class T>
  675. void laplace_poten_(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  676. //void laplace_poten(T* r_src_, int src_cnt, T* v_src_, int dof, T* r_trg_, int trg_cnt, T* k_out_){
  677. // int dim=3; //Only supporting 3D
  678. // T* r_src=mem::aligned_malloc<T>(src_cnt*dim);
  679. // T* r_trg=mem::aligned_malloc<T>(trg_cnt*dim);
  680. // T* v_src=mem::aligned_malloc<T>(src_cnt );
  681. // T* k_out=mem::aligned_malloc<T>(trg_cnt );
  682. // mem::memcopy(r_src,r_src_,src_cnt*dim*sizeof(T));
  683. // mem::memcopy(r_trg,r_trg_,trg_cnt*dim*sizeof(T));
  684. // mem::memcopy(v_src,v_src_,src_cnt *sizeof(T));
  685. // mem::memcopy(k_out,k_out_,trg_cnt *sizeof(T));
  686. #define EVAL_BLKSZ 32
  687. #define MAX_DOF 100
  688. //Compute source to target interactions.
  689. const T OOFP = 1.0/(4.0*const_pi<T>());
  690. if(dof==1){
  691. for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
  692. for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
  693. int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
  694. int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
  695. for(int t=t_;t<trg_blk;t++){
  696. T p=0;
  697. for(int s=s_;s<src_blk;s++){
  698. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  699. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  700. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  701. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  702. if (invR!=0) invR = 1.0/sqrt(invR);
  703. p += v_src[s]*invR;
  704. }
  705. k_out[t] += p*OOFP;
  706. }
  707. }
  708. }else if(dof==2){
  709. T p[MAX_DOF];
  710. for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
  711. for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
  712. int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
  713. int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
  714. for(int t=t_;t<trg_blk;t++){
  715. p[0]=0; p[1]=0;
  716. for(int s=s_;s<src_blk;s++){
  717. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  718. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  719. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  720. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  721. if (invR!=0) invR = 1.0/sqrt(invR);
  722. p[0] += v_src[s*dof+0]*invR;
  723. p[1] += v_src[s*dof+1]*invR;
  724. }
  725. k_out[t*dof+0] += p[0]*OOFP;
  726. k_out[t*dof+1] += p[1]*OOFP;
  727. }
  728. }
  729. }else if(dof==3){
  730. T p[MAX_DOF];
  731. for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
  732. for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
  733. int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
  734. int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
  735. for(int t=t_;t<trg_blk;t++){
  736. p[0]=0; p[1]=0; p[2]=0;
  737. for(int s=s_;s<src_blk;s++){
  738. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  739. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  740. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  741. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  742. if (invR!=0) invR = 1.0/sqrt(invR);
  743. p[0] += v_src[s*dof+0]*invR;
  744. p[1] += v_src[s*dof+1]*invR;
  745. p[2] += v_src[s*dof+2]*invR;
  746. }
  747. k_out[t*dof+0] += p[0]*OOFP;
  748. k_out[t*dof+1] += p[1]*OOFP;
  749. k_out[t*dof+2] += p[2]*OOFP;
  750. }
  751. }
  752. }else{
  753. T p[MAX_DOF];
  754. for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
  755. for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
  756. int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
  757. int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
  758. for(int t=t_;t<trg_blk;t++){
  759. for(int i=0;i<dof;i++) p[i]=0;
  760. for(int s=s_;s<src_blk;s++){
  761. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  762. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  763. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  764. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  765. if (invR!=0) invR = 1.0/sqrt(invR);
  766. for(int i=0;i<dof;i++)
  767. p[i] += v_src[s*dof+i]*invR;
  768. }
  769. for(int i=0;i<dof;i++)
  770. k_out[t*dof+i] += p[i]*OOFP;
  771. }
  772. }
  773. }
  774. #ifndef __MIC__
  775. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+2*dof));
  776. #endif
  777. #undef MAX_DOF
  778. #undef EVAL_BLKSZ
  779. // for (int t=0; t<trg_cnt; t++)
  780. // k_out_[t] += k_out[t];
  781. // mem::aligned_free(r_src);
  782. // mem::aligned_free(r_trg);
  783. // mem::aligned_free(v_src);
  784. // mem::aligned_free(k_out);
  785. }
  786. // Laplace double layer potential.
  787. template <class T>
  788. void laplace_dbl_poten(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  789. #ifndef __MIC__
  790. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(19*dof));
  791. #endif
  792. const T OOFP = -1.0/(4.0*const_pi<T>());
  793. for(int t=0;t<trg_cnt;t++){
  794. for(int i=0;i<dof;i++){
  795. T p=0;
  796. for(int s=0;s<src_cnt;s++){
  797. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  798. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  799. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  800. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  801. if (invR!=0) invR = 1.0/sqrt(invR);
  802. p = v_src[(s*dof+i)*4+3]*invR*invR*invR;
  803. k_out[t*dof+i] += p*OOFP*( dX_reg*v_src[(s*dof+i)*4+0] +
  804. dY_reg*v_src[(s*dof+i)*4+1] +
  805. dZ_reg*v_src[(s*dof+i)*4+2] );
  806. }
  807. }
  808. }
  809. }
  810. // Laplace grdient kernel.
  811. template <class T>
  812. void laplace_grad(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  813. #ifndef __MIC__
  814. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+12*dof));
  815. #endif
  816. const T OOFP = -1.0/(4.0*const_pi<T>());
  817. if(dof==1){
  818. for(int t=0;t<trg_cnt;t++){
  819. T p=0;
  820. for(int s=0;s<src_cnt;s++){
  821. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  822. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  823. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  824. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  825. if (invR!=0) invR = 1.0/sqrt(invR);
  826. p = v_src[s]*invR*invR*invR;
  827. k_out[(t)*3+0] += p*OOFP*dX_reg;
  828. k_out[(t)*3+1] += p*OOFP*dY_reg;
  829. k_out[(t)*3+2] += p*OOFP*dZ_reg;
  830. }
  831. }
  832. }else if(dof==2){
  833. for(int t=0;t<trg_cnt;t++){
  834. T p=0;
  835. for(int s=0;s<src_cnt;s++){
  836. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  837. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  838. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  839. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  840. if (invR!=0) invR = 1.0/sqrt(invR);
  841. p = v_src[s*dof+0]*invR*invR*invR;
  842. k_out[(t*dof+0)*3+0] += p*OOFP*dX_reg;
  843. k_out[(t*dof+0)*3+1] += p*OOFP*dY_reg;
  844. k_out[(t*dof+0)*3+2] += p*OOFP*dZ_reg;
  845. p = v_src[s*dof+1]*invR*invR*invR;
  846. k_out[(t*dof+1)*3+0] += p*OOFP*dX_reg;
  847. k_out[(t*dof+1)*3+1] += p*OOFP*dY_reg;
  848. k_out[(t*dof+1)*3+2] += p*OOFP*dZ_reg;
  849. }
  850. }
  851. }else if(dof==3){
  852. for(int t=0;t<trg_cnt;t++){
  853. T p=0;
  854. for(int s=0;s<src_cnt;s++){
  855. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  856. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  857. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  858. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  859. if (invR!=0) invR = 1.0/sqrt(invR);
  860. p = v_src[s*dof+0]*invR*invR*invR;
  861. k_out[(t*dof+0)*3+0] += p*OOFP*dX_reg;
  862. k_out[(t*dof+0)*3+1] += p*OOFP*dY_reg;
  863. k_out[(t*dof+0)*3+2] += p*OOFP*dZ_reg;
  864. p = v_src[s*dof+1]*invR*invR*invR;
  865. k_out[(t*dof+1)*3+0] += p*OOFP*dX_reg;
  866. k_out[(t*dof+1)*3+1] += p*OOFP*dY_reg;
  867. k_out[(t*dof+1)*3+2] += p*OOFP*dZ_reg;
  868. p = v_src[s*dof+2]*invR*invR*invR;
  869. k_out[(t*dof+2)*3+0] += p*OOFP*dX_reg;
  870. k_out[(t*dof+2)*3+1] += p*OOFP*dY_reg;
  871. k_out[(t*dof+2)*3+2] += p*OOFP*dZ_reg;
  872. }
  873. }
  874. }else{
  875. for(int t=0;t<trg_cnt;t++){
  876. for(int i=0;i<dof;i++){
  877. T p=0;
  878. for(int s=0;s<src_cnt;s++){
  879. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  880. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  881. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  882. T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  883. if (invR!=0) invR = 1.0/sqrt(invR);
  884. p = v_src[s*dof+i]*invR*invR*invR;
  885. k_out[(t*dof+i)*3+0] += p*OOFP*dX_reg;
  886. k_out[(t*dof+i)*3+1] += p*OOFP*dY_reg;
  887. k_out[(t*dof+i)*3+2] += p*OOFP*dZ_reg;
  888. }
  889. }
  890. }
  891. }
  892. }
  893. #ifndef __MIC__
  894. #ifdef USE_SSE
  895. namespace
  896. {
  897. #define IDEAL_ALIGNMENT 16
  898. #define SIMD_LEN (int)(IDEAL_ALIGNMENT / sizeof(double))
  899. #define DECL_SIMD_ALIGNED __declspec(align(IDEAL_ALIGNMENT))
  900. void laplaceSSE(
  901. const int ns,
  902. const int nt,
  903. const double *sx,
  904. const double *sy,
  905. const double *sz,
  906. const double *tx,
  907. const double *ty,
  908. const double *tz,
  909. const double *srcDen,
  910. double *trgVal)
  911. {
  912. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  913. abort();
  914. double OOFP = 1.0/(4.0*const_pi<double>());
  915. __m128d temp;
  916. double aux_arr[SIMD_LEN+1];
  917. double *tempval;
  918. // if aux_arr is misaligned
  919. if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempval = aux_arr + 1;
  920. else tempval = aux_arr;
  921. if (size_t(tempval)%IDEAL_ALIGNMENT) abort();
  922. /*! One over four pi */
  923. __m128d oofp = _mm_set1_pd (OOFP);
  924. __m128d half = _mm_set1_pd (0.5);
  925. __m128d opf = _mm_set1_pd (1.5);
  926. __m128d zero = _mm_setzero_pd ();
  927. // loop over sources
  928. int i = 0;
  929. for (; i < nt; i++) {
  930. temp = _mm_setzero_pd();
  931. __m128d txi = _mm_load1_pd (&tx[i]);
  932. __m128d tyi = _mm_load1_pd (&ty[i]);
  933. __m128d tzi = _mm_load1_pd (&tz[i]);
  934. int j = 0;
  935. // Load and calculate in groups of SIMD_LEN
  936. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  937. __m128d sxj = _mm_load_pd (&sx[j]);
  938. __m128d syj = _mm_load_pd (&sy[j]);
  939. __m128d szj = _mm_load_pd (&sz[j]);
  940. __m128d sden = _mm_set_pd (srcDen[j+1], srcDen[j]);
  941. __m128d dX, dY, dZ;
  942. __m128d dR2;
  943. __m128d S;
  944. dX = _mm_sub_pd(txi , sxj);
  945. dY = _mm_sub_pd(tyi , syj);
  946. dZ = _mm_sub_pd(tzi , szj);
  947. sxj = _mm_mul_pd(dX, dX);
  948. syj = _mm_mul_pd(dY, dY);
  949. szj = _mm_mul_pd(dZ, dZ);
  950. dR2 = _mm_add_pd(sxj, syj);
  951. dR2 = _mm_add_pd(szj, dR2);
  952. __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
  953. __m128d xhalf = _mm_mul_pd (half, dR2);
  954. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  955. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  956. __m128d S_d = _mm_cvtps_pd(S_s);
  957. // To handle the condition when src and trg coincide
  958. S_d = _mm_andnot_pd (reqzero, S_d);
  959. S = _mm_mul_pd (S_d, S_d);
  960. S = _mm_mul_pd (S, xhalf);
  961. S = _mm_sub_pd (opf, S);
  962. S = _mm_mul_pd (S, S_d);
  963. sden = _mm_mul_pd (sden, S);
  964. temp = _mm_add_pd (sden, temp);
  965. }
  966. temp = _mm_mul_pd (temp, oofp);
  967. _mm_store_pd(tempval, temp);
  968. for (int k = 0; k < SIMD_LEN; k++) {
  969. trgVal[i] += tempval[k];
  970. }
  971. for (; j < ns; j++) {
  972. double x = tx[i] - sx[j];
  973. double y = ty[i] - sy[j];
  974. double z = tz[i] - sz[j];
  975. double r2 = x*x + y*y + z*z;
  976. double r = sqrt(r2);
  977. double invdr;
  978. if (r == 0)
  979. invdr = 0;
  980. else
  981. invdr = 1/r;
  982. double den = srcDen[j];
  983. trgVal[i] += den*invdr*OOFP;
  984. }
  985. }
  986. return;
  987. }
  988. void laplaceDblSSE(
  989. const int ns,
  990. const int nt,
  991. const double *sx,
  992. const double *sy,
  993. const double *sz,
  994. const double *tx,
  995. const double *ty,
  996. const double *tz,
  997. const double *srcDen,
  998. double *trgVal)
  999. {
  1000. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1001. abort();
  1002. double OOFP = 1.0/(4.0*const_pi<double>());
  1003. __m128d temp;
  1004. double aux_arr[SIMD_LEN+1];
  1005. double *tempval;
  1006. // if aux_arr is misaligned
  1007. if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempval = aux_arr + 1;
  1008. else tempval = aux_arr;
  1009. if (size_t(tempval)%IDEAL_ALIGNMENT) abort();
  1010. /*! One over four pi */
  1011. __m128d oofp = _mm_set1_pd (OOFP);
  1012. __m128d half = _mm_set1_pd (0.5);
  1013. __m128d opf = _mm_set1_pd (1.5);
  1014. __m128d zero = _mm_setzero_pd ();
  1015. // loop over sources
  1016. int i = 0;
  1017. for (; i < nt; i++) {
  1018. temp = _mm_setzero_pd();
  1019. __m128d txi = _mm_load1_pd (&tx[i]);
  1020. __m128d tyi = _mm_load1_pd (&ty[i]);
  1021. __m128d tzi = _mm_load1_pd (&tz[i]);
  1022. int j = 0;
  1023. // Load and calculate in groups of SIMD_LEN
  1024. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1025. __m128d sxj = _mm_load_pd (&sx[j]);
  1026. __m128d syj = _mm_load_pd (&sy[j]);
  1027. __m128d szj = _mm_load_pd (&sz[j]);
  1028. __m128d snormx = _mm_set_pd (srcDen[(j+1)*4+0], srcDen[j*4+0]);
  1029. __m128d snormy = _mm_set_pd (srcDen[(j+1)*4+1], srcDen[j*4+1]);
  1030. __m128d snormz = _mm_set_pd (srcDen[(j+1)*4+2], srcDen[j*4+2]);
  1031. __m128d sden = _mm_set_pd (srcDen[(j+1)*4+3], srcDen[j*4+3]);
  1032. __m128d dX, dY, dZ;
  1033. __m128d dR2;
  1034. __m128d S;
  1035. __m128d S2;
  1036. __m128d S3;
  1037. dX = _mm_sub_pd(txi , sxj);
  1038. dY = _mm_sub_pd(tyi , syj);
  1039. dZ = _mm_sub_pd(tzi , szj);
  1040. sxj = _mm_mul_pd(dX, dX);
  1041. syj = _mm_mul_pd(dY, dY);
  1042. szj = _mm_mul_pd(dZ, dZ);
  1043. dR2 = _mm_add_pd(sxj, syj);
  1044. dR2 = _mm_add_pd(szj, dR2);
  1045. __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
  1046. __m128d xhalf = _mm_mul_pd (half, dR2);
  1047. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  1048. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  1049. __m128d S_d = _mm_cvtps_pd(S_s);
  1050. // To handle the condition when src and trg coincide
  1051. S_d = _mm_andnot_pd (reqzero, S_d);
  1052. S = _mm_mul_pd (S_d, S_d);
  1053. S = _mm_mul_pd (S, xhalf);
  1054. S = _mm_sub_pd (opf, S);
  1055. S = _mm_mul_pd (S, S_d);
  1056. S2 = _mm_mul_pd (S, S);
  1057. S3 = _mm_mul_pd (S2, S);
  1058. __m128d S3_sden=_mm_mul_pd(S3, sden);
  1059. __m128d dot_sum = _mm_add_pd(_mm_mul_pd(snormx,dX),_mm_mul_pd(snormy,dY));
  1060. dot_sum = _mm_add_pd(dot_sum,_mm_mul_pd(snormz,dZ));
  1061. temp = _mm_add_pd(_mm_mul_pd(S3_sden,dot_sum),temp);
  1062. }
  1063. temp = _mm_mul_pd (temp, oofp);
  1064. _mm_store_pd(tempval, temp);
  1065. for (int k = 0; k < SIMD_LEN; k++) {
  1066. trgVal[i] += tempval[k];
  1067. }
  1068. for (; j < ns; j++) {
  1069. double x = tx[i] - sx[j];
  1070. double y = ty[i] - sy[j];
  1071. double z = tz[i] - sz[j];
  1072. double r2 = x*x + y*y + z*z;
  1073. double r = sqrt(r2);
  1074. double invdr;
  1075. if (r == 0)
  1076. invdr = 0;
  1077. else
  1078. invdr = 1/r;
  1079. double invdr2=invdr*invdr;
  1080. double invdr3=invdr2*invdr;
  1081. double dot_sum = x*srcDen[j*4+0] + y*srcDen[j*4+1] + z*srcDen[j*4+2];
  1082. trgVal[i] += OOFP*invdr3*x*srcDen[j*4+3]*dot_sum;
  1083. }
  1084. }
  1085. return;
  1086. }
  1087. void laplaceGradSSE(
  1088. const int ns,
  1089. const int nt,
  1090. const double *sx,
  1091. const double *sy,
  1092. const double *sz,
  1093. const double *tx,
  1094. const double *ty,
  1095. const double *tz,
  1096. const double *srcDen,
  1097. double *trgVal)
  1098. {
  1099. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1100. abort();
  1101. double OOFP = 1.0/(4.0*const_pi<double>());
  1102. __m128d tempx; __m128d tempy; __m128d tempz;
  1103. double aux_arr[3*SIMD_LEN+1];
  1104. double *tempvalx, *tempvaly, *tempvalz;
  1105. // if aux_arr is misaligned
  1106. if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempvalx = aux_arr + 1;
  1107. else tempvalx = aux_arr;
  1108. if (size_t(tempvalx)%IDEAL_ALIGNMENT) abort();
  1109. tempvaly=tempvalx+SIMD_LEN;
  1110. tempvalz=tempvaly+SIMD_LEN;
  1111. /*! One over four pi */
  1112. __m128d oofp = _mm_set1_pd (OOFP);
  1113. __m128d half = _mm_set1_pd (0.5);
  1114. __m128d opf = _mm_set1_pd (1.5);
  1115. __m128d zero = _mm_setzero_pd ();
  1116. // loop over sources
  1117. int i = 0;
  1118. for (; i < nt; i++) {
  1119. tempx = _mm_setzero_pd();
  1120. tempy = _mm_setzero_pd();
  1121. tempz = _mm_setzero_pd();
  1122. __m128d txi = _mm_load1_pd (&tx[i]);
  1123. __m128d tyi = _mm_load1_pd (&ty[i]);
  1124. __m128d tzi = _mm_load1_pd (&tz[i]);
  1125. int j = 0;
  1126. // Load and calculate in groups of SIMD_LEN
  1127. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1128. __m128d sxj = _mm_load_pd (&sx[j]);
  1129. __m128d syj = _mm_load_pd (&sy[j]);
  1130. __m128d szj = _mm_load_pd (&sz[j]);
  1131. __m128d sden = _mm_set_pd (srcDen[j+1], srcDen[j]);
  1132. __m128d dX, dY, dZ;
  1133. __m128d dR2;
  1134. __m128d S;
  1135. __m128d S2;
  1136. __m128d S3;
  1137. dX = _mm_sub_pd(txi , sxj);
  1138. dY = _mm_sub_pd(tyi , syj);
  1139. dZ = _mm_sub_pd(tzi , szj);
  1140. sxj = _mm_mul_pd(dX, dX);
  1141. syj = _mm_mul_pd(dY, dY);
  1142. szj = _mm_mul_pd(dZ, dZ);
  1143. dR2 = _mm_add_pd(sxj, syj);
  1144. dR2 = _mm_add_pd(szj, dR2);
  1145. __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
  1146. __m128d xhalf = _mm_mul_pd (half, dR2);
  1147. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  1148. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  1149. __m128d S_d = _mm_cvtps_pd(S_s);
  1150. // To handle the condition when src and trg coincide
  1151. S_d = _mm_andnot_pd (reqzero, S_d);
  1152. S = _mm_mul_pd (S_d, S_d);
  1153. S = _mm_mul_pd (S, xhalf);
  1154. S = _mm_sub_pd (opf, S);
  1155. S = _mm_mul_pd (S, S_d);
  1156. S2 = _mm_mul_pd (S, S);
  1157. S3 = _mm_mul_pd (S2, S);
  1158. __m128d S3_sden=_mm_mul_pd(S3, sden);
  1159. tempx = _mm_add_pd(_mm_mul_pd(S3_sden,dX),tempx);
  1160. tempy = _mm_add_pd(_mm_mul_pd(S3_sden,dY),tempy);
  1161. tempz = _mm_add_pd(_mm_mul_pd(S3_sden,dZ),tempz);
  1162. }
  1163. tempx = _mm_mul_pd (tempx, oofp);
  1164. tempy = _mm_mul_pd (tempy, oofp);
  1165. tempz = _mm_mul_pd (tempz, oofp);
  1166. _mm_store_pd(tempvalx, tempx);
  1167. _mm_store_pd(tempvaly, tempy);
  1168. _mm_store_pd(tempvalz, tempz);
  1169. for (int k = 0; k < SIMD_LEN; k++) {
  1170. trgVal[i*3 ] += tempvalx[k];
  1171. trgVal[i*3+1] += tempvaly[k];
  1172. trgVal[i*3+2] += tempvalz[k];
  1173. }
  1174. for (; j < ns; j++) {
  1175. double x = tx[i] - sx[j];
  1176. double y = ty[i] - sy[j];
  1177. double z = tz[i] - sz[j];
  1178. double r2 = x*x + y*y + z*z;
  1179. double r = sqrt(r2);
  1180. double invdr;
  1181. if (r == 0)
  1182. invdr = 0;
  1183. else
  1184. invdr = 1/r;
  1185. double invdr2=invdr*invdr;
  1186. double invdr3=invdr2*invdr;
  1187. trgVal[i*3 ] += OOFP*invdr3*x*srcDen[j];
  1188. trgVal[i*3+1] += OOFP*invdr3*y*srcDen[j];
  1189. trgVal[i*3+2] += OOFP*invdr3*z*srcDen[j];
  1190. }
  1191. }
  1192. return;
  1193. }
  1194. #undef SIMD_LEN
  1195. #define X(s,k) (s)[(k)*COORD_DIM]
  1196. #define Y(s,k) (s)[(k)*COORD_DIM+1]
  1197. #define Z(s,k) (s)[(k)*COORD_DIM+2]
  1198. void laplaceSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
  1199. {
  1200. // TODO
  1201. }
  1202. void laplaceSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
  1203. {
  1204. double* buff=NULL;
  1205. buff=mem::aligned_new<double>((ns+1+nt)*3,mem_mgr);
  1206. double* buff_=buff;
  1207. pvfmm::Vector<double> xs(ns+1,buff_,false); buff_+=ns+1;
  1208. pvfmm::Vector<double> ys(ns+1,buff_,false); buff_+=ns+1;
  1209. pvfmm::Vector<double> zs(ns+1,buff_,false); buff_+=ns+1;
  1210. pvfmm::Vector<double> xt(nt ,buff_,false); buff_+=nt ;
  1211. pvfmm::Vector<double> yt(nt ,buff_,false); buff_+=nt ;
  1212. pvfmm::Vector<double> zt(nt ,buff_,false); buff_+=nt ;
  1213. //std::vector<double> xs(ns+1);
  1214. //std::vector<double> ys(ns+1);
  1215. //std::vector<double> zs(ns+1);
  1216. //std::vector<double> xt(nt );
  1217. //std::vector<double> yt(nt );
  1218. //std::vector<double> zt(nt );
  1219. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1220. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  1221. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1222. //1. reshuffle memory
  1223. for (int k =0;k<ns;k++){
  1224. xs[k+x_shift]=X(src,k);
  1225. ys[k+y_shift]=Y(src,k);
  1226. zs[k+z_shift]=Z(src,k);
  1227. }
  1228. for (int k=0;k<nt;k++){
  1229. xt[k]=X(trg,k);
  1230. yt[k]=Y(trg,k);
  1231. zt[k]=Z(trg,k);
  1232. }
  1233. //2. perform caclulation
  1234. laplaceSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
  1235. mem::aligned_delete<double>(buff,mem_mgr);
  1236. return;
  1237. }
  1238. void laplaceDblSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
  1239. {
  1240. // TODO
  1241. }
  1242. void laplaceDblSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
  1243. {
  1244. std::vector<double> xs(ns+1); std::vector<double> xt(nt);
  1245. std::vector<double> ys(ns+1); std::vector<double> yt(nt);
  1246. std::vector<double> zs(ns+1); std::vector<double> zt(nt);
  1247. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1248. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  1249. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1250. //1. reshuffle memory
  1251. for (int k =0;k<ns;k++){
  1252. xs[k+x_shift]=X(src,k);
  1253. ys[k+y_shift]=Y(src,k);
  1254. zs[k+z_shift]=Z(src,k);
  1255. }
  1256. for (int k=0;k<nt;k++){
  1257. xt[k]=X(trg,k);
  1258. yt[k]=Y(trg,k);
  1259. zt[k]=Z(trg,k);
  1260. }
  1261. //2. perform caclulation
  1262. laplaceDblSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
  1263. return;
  1264. }
  1265. void laplaceGradSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
  1266. {
  1267. // TODO
  1268. }
  1269. void laplaceGradSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
  1270. {
  1271. int tid=omp_get_thread_num();
  1272. static std::vector<std::vector<double> > xs_(100); static std::vector<std::vector<double> > xt_(100);
  1273. static std::vector<std::vector<double> > ys_(100); static std::vector<std::vector<double> > yt_(100);
  1274. static std::vector<std::vector<double> > zs_(100); static std::vector<std::vector<double> > zt_(100);
  1275. std::vector<double>& xs=xs_[tid]; std::vector<double>& xt=xt_[tid];
  1276. std::vector<double>& ys=ys_[tid]; std::vector<double>& yt=yt_[tid];
  1277. std::vector<double>& zs=zs_[tid]; std::vector<double>& zt=zt_[tid];
  1278. xs.resize(ns+1); xt.resize(nt);
  1279. ys.resize(ns+1); yt.resize(nt);
  1280. zs.resize(ns+1); zt.resize(nt);
  1281. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1282. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  1283. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  1284. //1. reshuffle memory
  1285. for (int k =0;k<ns;k++){
  1286. xs[k+x_shift]=X(src,k);
  1287. ys[k+y_shift]=Y(src,k);
  1288. zs[k+z_shift]=Z(src,k);
  1289. }
  1290. for (int k=0;k<nt;k++){
  1291. xt[k]=X(trg,k);
  1292. yt[k]=Y(trg,k);
  1293. zt[k]=Z(trg,k);
  1294. }
  1295. //2. perform caclulation
  1296. laplaceGradSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
  1297. return;
  1298. }
  1299. #undef X
  1300. #undef Y
  1301. #undef Z
  1302. #undef IDEAL_ALIGNMENT
  1303. #undef DECL_SIMD_ALIGNED
  1304. }
  1305. template <>
  1306. void laplace_poten<double>(double* r_src, int src_cnt, double* v_src, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  1307. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(12*dof));
  1308. if(dof==1){
  1309. laplaceSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
  1310. return;
  1311. }
  1312. }
  1313. template <>
  1314. void laplace_dbl_poten<double>(double* r_src, int src_cnt, double* v_src, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  1315. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(19*dof));
  1316. if(dof==1){
  1317. laplaceDblSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
  1318. return;
  1319. }
  1320. }
  1321. template <>
  1322. void laplace_grad<double>(double* r_src, int src_cnt, double* v_src, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  1323. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+12*dof));
  1324. if(dof==1){
  1325. laplaceGradSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
  1326. return;
  1327. }
  1328. }
  1329. #endif
  1330. #endif
  1331. ////////////////////////////////////////////////////////////////////////////////
  1332. //////// STOKES KERNEL ////////
  1333. ////////////////////////////////////////////////////////////////////////////////
  1334. /**
  1335. * \brief Green's function for the Stokes's equation. Kernel tensor
  1336. * dimension = 3x3.
  1337. */
  1338. template <class T>
  1339. void stokes_vel(T* r_src, int src_cnt, T* v_src_, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  1340. #ifndef __MIC__
  1341. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(28*dof));
  1342. #endif
  1343. const T mu=1.0;
  1344. const T OOEPMU = 1.0/(8.0*const_pi<T>()*mu);
  1345. for(int t=0;t<trg_cnt;t++){
  1346. for(int i=0;i<dof;i++){
  1347. T p[3]={0,0,0};
  1348. for(int s=0;s<src_cnt;s++){
  1349. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  1350. r_trg[3*t+1]-r_src[3*s+1],
  1351. r_trg[3*t+2]-r_src[3*s+2]};
  1352. T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  1353. if (R!=0){
  1354. T invR2=1.0/R;
  1355. T invR=sqrt(invR2);
  1356. T v_src[3]={v_src_[(s*dof+i)*3 ],
  1357. v_src_[(s*dof+i)*3+1],
  1358. v_src_[(s*dof+i)*3+2]};
  1359. T inner_prod=(v_src[0]*dR[0] +
  1360. v_src[1]*dR[1] +
  1361. v_src[2]*dR[2])* invR2;
  1362. p[0] += (v_src[0] + dR[0]*inner_prod)*invR;
  1363. p[1] += (v_src[1] + dR[1]*inner_prod)*invR;
  1364. p[2] += (v_src[2] + dR[2]*inner_prod)*invR;
  1365. }
  1366. }
  1367. k_out[(t*dof+i)*3+0] += p[0]*OOEPMU;
  1368. k_out[(t*dof+i)*3+1] += p[1]*OOEPMU;
  1369. k_out[(t*dof+i)*3+2] += p[2]*OOEPMU;
  1370. }
  1371. }
  1372. }
  1373. template <class T>
  1374. void stokes_sym_dip(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  1375. #ifndef __MIC__
  1376. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(47*dof));
  1377. #endif
  1378. const T mu=1.0;
  1379. const T OOEPMU = -1.0/(8.0*const_pi<T>()*mu);
  1380. for(int t=0;t<trg_cnt;t++){
  1381. for(int i=0;i<dof;i++){
  1382. T p[3]={0,0,0};
  1383. for(int s=0;s<src_cnt;s++){
  1384. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  1385. r_trg[3*t+1]-r_src[3*s+1],
  1386. r_trg[3*t+2]-r_src[3*s+2]};
  1387. T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  1388. if (R!=0){
  1389. T invR2=1.0/R;
  1390. T invR=sqrt(invR2);
  1391. T invR3=invR2*invR;
  1392. T* f=&v_src[(s*dof+i)*6+0];
  1393. T* n=&v_src[(s*dof+i)*6+3];
  1394. T r_dot_n=(n[0]*dR[0]+n[1]*dR[1]+n[2]*dR[2]);
  1395. T r_dot_f=(f[0]*dR[0]+f[1]*dR[1]+f[2]*dR[2]);
  1396. T n_dot_f=(f[0]* n[0]+f[1]* n[1]+f[2]* n[2]);
  1397. p[0] += dR[0]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
  1398. p[1] += dR[1]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
  1399. p[2] += dR[2]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
  1400. }
  1401. }
  1402. k_out[(t*dof+i)*3+0] += p[0]*OOEPMU;
  1403. k_out[(t*dof+i)*3+1] += p[1]*OOEPMU;
  1404. k_out[(t*dof+i)*3+2] += p[2]*OOEPMU;
  1405. }
  1406. }
  1407. }
  1408. template <class T>
  1409. void stokes_press(T* r_src, int src_cnt, T* v_src_, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  1410. #ifndef __MIC__
  1411. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(17*dof));
  1412. #endif
  1413. const T OOFP = 1.0/(4.0*const_pi<T>());
  1414. for(int t=0;t<trg_cnt;t++){
  1415. for(int i=0;i<dof;i++){
  1416. T p=0;
  1417. for(int s=0;s<src_cnt;s++){
  1418. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  1419. r_trg[3*t+1]-r_src[3*s+1],
  1420. r_trg[3*t+2]-r_src[3*s+2]};
  1421. T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  1422. if (R!=0){
  1423. T invR2=1.0/R;
  1424. T invR=sqrt(invR2);
  1425. T invR3=invR2*invR;
  1426. T v_src[3]={v_src_[(s*dof+i)*3 ],
  1427. v_src_[(s*dof+i)*3+1],
  1428. v_src_[(s*dof+i)*3+2]};
  1429. T inner_prod=(v_src[0]*dR[0] +
  1430. v_src[1]*dR[1] +
  1431. v_src[2]*dR[2])* invR3;
  1432. p += inner_prod;
  1433. }
  1434. }
  1435. k_out[t*dof+i] += p*OOFP;
  1436. }
  1437. }
  1438. }
  1439. template <class T>
  1440. void stokes_stress(T* r_src, int src_cnt, T* v_src_, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  1441. #ifndef __MIC__
  1442. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(45*dof));
  1443. #endif
  1444. const T TOFP = -3.0/(4.0*const_pi<T>());
  1445. for(int t=0;t<trg_cnt;t++){
  1446. for(int i=0;i<dof;i++){
  1447. T p[9]={0,0,0,
  1448. 0,0,0,
  1449. 0,0,0};
  1450. for(int s=0;s<src_cnt;s++){
  1451. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  1452. r_trg[3*t+1]-r_src[3*s+1],
  1453. r_trg[3*t+2]-r_src[3*s+2]};
  1454. T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  1455. if (R!=0){
  1456. T invR2=1.0/R;
  1457. T invR=sqrt(invR2);
  1458. T invR3=invR2*invR;
  1459. T invR5=invR3*invR2;
  1460. T v_src[3]={v_src_[(s*dof+i)*3 ],
  1461. v_src_[(s*dof+i)*3+1],
  1462. v_src_[(s*dof+i)*3+2]};
  1463. T inner_prod=(v_src[0]*dR[0] +
  1464. v_src[1]*dR[1] +
  1465. v_src[2]*dR[2])* invR5;
  1466. p[0] += inner_prod*dR[0]*dR[0]; p[1] += inner_prod*dR[1]*dR[0]; p[2] += inner_prod*dR[2]*dR[0];
  1467. p[3] += inner_prod*dR[0]*dR[1]; p[4] += inner_prod*dR[1]*dR[1]; p[5] += inner_prod*dR[2]*dR[1];
  1468. p[6] += inner_prod*dR[0]*dR[2]; p[7] += inner_prod*dR[1]*dR[2]; p[8] += inner_prod*dR[2]*dR[2];
  1469. }
  1470. }
  1471. k_out[(t*dof+i)*9+0] += p[0]*TOFP;
  1472. k_out[(t*dof+i)*9+1] += p[1]*TOFP;
  1473. k_out[(t*dof+i)*9+2] += p[2]*TOFP;
  1474. k_out[(t*dof+i)*9+3] += p[3]*TOFP;
  1475. k_out[(t*dof+i)*9+4] += p[4]*TOFP;
  1476. k_out[(t*dof+i)*9+5] += p[5]*TOFP;
  1477. k_out[(t*dof+i)*9+6] += p[6]*TOFP;
  1478. k_out[(t*dof+i)*9+7] += p[7]*TOFP;
  1479. k_out[(t*dof+i)*9+8] += p[8]*TOFP;
  1480. }
  1481. }
  1482. }
  1483. template <class T>
  1484. void stokes_grad(T* r_src, int src_cnt, T* v_src_, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  1485. #ifndef __MIC__
  1486. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(89*dof));
  1487. #endif
  1488. const T mu=1.0;
  1489. const T OOEPMU = 1.0/(8.0*const_pi<T>()*mu);
  1490. for(int t=0;t<trg_cnt;t++){
  1491. for(int i=0;i<dof;i++){
  1492. T p[9]={0,0,0,
  1493. 0,0,0,
  1494. 0,0,0};
  1495. for(int s=0;s<src_cnt;s++){
  1496. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  1497. r_trg[3*t+1]-r_src[3*s+1],
  1498. r_trg[3*t+2]-r_src[3*s+2]};
  1499. T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  1500. if (R!=0){
  1501. T invR2=1.0/R;
  1502. T invR=sqrt(invR2);
  1503. T invR3=invR2*invR;
  1504. T v_src[3]={v_src_[(s*dof+i)*3 ],
  1505. v_src_[(s*dof+i)*3+1],
  1506. v_src_[(s*dof+i)*3+2]};
  1507. T inner_prod=(v_src[0]*dR[0] +
  1508. v_src[1]*dR[1] +
  1509. v_src[2]*dR[2]);
  1510. p[0] += ( inner_prod*(1-3*dR[0]*dR[0]*invR2))*invR3; //6
  1511. p[1] += (dR[1]*v_src[0]-v_src[1]*dR[0]+inner_prod*( -3*dR[1]*dR[0]*invR2))*invR3; //9
  1512. p[2] += (dR[2]*v_src[0]-v_src[2]*dR[0]+inner_prod*( -3*dR[2]*dR[0]*invR2))*invR3;
  1513. p[3] += (dR[0]*v_src[1]-v_src[0]*dR[1]+inner_prod*( -3*dR[0]*dR[1]*invR2))*invR3;
  1514. p[4] += ( inner_prod*(1-3*dR[1]*dR[1]*invR2))*invR3;
  1515. p[5] += (dR[2]*v_src[1]-v_src[2]*dR[1]+inner_prod*( -3*dR[2]*dR[1]*invR2))*invR3;
  1516. p[6] += (dR[0]*v_src[2]-v_src[0]*dR[2]+inner_prod*( -3*dR[0]*dR[2]*invR2))*invR3;
  1517. p[7] += (dR[1]*v_src[2]-v_src[1]*dR[2]+inner_prod*( -3*dR[1]*dR[2]*invR2))*invR3;
  1518. p[8] += ( inner_prod*(1-3*dR[2]*dR[2]*invR2))*invR3;
  1519. }
  1520. }
  1521. k_out[(t*dof+i)*9+0] += p[0]*OOEPMU;
  1522. k_out[(t*dof+i)*9+1] += p[1]*OOEPMU;
  1523. k_out[(t*dof+i)*9+2] += p[2]*OOEPMU;
  1524. k_out[(t*dof+i)*9+3] += p[3]*OOEPMU;
  1525. k_out[(t*dof+i)*9+4] += p[4]*OOEPMU;
  1526. k_out[(t*dof+i)*9+5] += p[5]*OOEPMU;
  1527. k_out[(t*dof+i)*9+6] += p[6]*OOEPMU;
  1528. k_out[(t*dof+i)*9+7] += p[7]*OOEPMU;
  1529. k_out[(t*dof+i)*9+8] += p[8]*OOEPMU;
  1530. }
  1531. }
  1532. }
  1533. #ifndef __MIC__
  1534. #ifdef USE_SSE
  1535. namespace
  1536. {
  1537. #define IDEAL_ALIGNMENT 16
  1538. #define SIMD_LEN (int)(IDEAL_ALIGNMENT / sizeof(double))
  1539. #define DECL_SIMD_ALIGNED __declspec(align(IDEAL_ALIGNMENT))
  1540. void stokesDirectVecSSE(
  1541. const int ns,
  1542. const int nt,
  1543. const double *sx,
  1544. const double *sy,
  1545. const double *sz,
  1546. const double *tx,
  1547. const double *ty,
  1548. const double *tz,
  1549. const double *srcDen,
  1550. double *trgVal,
  1551. const double cof )
  1552. {
  1553. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1554. abort();
  1555. double mu = cof;
  1556. double OOEP = 1.0/(8.0*const_pi<double>());
  1557. __m128d tempx;
  1558. __m128d tempy;
  1559. __m128d tempz;
  1560. double oomeu = 1/mu;
  1561. double aux_arr[3*SIMD_LEN+1];
  1562. double *tempvalx;
  1563. double *tempvaly;
  1564. double *tempvalz;
  1565. if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
  1566. {
  1567. tempvalx = aux_arr + 1;
  1568. if (size_t(tempvalx)%IDEAL_ALIGNMENT)
  1569. abort();
  1570. }
  1571. else
  1572. tempvalx = aux_arr;
  1573. tempvaly=tempvalx+SIMD_LEN;
  1574. tempvalz=tempvaly+SIMD_LEN;
  1575. /*! One over eight pi */
  1576. __m128d ooep = _mm_set1_pd (OOEP);
  1577. __m128d half = _mm_set1_pd (0.5);
  1578. __m128d opf = _mm_set1_pd (1.5);
  1579. __m128d zero = _mm_setzero_pd ();
  1580. __m128d oomu = _mm_set1_pd (1/mu);
  1581. // loop over sources
  1582. int i = 0;
  1583. for (; i < nt; i++) {
  1584. tempx = _mm_setzero_pd();
  1585. tempy = _mm_setzero_pd();
  1586. tempz = _mm_setzero_pd();
  1587. __m128d txi = _mm_load1_pd (&tx[i]);
  1588. __m128d tyi = _mm_load1_pd (&ty[i]);
  1589. __m128d tzi = _mm_load1_pd (&tz[i]);
  1590. int j = 0;
  1591. // Load and calculate in groups of SIMD_LEN
  1592. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1593. __m128d sxj = _mm_load_pd (&sx[j]);
  1594. __m128d syj = _mm_load_pd (&sy[j]);
  1595. __m128d szj = _mm_load_pd (&sz[j]);
  1596. __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
  1597. __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
  1598. __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
  1599. __m128d dX, dY, dZ;
  1600. __m128d dR2;
  1601. __m128d S;
  1602. dX = _mm_sub_pd(txi , sxj);
  1603. dY = _mm_sub_pd(tyi , syj);
  1604. dZ = _mm_sub_pd(tzi , szj);
  1605. sxj = _mm_mul_pd(dX, dX);
  1606. syj = _mm_mul_pd(dY, dY);
  1607. szj = _mm_mul_pd(dZ, dZ);
  1608. dR2 = _mm_add_pd(sxj, syj);
  1609. dR2 = _mm_add_pd(szj, dR2);
  1610. __m128d temp = _mm_cmpeq_pd (dR2, zero);
  1611. __m128d xhalf = _mm_mul_pd (half, dR2);
  1612. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  1613. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  1614. __m128d S_d = _mm_cvtps_pd(S_s);
  1615. // To handle the condition when src and trg coincide
  1616. S_d = _mm_andnot_pd (temp, S_d);
  1617. S = _mm_mul_pd (S_d, S_d);
  1618. S = _mm_mul_pd (S, xhalf);
  1619. S = _mm_sub_pd (opf, S);
  1620. S = _mm_mul_pd (S, S_d);
  1621. __m128d dotx = _mm_mul_pd (dX, sdenx);
  1622. __m128d doty = _mm_mul_pd (dY, sdeny);
  1623. __m128d dotz = _mm_mul_pd (dZ, sdenz);
  1624. __m128d dot_sum = _mm_add_pd (dotx, doty);
  1625. dot_sum = _mm_add_pd (dot_sum, dotz);
  1626. dot_sum = _mm_mul_pd (dot_sum, S);
  1627. dot_sum = _mm_mul_pd (dot_sum, S);
  1628. dotx = _mm_mul_pd (dot_sum, dX);
  1629. doty = _mm_mul_pd (dot_sum, dY);
  1630. dotz = _mm_mul_pd (dot_sum, dZ);
  1631. sdenx = _mm_add_pd (sdenx, dotx);
  1632. sdeny = _mm_add_pd (sdeny, doty);
  1633. sdenz = _mm_add_pd (sdenz, dotz);
  1634. sdenx = _mm_mul_pd (sdenx, S);
  1635. sdeny = _mm_mul_pd (sdeny, S);
  1636. sdenz = _mm_mul_pd (sdenz, S);
  1637. tempx = _mm_add_pd (sdenx, tempx);
  1638. tempy = _mm_add_pd (sdeny, tempy);
  1639. tempz = _mm_add_pd (sdenz, tempz);
  1640. }
  1641. tempx = _mm_mul_pd (tempx, ooep);
  1642. tempy = _mm_mul_pd (tempy, ooep);
  1643. tempz = _mm_mul_pd (tempz, ooep);
  1644. tempx = _mm_mul_pd (tempx, oomu);
  1645. tempy = _mm_mul_pd (tempy, oomu);
  1646. tempz = _mm_mul_pd (tempz, oomu);
  1647. _mm_store_pd(tempvalx, tempx);
  1648. _mm_store_pd(tempvaly, tempy);
  1649. _mm_store_pd(tempvalz, tempz);
  1650. for (int k = 0; k < SIMD_LEN; k++) {
  1651. trgVal[i*3] += tempvalx[k];
  1652. trgVal[i*3+1] += tempvaly[k];
  1653. trgVal[i*3+2] += tempvalz[k];
  1654. }
  1655. for (; j < ns; j++) {
  1656. double x = tx[i] - sx[j];
  1657. double y = ty[i] - sy[j];
  1658. double z = tz[i] - sz[j];
  1659. double r2 = x*x + y*y + z*z;
  1660. double r = sqrt(r2);
  1661. double invdr;
  1662. if (r == 0)
  1663. invdr = 0;
  1664. else
  1665. invdr = 1/r;
  1666. double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr * invdr;
  1667. double denx = srcDen[j*3] + dot*x;
  1668. double deny = srcDen[j*3+1] + dot*y;
  1669. double denz = srcDen[j*3+2] + dot*z;
  1670. trgVal[i*3] += denx*invdr*OOEP*oomeu;
  1671. trgVal[i*3+1] += deny*invdr*OOEP*oomeu;
  1672. trgVal[i*3+2] += denz*invdr*OOEP*oomeu;
  1673. }
  1674. }
  1675. return;
  1676. }
  1677. void stokesPressureSSE(
  1678. const int ns,
  1679. const int nt,
  1680. const double *sx,
  1681. const double *sy,
  1682. const double *sz,
  1683. const double *tx,
  1684. const double *ty,
  1685. const double *tz,
  1686. const double *srcDen,
  1687. double *trgVal)
  1688. {
  1689. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1690. abort();
  1691. double OOFP = 1.0/(4.0*const_pi<double>());
  1692. __m128d temp_press;
  1693. double aux_arr[SIMD_LEN+1];
  1694. double *tempval_press;
  1695. if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
  1696. {
  1697. tempval_press = aux_arr + 1;
  1698. if (size_t(tempval_press)%IDEAL_ALIGNMENT)
  1699. abort();
  1700. }
  1701. else
  1702. tempval_press = aux_arr;
  1703. /*! One over eight pi */
  1704. __m128d oofp = _mm_set1_pd (OOFP);
  1705. __m128d half = _mm_set1_pd (0.5);
  1706. __m128d opf = _mm_set1_pd (1.5);
  1707. __m128d zero = _mm_setzero_pd ();
  1708. // loop over sources
  1709. int i = 0;
  1710. for (; i < nt; i++) {
  1711. temp_press = _mm_setzero_pd();
  1712. __m128d txi = _mm_load1_pd (&tx[i]);
  1713. __m128d tyi = _mm_load1_pd (&ty[i]);
  1714. __m128d tzi = _mm_load1_pd (&tz[i]);
  1715. int j = 0;
  1716. // Load and calculate in groups of SIMD_LEN
  1717. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1718. __m128d sxj = _mm_load_pd (&sx[j]);
  1719. __m128d syj = _mm_load_pd (&sy[j]);
  1720. __m128d szj = _mm_load_pd (&sz[j]);
  1721. __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
  1722. __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
  1723. __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
  1724. __m128d dX, dY, dZ;
  1725. __m128d dR2;
  1726. __m128d S;
  1727. dX = _mm_sub_pd(txi , sxj);
  1728. dY = _mm_sub_pd(tyi , syj);
  1729. dZ = _mm_sub_pd(tzi , szj);
  1730. sxj = _mm_mul_pd(dX, dX);
  1731. syj = _mm_mul_pd(dY, dY);
  1732. szj = _mm_mul_pd(dZ, dZ);
  1733. dR2 = _mm_add_pd(sxj, syj);
  1734. dR2 = _mm_add_pd(szj, dR2);
  1735. __m128d temp = _mm_cmpeq_pd (dR2, zero);
  1736. __m128d xhalf = _mm_mul_pd (half, dR2);
  1737. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  1738. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  1739. __m128d S_d = _mm_cvtps_pd(S_s);
  1740. // To handle the condition when src and trg coincide
  1741. S_d = _mm_andnot_pd (temp, S_d);
  1742. S = _mm_mul_pd (S_d, S_d);
  1743. S = _mm_mul_pd (S, xhalf);
  1744. S = _mm_sub_pd (opf, S);
  1745. S = _mm_mul_pd (S, S_d);
  1746. __m128d dotx = _mm_mul_pd (dX, sdenx);
  1747. __m128d doty = _mm_mul_pd (dY, sdeny);
  1748. __m128d dotz = _mm_mul_pd (dZ, sdenz);
  1749. __m128d dot_sum = _mm_add_pd (dotx, doty);
  1750. dot_sum = _mm_add_pd (dot_sum, dotz);
  1751. dot_sum = _mm_mul_pd (dot_sum, S);
  1752. dot_sum = _mm_mul_pd (dot_sum, S);
  1753. dot_sum = _mm_mul_pd (dot_sum, S);
  1754. temp_press = _mm_add_pd (dot_sum, temp_press);
  1755. }
  1756. temp_press = _mm_mul_pd (temp_press, oofp);
  1757. _mm_store_pd(tempval_press, temp_press);
  1758. for (int k = 0; k < SIMD_LEN; k++) {
  1759. trgVal[i] += tempval_press[k];
  1760. }
  1761. for (; j < ns; j++) {
  1762. double x = tx[i] - sx[j];
  1763. double y = ty[i] - sy[j];
  1764. double z = tz[i] - sz[j];
  1765. double r2 = x*x + y*y + z*z;
  1766. double r = sqrt(r2);
  1767. double invdr;
  1768. if (r == 0)
  1769. invdr = 0;
  1770. else
  1771. invdr = 1/r;
  1772. double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr * invdr * invdr;
  1773. trgVal[i] += dot*OOFP;
  1774. }
  1775. }
  1776. return;
  1777. }
  1778. void stokesStressSSE(
  1779. const int ns,
  1780. const int nt,
  1781. const double *sx,
  1782. const double *sy,
  1783. const double *sz,
  1784. const double *tx,
  1785. const double *ty,
  1786. const double *tz,
  1787. const double *srcDen,
  1788. double *trgVal)
  1789. {
  1790. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1791. abort();
  1792. double TOFP = -3.0/(4.0*const_pi<double>());
  1793. __m128d tempxx; __m128d tempxy; __m128d tempxz;
  1794. __m128d tempyx; __m128d tempyy; __m128d tempyz;
  1795. __m128d tempzx; __m128d tempzy; __m128d tempzz;
  1796. double aux_arr[9*SIMD_LEN+1];
  1797. double *tempvalxx, *tempvalxy, *tempvalxz;
  1798. double *tempvalyx, *tempvalyy, *tempvalyz;
  1799. double *tempvalzx, *tempvalzy, *tempvalzz;
  1800. if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
  1801. {
  1802. tempvalxx = aux_arr + 1;
  1803. if (size_t(tempvalxx)%IDEAL_ALIGNMENT)
  1804. abort();
  1805. }
  1806. else
  1807. tempvalxx = aux_arr;
  1808. tempvalxy=tempvalxx+SIMD_LEN;
  1809. tempvalxz=tempvalxy+SIMD_LEN;
  1810. tempvalyx=tempvalxz+SIMD_LEN;
  1811. tempvalyy=tempvalyx+SIMD_LEN;
  1812. tempvalyz=tempvalyy+SIMD_LEN;
  1813. tempvalzx=tempvalyz+SIMD_LEN;
  1814. tempvalzy=tempvalzx+SIMD_LEN;
  1815. tempvalzz=tempvalzy+SIMD_LEN;
  1816. /*! One over eight pi */
  1817. __m128d tofp = _mm_set1_pd (TOFP);
  1818. __m128d half = _mm_set1_pd (0.5);
  1819. __m128d opf = _mm_set1_pd (1.5);
  1820. __m128d zero = _mm_setzero_pd ();
  1821. // loop over sources
  1822. int i = 0;
  1823. for (; i < nt; i++) {
  1824. tempxx = _mm_setzero_pd(); tempxy = _mm_setzero_pd(); tempxz = _mm_setzero_pd();
  1825. tempyx = _mm_setzero_pd(); tempyy = _mm_setzero_pd(); tempyz = _mm_setzero_pd();
  1826. tempzx = _mm_setzero_pd(); tempzy = _mm_setzero_pd(); tempzz = _mm_setzero_pd();
  1827. __m128d txi = _mm_load1_pd (&tx[i]);
  1828. __m128d tyi = _mm_load1_pd (&ty[i]);
  1829. __m128d tzi = _mm_load1_pd (&tz[i]);
  1830. int j = 0;
  1831. // Load and calculate in groups of SIMD_LEN
  1832. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1833. __m128d sxj = _mm_load_pd (&sx[j]);
  1834. __m128d syj = _mm_load_pd (&sy[j]);
  1835. __m128d szj = _mm_load_pd (&sz[j]);
  1836. __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
  1837. __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
  1838. __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
  1839. __m128d dX, dY, dZ;
  1840. __m128d dR2;
  1841. __m128d S;
  1842. __m128d S2;
  1843. dX = _mm_sub_pd(txi , sxj);
  1844. dY = _mm_sub_pd(tyi , syj);
  1845. dZ = _mm_sub_pd(tzi , szj);
  1846. sxj = _mm_mul_pd(dX, dX);
  1847. syj = _mm_mul_pd(dY, dY);
  1848. szj = _mm_mul_pd(dZ, dZ);
  1849. dR2 = _mm_add_pd(sxj, syj);
  1850. dR2 = _mm_add_pd(szj, dR2);
  1851. __m128d temp = _mm_cmpeq_pd (dR2, zero);
  1852. __m128d xhalf = _mm_mul_pd (half, dR2);
  1853. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  1854. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  1855. __m128d S_d = _mm_cvtps_pd(S_s);
  1856. // To handle the condition when src and trg coincide
  1857. S_d = _mm_andnot_pd (temp, S_d);
  1858. S = _mm_mul_pd (S_d, S_d);
  1859. S = _mm_mul_pd (S, xhalf);
  1860. S = _mm_sub_pd (opf, S);
  1861. S = _mm_mul_pd (S, S_d);
  1862. S2 = _mm_mul_pd (S, S);
  1863. __m128d dotx = _mm_mul_pd (dX, sdenx);
  1864. __m128d doty = _mm_mul_pd (dY, sdeny);
  1865. __m128d dotz = _mm_mul_pd (dZ, sdenz);
  1866. __m128d dot_sum = _mm_add_pd (dotx, doty);
  1867. dot_sum = _mm_add_pd (dot_sum, dotz);
  1868. dot_sum = _mm_mul_pd (dot_sum, S);
  1869. dot_sum = _mm_mul_pd (dot_sum, S2);
  1870. dot_sum = _mm_mul_pd (dot_sum, S2);
  1871. dotx = _mm_mul_pd (dot_sum, dX);
  1872. doty = _mm_mul_pd (dot_sum, dY);
  1873. dotz = _mm_mul_pd (dot_sum, dZ);
  1874. tempxx = _mm_add_pd (_mm_mul_pd(dotx,dX), tempxx);
  1875. tempxy = _mm_add_pd (_mm_mul_pd(dotx,dY), tempxy);
  1876. tempxz = _mm_add_pd (_mm_mul_pd(dotx,dZ), tempxz);
  1877. tempyx = _mm_add_pd (_mm_mul_pd(doty,dX), tempyx);
  1878. tempyy = _mm_add_pd (_mm_mul_pd(doty,dY), tempyy);
  1879. tempyz = _mm_add_pd (_mm_mul_pd(doty,dZ), tempyz);
  1880. tempzx = _mm_add_pd (_mm_mul_pd(dotz,dX), tempzx);
  1881. tempzy = _mm_add_pd (_mm_mul_pd(dotz,dY), tempzy);
  1882. tempzz = _mm_add_pd (_mm_mul_pd(dotz,dZ), tempzz);
  1883. }
  1884. tempxx = _mm_mul_pd (tempxx, tofp);
  1885. tempxy = _mm_mul_pd (tempxy, tofp);
  1886. tempxz = _mm_mul_pd (tempxz, tofp);
  1887. tempyx = _mm_mul_pd (tempyx, tofp);
  1888. tempyy = _mm_mul_pd (tempyy, tofp);
  1889. tempyz = _mm_mul_pd (tempyz, tofp);
  1890. tempzx = _mm_mul_pd (tempzx, tofp);
  1891. tempzy = _mm_mul_pd (tempzy, tofp);
  1892. tempzz = _mm_mul_pd (tempzz, tofp);
  1893. _mm_store_pd(tempvalxx, tempxx); _mm_store_pd(tempvalxy, tempxy); _mm_store_pd(tempvalxz, tempxz);
  1894. _mm_store_pd(tempvalyx, tempyx); _mm_store_pd(tempvalyy, tempyy); _mm_store_pd(tempvalyz, tempyz);
  1895. _mm_store_pd(tempvalzx, tempzx); _mm_store_pd(tempvalzy, tempzy); _mm_store_pd(tempvalzz, tempzz);
  1896. for (int k = 0; k < SIMD_LEN; k++) {
  1897. trgVal[i*9 ] += tempvalxx[k];
  1898. trgVal[i*9+1] += tempvalxy[k];
  1899. trgVal[i*9+2] += tempvalxz[k];
  1900. trgVal[i*9+3] += tempvalyx[k];
  1901. trgVal[i*9+4] += tempvalyy[k];
  1902. trgVal[i*9+5] += tempvalyz[k];
  1903. trgVal[i*9+6] += tempvalzx[k];
  1904. trgVal[i*9+7] += tempvalzy[k];
  1905. trgVal[i*9+8] += tempvalzz[k];
  1906. }
  1907. for (; j < ns; j++) {
  1908. double x = tx[i] - sx[j];
  1909. double y = ty[i] - sy[j];
  1910. double z = tz[i] - sz[j];
  1911. double r2 = x*x + y*y + z*z;
  1912. double r = sqrt(r2);
  1913. double invdr;
  1914. if (r == 0)
  1915. invdr = 0;
  1916. else
  1917. invdr = 1/r;
  1918. double invdr2=invdr*invdr;
  1919. double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr2 * invdr2 * invdr;
  1920. double denx = dot*x;
  1921. double deny = dot*y;
  1922. double denz = dot*z;
  1923. trgVal[i*9 ] += denx*x*TOFP;
  1924. trgVal[i*9+1] += denx*y*TOFP;
  1925. trgVal[i*9+2] += denx*z*TOFP;
  1926. trgVal[i*9+3] += deny*x*TOFP;
  1927. trgVal[i*9+4] += deny*y*TOFP;
  1928. trgVal[i*9+5] += deny*z*TOFP;
  1929. trgVal[i*9+6] += denz*x*TOFP;
  1930. trgVal[i*9+7] += denz*y*TOFP;
  1931. trgVal[i*9+8] += denz*z*TOFP;
  1932. }
  1933. }
  1934. return;
  1935. }
  1936. void stokesGradSSE(
  1937. const int ns,
  1938. const int nt,
  1939. const double *sx,
  1940. const double *sy,
  1941. const double *sz,
  1942. const double *tx,
  1943. const double *ty,
  1944. const double *tz,
  1945. const double *srcDen,
  1946. double *trgVal,
  1947. const double cof )
  1948. {
  1949. if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
  1950. abort();
  1951. double mu = cof;
  1952. double OOEP = 1.0/(8.0*const_pi<double>());
  1953. __m128d tempxx; __m128d tempxy; __m128d tempxz;
  1954. __m128d tempyx; __m128d tempyy; __m128d tempyz;
  1955. __m128d tempzx; __m128d tempzy; __m128d tempzz;
  1956. double oomeu = 1/mu;
  1957. double aux_arr[9*SIMD_LEN+1];
  1958. double *tempvalxx, *tempvalxy, *tempvalxz;
  1959. double *tempvalyx, *tempvalyy, *tempvalyz;
  1960. double *tempvalzx, *tempvalzy, *tempvalzz;
  1961. if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
  1962. {
  1963. tempvalxx = aux_arr + 1;
  1964. if (size_t(tempvalxx)%IDEAL_ALIGNMENT)
  1965. abort();
  1966. }
  1967. else
  1968. tempvalxx = aux_arr;
  1969. tempvalxy=tempvalxx+SIMD_LEN;
  1970. tempvalxz=tempvalxy+SIMD_LEN;
  1971. tempvalyx=tempvalxz+SIMD_LEN;
  1972. tempvalyy=tempvalyx+SIMD_LEN;
  1973. tempvalyz=tempvalyy+SIMD_LEN;
  1974. tempvalzx=tempvalyz+SIMD_LEN;
  1975. tempvalzy=tempvalzx+SIMD_LEN;
  1976. tempvalzz=tempvalzy+SIMD_LEN;
  1977. /*! One over eight pi */
  1978. __m128d ooep = _mm_set1_pd (OOEP);
  1979. __m128d half = _mm_set1_pd (0.5);
  1980. __m128d opf = _mm_set1_pd (1.5);
  1981. __m128d three = _mm_set1_pd (3.0);
  1982. __m128d zero = _mm_setzero_pd ();
  1983. __m128d oomu = _mm_set1_pd (1/mu);
  1984. __m128d ooepmu = _mm_mul_pd(ooep,oomu);
  1985. // loop over sources
  1986. int i = 0;
  1987. for (; i < nt; i++) {
  1988. tempxx = _mm_setzero_pd(); tempxy = _mm_setzero_pd(); tempxz = _mm_setzero_pd();
  1989. tempyx = _mm_setzero_pd(); tempyy = _mm_setzero_pd(); tempyz = _mm_setzero_pd();
  1990. tempzx = _mm_setzero_pd(); tempzy = _mm_setzero_pd(); tempzz = _mm_setzero_pd();
  1991. __m128d txi = _mm_load1_pd (&tx[i]);
  1992. __m128d tyi = _mm_load1_pd (&ty[i]);
  1993. __m128d tzi = _mm_load1_pd (&tz[i]);
  1994. int j = 0;
  1995. // Load and calculate in groups of SIMD_LEN
  1996. for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
  1997. __m128d sxj = _mm_load_pd (&sx[j]);
  1998. __m128d syj = _mm_load_pd (&sy[j]);
  1999. __m128d szj = _mm_load_pd (&sz[j]);
  2000. __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
  2001. __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
  2002. __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
  2003. __m128d dX, dY, dZ;
  2004. __m128d dR2;
  2005. __m128d S;
  2006. __m128d S2;
  2007. __m128d S3;
  2008. dX = _mm_sub_pd(txi , sxj);
  2009. dY = _mm_sub_pd(tyi , syj);
  2010. dZ = _mm_sub_pd(tzi , szj);
  2011. sxj = _mm_mul_pd(dX, dX);
  2012. syj = _mm_mul_pd(dY, dY);
  2013. szj = _mm_mul_pd(dZ, dZ);
  2014. dR2 = _mm_add_pd(sxj, syj);
  2015. dR2 = _mm_add_pd(szj, dR2);
  2016. __m128d temp = _mm_cmpeq_pd (dR2, zero);
  2017. __m128d xhalf = _mm_mul_pd (half, dR2);
  2018. __m128 dR2_s = _mm_cvtpd_ps(dR2);
  2019. __m128 S_s = _mm_rsqrt_ps(dR2_s);
  2020. __m128d S_d = _mm_cvtps_pd(S_s);
  2021. // To handle the condition when src and trg coincide
  2022. S_d = _mm_andnot_pd (temp, S_d);
  2023. S = _mm_mul_pd (S_d, S_d);
  2024. S = _mm_mul_pd (S, xhalf);
  2025. S = _mm_sub_pd (opf, S);
  2026. S = _mm_mul_pd (S, S_d);
  2027. S2 = _mm_mul_pd (S, S);
  2028. S3 = _mm_mul_pd (S2, S);
  2029. __m128d dotx = _mm_mul_pd (dX, sdenx);
  2030. __m128d doty = _mm_mul_pd (dY, sdeny);
  2031. __m128d dotz = _mm_mul_pd (dZ, sdenz);
  2032. __m128d dot_sum = _mm_add_pd (dotx, doty);
  2033. dot_sum = _mm_add_pd (dot_sum, dotz);
  2034. dot_sum = _mm_mul_pd (dot_sum, S2);
  2035. tempxx = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dX, sdenx), _mm_mul_pd(sdenx, dX)), _mm_mul_pd(dot_sum, _mm_sub_pd(dR2 , _mm_mul_pd(three, _mm_mul_pd(dX, dX)))))),tempxx);
  2036. tempxy = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dY, sdenx), _mm_mul_pd(sdeny, dX)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dY, dX)))))),tempxy);
  2037. tempxz = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dZ, sdenx), _mm_mul_pd(sdenz, dX)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dZ, dX)))))),tempxz);
  2038. tempyx = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dX, sdeny), _mm_mul_pd(sdenx, dY)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dX, dY)))))),tempyx);
  2039. tempyy = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dY, sdeny), _mm_mul_pd(sdeny, dY)), _mm_mul_pd(dot_sum, _mm_sub_pd(dR2 , _mm_mul_pd(three, _mm_mul_pd(dY, dY)))))),tempyy);
  2040. tempyz = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dZ, sdeny), _mm_mul_pd(sdenz, dY)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dZ, dY)))))),tempyz);
  2041. tempzx = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dX, sdenz), _mm_mul_pd(sdenx, dZ)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dX, dZ)))))),tempzx);
  2042. tempzy = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dY, sdenz), _mm_mul_pd(sdeny, dZ)), _mm_mul_pd(dot_sum, _mm_sub_pd(zero, _mm_mul_pd(three, _mm_mul_pd(dY, dZ)))))),tempzy);
  2043. tempzz = _mm_add_pd(_mm_mul_pd(S3,_mm_add_pd(_mm_sub_pd(_mm_mul_pd(dZ, sdenz), _mm_mul_pd(sdenz, dZ)), _mm_mul_pd(dot_sum, _mm_sub_pd(dR2 , _mm_mul_pd(three, _mm_mul_pd(dZ, dZ)))))),tempzz);
  2044. }
  2045. tempxx = _mm_mul_pd (tempxx, ooepmu);
  2046. tempxy = _mm_mul_pd (tempxy, ooepmu);
  2047. tempxz = _mm_mul_pd (tempxz, ooepmu);
  2048. tempyx = _mm_mul_pd (tempyx, ooepmu);
  2049. tempyy = _mm_mul_pd (tempyy, ooepmu);
  2050. tempyz = _mm_mul_pd (tempyz, ooepmu);
  2051. tempzx = _mm_mul_pd (tempzx, ooepmu);
  2052. tempzy = _mm_mul_pd (tempzy, ooepmu);
  2053. tempzz = _mm_mul_pd (tempzz, ooepmu);
  2054. _mm_store_pd(tempvalxx, tempxx); _mm_store_pd(tempvalxy, tempxy); _mm_store_pd(tempvalxz, tempxz);
  2055. _mm_store_pd(tempvalyx, tempyx); _mm_store_pd(tempvalyy, tempyy); _mm_store_pd(tempvalyz, tempyz);
  2056. _mm_store_pd(tempvalzx, tempzx); _mm_store_pd(tempvalzy, tempzy); _mm_store_pd(tempvalzz, tempzz);
  2057. for (int k = 0; k < SIMD_LEN; k++) {
  2058. trgVal[i*9 ] += tempvalxx[k];
  2059. trgVal[i*9+1] += tempvalxy[k];
  2060. trgVal[i*9+2] += tempvalxz[k];
  2061. trgVal[i*9+3] += tempvalyx[k];
  2062. trgVal[i*9+4] += tempvalyy[k];
  2063. trgVal[i*9+5] += tempvalyz[k];
  2064. trgVal[i*9+6] += tempvalzx[k];
  2065. trgVal[i*9+7] += tempvalzy[k];
  2066. trgVal[i*9+8] += tempvalzz[k];
  2067. }
  2068. for (; j < ns; j++) {
  2069. double x = tx[i] - sx[j];
  2070. double y = ty[i] - sy[j];
  2071. double z = tz[i] - sz[j];
  2072. double r2 = x*x + y*y + z*z;
  2073. double r = sqrt(r2);
  2074. double invdr;
  2075. if (r == 0)
  2076. invdr = 0;
  2077. else
  2078. invdr = 1/r;
  2079. double invdr2=invdr*invdr;
  2080. double invdr3=invdr2*invdr;
  2081. double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]);
  2082. trgVal[i*9 ] += OOEP*oomeu*invdr3*( x*srcDen[j*3 ] - srcDen[j*3 ]*x + dot*(1-3*x*x*invdr2) );
  2083. trgVal[i*9+1] += OOEP*oomeu*invdr3*( y*srcDen[j*3 ] - srcDen[j*3+1]*x + dot*(0-3*y*x*invdr2) );
  2084. trgVal[i*9+2] += OOEP*oomeu*invdr3*( z*srcDen[j*3 ] - srcDen[j*3+2]*x + dot*(0-3*z*x*invdr2) );
  2085. trgVal[i*9+3] += OOEP*oomeu*invdr3*( x*srcDen[j*3+1] - srcDen[j*3 ]*y + dot*(0-3*x*y*invdr2) );
  2086. trgVal[i*9+4] += OOEP*oomeu*invdr3*( y*srcDen[j*3+1] - srcDen[j*3+1]*y + dot*(1-3*y*y*invdr2) );
  2087. trgVal[i*9+5] += OOEP*oomeu*invdr3*( z*srcDen[j*3+1] - srcDen[j*3+2]*y + dot*(0-3*z*y*invdr2) );
  2088. trgVal[i*9+6] += OOEP*oomeu*invdr3*( x*srcDen[j*3+2] - srcDen[j*3 ]*z + dot*(0-3*x*z*invdr2) );
  2089. trgVal[i*9+7] += OOEP*oomeu*invdr3*( y*srcDen[j*3+2] - srcDen[j*3+1]*z + dot*(0-3*y*z*invdr2) );
  2090. trgVal[i*9+8] += OOEP*oomeu*invdr3*( z*srcDen[j*3+2] - srcDen[j*3+2]*z + dot*(1-3*z*z*invdr2) );
  2091. }
  2092. }
  2093. return;
  2094. }
  2095. #undef SIMD_LEN
  2096. #define X(s,k) (s)[(k)*COORD_DIM]
  2097. #define Y(s,k) (s)[(k)*COORD_DIM+1]
  2098. #define Z(s,k) (s)[(k)*COORD_DIM+2]
  2099. void stokesDirectSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], const double kernel_coef, mem::MemoryManager* mem_mgr=NULL)
  2100. {
  2101. std::vector<double> xs(ns+1); std::vector<double> xt(nt);
  2102. std::vector<double> ys(ns+1); std::vector<double> yt(nt);
  2103. std::vector<double> zs(ns+1); std::vector<double> zt(nt);
  2104. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2105. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  2106. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2107. //1. reshuffle memory
  2108. for (int k =0;k<ns;k++){
  2109. xs[k+x_shift]=X(src,k);
  2110. ys[k+y_shift]=Y(src,k);
  2111. zs[k+z_shift]=Z(src,k);
  2112. }
  2113. for (int k=0;k<nt;k++){
  2114. xt[k]=X(trg,k);
  2115. yt[k]=Y(trg,k);
  2116. zt[k]=Z(trg,k);
  2117. }
  2118. //2. perform caclulation
  2119. stokesDirectVecSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot,kernel_coef);
  2120. return;
  2121. }
  2122. void stokesPressureSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
  2123. {
  2124. std::vector<double> xs(ns+1); std::vector<double> xt(nt);
  2125. std::vector<double> ys(ns+1); std::vector<double> yt(nt);
  2126. std::vector<double> zs(ns+1); std::vector<double> zt(nt);
  2127. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2128. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  2129. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2130. //1. reshuffle memory
  2131. for (int k =0;k<ns;k++){
  2132. xs[k+x_shift]=X(src,k);
  2133. ys[k+y_shift]=Y(src,k);
  2134. zs[k+z_shift]=Z(src,k);
  2135. }
  2136. for (int k=0;k<nt;k++){
  2137. xt[k]=X(trg,k);
  2138. yt[k]=Y(trg,k);
  2139. zt[k]=Z(trg,k);
  2140. }
  2141. //2. perform caclulation
  2142. stokesPressureSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
  2143. return;
  2144. }
  2145. void stokesStressSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
  2146. {
  2147. std::vector<double> xs(ns+1); std::vector<double> xt(nt);
  2148. std::vector<double> ys(ns+1); std::vector<double> yt(nt);
  2149. std::vector<double> zs(ns+1); std::vector<double> zt(nt);
  2150. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2151. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  2152. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2153. //1. reshuffle memory
  2154. for (int k =0;k<ns;k++){
  2155. xs[k+x_shift]=X(src,k);
  2156. ys[k+y_shift]=Y(src,k);
  2157. zs[k+z_shift]=Z(src,k);
  2158. }
  2159. for (int k=0;k<nt;k++){
  2160. xt[k]=X(trg,k);
  2161. yt[k]=Y(trg,k);
  2162. zt[k]=Z(trg,k);
  2163. }
  2164. //2. perform caclulation
  2165. stokesStressSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
  2166. return;
  2167. }
  2168. void stokesGradSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], const double kernel_coef, mem::MemoryManager* mem_mgr=NULL)
  2169. {
  2170. std::vector<double> xs(ns+1); std::vector<double> xt(nt);
  2171. std::vector<double> ys(ns+1); std::vector<double> yt(nt);
  2172. std::vector<double> zs(ns+1); std::vector<double> zt(nt);
  2173. int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2174. int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
  2175. int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
  2176. //1. reshuffle memory
  2177. for (int k =0;k<ns;k++){
  2178. xs[k+x_shift]=X(src,k);
  2179. ys[k+y_shift]=Y(src,k);
  2180. zs[k+z_shift]=Z(src,k);
  2181. }
  2182. for (int k=0;k<nt;k++){
  2183. xt[k]=X(trg,k);
  2184. yt[k]=Y(trg,k);
  2185. zt[k]=Z(trg,k);
  2186. }
  2187. //2. perform caclulation
  2188. stokesGradSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot,kernel_coef);
  2189. return;
  2190. }
  2191. #undef X
  2192. #undef Y
  2193. #undef Z
  2194. #undef IDEAL_ALIGNMENT
  2195. #undef DECL_SIMD_ALIGNED
  2196. }
  2197. template <>
  2198. void stokes_vel<double>(double* r_src, int src_cnt, double* v_src_, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  2199. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(28*dof));
  2200. const double mu=1.0;
  2201. stokesDirectSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mu, mem_mgr);
  2202. }
  2203. template <>
  2204. void stokes_press<double>(double* r_src, int src_cnt, double* v_src_, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  2205. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(17*dof));
  2206. stokesPressureSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mem_mgr);
  2207. return;
  2208. }
  2209. template <>
  2210. void stokes_stress<double>(double* r_src, int src_cnt, double* v_src_, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  2211. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(45*dof));
  2212. stokesStressSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mem_mgr);
  2213. }
  2214. template <>
  2215. void stokes_grad<double>(double* r_src, int src_cnt, double* v_src_, int dof, double* r_trg, int trg_cnt, double* k_out, mem::MemoryManager* mem_mgr){
  2216. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(89*dof));
  2217. const double mu=1.0;
  2218. stokesGradSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mu, mem_mgr);
  2219. }
  2220. #endif
  2221. #endif
  2222. ////////////////////////////////////////////////////////////////////////////////
  2223. //////// BIOT-SAVART KERNEL ////////
  2224. ////////////////////////////////////////////////////////////////////////////////
  2225. template <class T>
  2226. void biot_savart(T* r_src, int src_cnt, T* v_src_, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  2227. #ifndef __MIC__
  2228. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(26*dof));
  2229. #endif
  2230. const T OOFP = -1.0/(4.0*const_pi<T>());
  2231. for(int t=0;t<trg_cnt;t++){
  2232. for(int i=0;i<dof;i++){
  2233. T p[3]={0,0,0};
  2234. for(int s=0;s<src_cnt;s++){
  2235. T dR[3]={r_trg[3*t ]-r_src[3*s ],
  2236. r_trg[3*t+1]-r_src[3*s+1],
  2237. r_trg[3*t+2]-r_src[3*s+2]};
  2238. T R2 = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
  2239. if (R2!=0){
  2240. T invR2=1.0/R2;
  2241. T invR=sqrt(invR2);
  2242. T invR3=invR*invR2;
  2243. T v_src[3]={v_src_[(s*dof+i)*3 ],
  2244. v_src_[(s*dof+i)*3+1],
  2245. v_src_[(s*dof+i)*3+2]};
  2246. p[0] -= (v_src[1]*dR[2]-v_src[2]*dR[1])*invR3;
  2247. p[1] -= (v_src[2]*dR[0]-v_src[0]*dR[2])*invR3;
  2248. p[2] -= (v_src[0]*dR[1]-v_src[1]*dR[0])*invR3;
  2249. }
  2250. }
  2251. k_out[(t*dof+i)*3+0] += p[0]*OOFP;
  2252. k_out[(t*dof+i)*3+1] += p[1]*OOFP;
  2253. k_out[(t*dof+i)*3+2] += p[2]*OOFP;
  2254. }
  2255. }
  2256. }
  2257. ////////////////////////////////////////////////////////////////////////////////
  2258. //////// HELMHOLTZ KERNEL ////////
  2259. ////////////////////////////////////////////////////////////////////////////////
  2260. /**
  2261. * \brief Green's function for the Helmholtz's equation. Kernel tensor
  2262. * dimension = 2x2.
  2263. */
  2264. template <class T>
  2265. void helmholtz_poten(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  2266. #ifndef __MIC__
  2267. Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(24*dof));
  2268. #endif
  2269. const T mu = (20.0*const_pi<T>());
  2270. for(int t=0;t<trg_cnt;t++){
  2271. for(int i=0;i<dof;i++){
  2272. T p[2]={0,0};
  2273. for(int s=0;s<src_cnt;s++){
  2274. T dX_reg=r_trg[3*t ]-r_src[3*s ];
  2275. T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
  2276. T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
  2277. T R = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
  2278. if (R!=0){
  2279. R = sqrt(R);
  2280. T invR=1.0/R;
  2281. T G[2]={cos(mu*R)*invR, sin(mu*R)*invR};
  2282. p[0] += v_src[(s*dof+i)*2+0]*G[0] - v_src[(s*dof+i)*2+1]*G[1];
  2283. p[1] += v_src[(s*dof+i)*2+0]*G[1] + v_src[(s*dof+i)*2+1]*G[0];
  2284. }
  2285. }
  2286. k_out[(t*dof+i)*2+0] += p[0];
  2287. k_out[(t*dof+i)*2+1] += p[1];
  2288. }
  2289. }
  2290. }
  2291. template <class T>
  2292. void helmholtz_grad(T* r_src, int src_cnt, T* v_src, int dof, T* r_trg, int trg_cnt, T* k_out, mem::MemoryManager* mem_mgr){
  2293. //TODO Implement this.
  2294. }
  2295. }//end namespace