kernel.txx 77 KB

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