kernel.txx 82 KB

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