123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986 |
- /**
- * \file kernel.txx
- * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
- * \date 12-20-2011
- * \brief This file contains the implementation of the struct Kernel and also the
- * implementation of various kernels for FMM.
- */
- #ifdef USE_SSE
- #include <emmintrin.h>
- #endif
- #include <math.h>
- #include <assert.h>
- #include <vector>
- #include <profile.hpp>
- namespace pvfmm{
- /**
- * \brief Constructor.
- */
- template <class T>
- Kernel<T>::Kernel(): dim(0){
- ker_dim[0]=0;
- ker_dim[1]=0;
- }
- /**
- * \brief Constructor.
- */
- template <class T>
- Kernel<T>::Kernel(Ker_t poten, Ker_t dbl_poten, const char* name, int dim_,
- const int (&k_dim)[2], bool homogen_, T ker_scale,
- size_t dev_poten, size_t dev_dbl_poten){
- dim=dim_;
- ker_dim[0]=k_dim[0];
- ker_dim[1]=k_dim[1];
- ker_poten=poten;
- dbl_layer_poten=dbl_poten;
- homogen=homogen_;
- poten_scale=ker_scale;
- ker_name=std::string(name);
- dev_ker_poten=dev_poten;
- dev_dbl_layer_poten=dev_dbl_poten;
- }
- /**
- * \brief Compute the transformation matrix (on the source strength vector)
- * to get potential at target coordinates due to sources at the given
- * coordinates.
- * \param[in] r_src Coordinates of source points.
- * \param[in] src_cnt Number of source points.
- * \param[in] r_trg Coordinates of target points.
- * \param[in] trg_cnt Number of target points.
- * \param[out] k_out Output array with potential values.
- */
- template <class T>
- void Kernel<T>::BuildMatrix(T* r_src, int src_cnt,
- T* r_trg, int trg_cnt, T* k_out){
- int dim=3; //Only supporting 3D
- memset(k_out, 0, src_cnt*ker_dim[0]*trg_cnt*ker_dim[1]*sizeof(T));
- for(int i=0;i<src_cnt;i++) //TODO Optimize this.
- for(int j=0;j<ker_dim[0];j++){
- std::vector<T> v_src(ker_dim[0],0);
- v_src[j]=1.0;
- ker_poten(&r_src[i*dim], 1, &v_src[0], 1, r_trg, trg_cnt,
- &k_out[(i*ker_dim[0]+j)*trg_cnt*ker_dim[1]], NULL);
- }
- }
- ////////////////////////////////////////////////////////////////////////////////
- //////// LAPLACE KERNEL ////////
- ////////////////////////////////////////////////////////////////////////////////
- /**
- * \brief Green's function for the Poisson's equation. Kernel tensor
- * dimension = 1x1.
- */
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(12*dof));
- #endif
- const T OOFP = 1.0/(4.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p += v_src[s*dof+i]*invR;
- }
- k_out[t*dof+i] += p*OOFP;
- }
- }
- }
- template <class T>
- 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){
- //void laplace_poten(T* r_src_, int src_cnt, T* v_src_, int dof, T* r_trg_, int trg_cnt, T* k_out_){
- // int dim=3; //Only supporting 3D
- // T* r_src=new T[src_cnt*dim];
- // T* r_trg=new T[trg_cnt*dim];
- // T* v_src=new T[src_cnt ];
- // T* k_out=new T[trg_cnt ];
- // mem::memcopy(r_src,r_src_,src_cnt*dim*sizeof(T));
- // mem::memcopy(r_trg,r_trg_,trg_cnt*dim*sizeof(T));
- // mem::memcopy(v_src,v_src_,src_cnt *sizeof(T));
- // mem::memcopy(k_out,k_out_,trg_cnt *sizeof(T));
- #define EVAL_BLKSZ 32
- #define MAX_DOF 100
- //Compute source to target interactions.
- const T OOFP = 1.0/(4.0*const_pi<T>());
- if(dof==1){
- for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
- for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
- int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
- int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
- for(int t=t_;t<trg_blk;t++){
- T p=0;
- for(int s=s_;s<src_blk;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p += v_src[s]*invR;
- }
- k_out[t] += p*OOFP;
- }
- }
- }else if(dof==2){
- T p[MAX_DOF];
- for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
- for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
- int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
- int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
- for(int t=t_;t<trg_blk;t++){
- p[0]=0; p[1]=0;
- for(int s=s_;s<src_blk;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p[0] += v_src[s*dof+0]*invR;
- p[1] += v_src[s*dof+1]*invR;
- }
- k_out[t*dof+0] += p[0]*OOFP;
- k_out[t*dof+1] += p[1]*OOFP;
- }
- }
- }else if(dof==3){
- T p[MAX_DOF];
- for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
- for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
- int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
- int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
- for(int t=t_;t<trg_blk;t++){
- p[0]=0; p[1]=0; p[2]=0;
- for(int s=s_;s<src_blk;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p[0] += v_src[s*dof+0]*invR;
- p[1] += v_src[s*dof+1]*invR;
- p[2] += v_src[s*dof+2]*invR;
- }
- k_out[t*dof+0] += p[0]*OOFP;
- k_out[t*dof+1] += p[1]*OOFP;
- k_out[t*dof+2] += p[2]*OOFP;
- }
- }
- }else{
- T p[MAX_DOF];
- for (int t_=0; t_<trg_cnt; t_+=EVAL_BLKSZ)
- for (int s_=0; s_<src_cnt; s_+=EVAL_BLKSZ){
- int src_blk=s_+EVAL_BLKSZ; src_blk=(src_blk>src_cnt?src_cnt:src_blk);
- int trg_blk=t_+EVAL_BLKSZ; trg_blk=(trg_blk>trg_cnt?trg_cnt:trg_blk);
- for(int t=t_;t<trg_blk;t++){
- for(int i=0;i<dof;i++) p[i]=0;
- for(int s=s_;s<src_blk;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- for(int i=0;i<dof;i++)
- p[i] += v_src[s*dof+i]*invR;
- }
- for(int i=0;i<dof;i++)
- k_out[t*dof+i] += p[i]*OOFP;
- }
- }
- }
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+2*dof));
- #endif
- #undef MAX_DOF
- #undef EVAL_BLKSZ
- // for (int t=0; t<trg_cnt; t++)
- // k_out_[t] += k_out[t];
- // delete[] r_src;
- // delete[] r_trg;
- // delete[] v_src;
- // delete[] k_out;
- }
- // Laplace double layer potential.
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(19*dof));
- #endif
- const T OOFP = -1.0/(4.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p = v_src[(s*dof+i)*4+3]*invR*invR*invR;
- k_out[t*dof+i] += p*OOFP*( dX_reg*v_src[(s*dof+i)*4+0] +
- dY_reg*v_src[(s*dof+i)*4+1] +
- dZ_reg*v_src[(s*dof+i)*4+2] );
- }
- }
- }
- }
- // Laplace grdient kernel.
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+12*dof));
- #endif
- const T OOFP = -1.0/(4.0*const_pi<T>());
- if(dof==1){
- for(int t=0;t<trg_cnt;t++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p = v_src[s]*invR*invR*invR;
- k_out[(t)*3+0] += p*OOFP*dX_reg;
- k_out[(t)*3+1] += p*OOFP*dY_reg;
- k_out[(t)*3+2] += p*OOFP*dZ_reg;
- }
- }
- }else if(dof==2){
- for(int t=0;t<trg_cnt;t++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p = v_src[s*dof+0]*invR*invR*invR;
- k_out[(t*dof+0)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+0)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+0)*3+2] += p*OOFP*dZ_reg;
- p = v_src[s*dof+1]*invR*invR*invR;
- k_out[(t*dof+1)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+1)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+1)*3+2] += p*OOFP*dZ_reg;
- }
- }
- }else if(dof==3){
- for(int t=0;t<trg_cnt;t++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p = v_src[s*dof+0]*invR*invR*invR;
- k_out[(t*dof+0)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+0)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+0)*3+2] += p*OOFP*dZ_reg;
- p = v_src[s*dof+1]*invR*invR*invR;
- k_out[(t*dof+1)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+1)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+1)*3+2] += p*OOFP*dZ_reg;
- p = v_src[s*dof+2]*invR*invR*invR;
- k_out[(t*dof+2)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+2)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+2)*3+2] += p*OOFP*dZ_reg;
- }
- }
- }else{
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T invR = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (invR!=0) invR = 1.0/sqrt(invR);
- p = v_src[s*dof+i]*invR*invR*invR;
- k_out[(t*dof+i)*3+0] += p*OOFP*dX_reg;
- k_out[(t*dof+i)*3+1] += p*OOFP*dY_reg;
- k_out[(t*dof+i)*3+2] += p*OOFP*dZ_reg;
- }
- }
- }
- }
- }
- #ifndef __MIC__
- #ifdef USE_SSE
- namespace
- {
- #define IDEAL_ALIGNMENT 16
- #define SIMD_LEN (int)(IDEAL_ALIGNMENT / sizeof(double))
- #define DECL_SIMD_ALIGNED __declspec(align(IDEAL_ALIGNMENT))
- void laplaceSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal)
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double OOFP = 1.0/(4.0*const_pi<double>());
- __m128d temp;
- double aux_arr[SIMD_LEN+1];
- double *tempval;
- // if aux_arr is misaligned
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempval = aux_arr + 1;
- else tempval = aux_arr;
- if (size_t(tempval)%IDEAL_ALIGNMENT) abort();
- /*! One over four pi */
- __m128d oofp = _mm_set1_pd (OOFP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- temp = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sden = _mm_set_pd (srcDen[j+1], srcDen[j]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (reqzero, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- sden = _mm_mul_pd (sden, S);
- temp = _mm_add_pd (sden, temp);
- }
- temp = _mm_mul_pd (temp, oofp);
- _mm_store_pd(tempval, temp);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i] += tempval[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double den = srcDen[j];
- trgVal[i] += den*invdr*OOFP;
- }
- }
- return;
- }
- void laplaceDblSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal)
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double OOFP = 1.0/(4.0*const_pi<double>());
- __m128d temp;
- double aux_arr[SIMD_LEN+1];
- double *tempval;
- // if aux_arr is misaligned
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempval = aux_arr + 1;
- else tempval = aux_arr;
- if (size_t(tempval)%IDEAL_ALIGNMENT) abort();
- /*! One over four pi */
- __m128d oofp = _mm_set1_pd (OOFP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- temp = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d snormx = _mm_set_pd (srcDen[(j+1)*4+0], srcDen[j*4+0]);
- __m128d snormy = _mm_set_pd (srcDen[(j+1)*4+1], srcDen[j*4+1]);
- __m128d snormz = _mm_set_pd (srcDen[(j+1)*4+2], srcDen[j*4+2]);
- __m128d sden = _mm_set_pd (srcDen[(j+1)*4+3], srcDen[j*4+3]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- __m128d S2;
- __m128d S3;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (reqzero, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- S2 = _mm_mul_pd (S, S);
- S3 = _mm_mul_pd (S2, S);
- __m128d S3_sden=_mm_mul_pd(S3, sden);
- __m128d dot_sum = _mm_add_pd(_mm_mul_pd(snormx,dX),_mm_mul_pd(snormy,dY));
- dot_sum = _mm_add_pd(dot_sum,_mm_mul_pd(snormz,dZ));
- temp = _mm_add_pd(_mm_mul_pd(S3_sden,dot_sum),temp);
- }
- temp = _mm_mul_pd (temp, oofp);
- _mm_store_pd(tempval, temp);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i] += tempval[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double invdr2=invdr*invdr;
- double invdr3=invdr2*invdr;
- double dot_sum = x*srcDen[j*4+0] + y*srcDen[j*4+1] + z*srcDen[j*4+2];
- trgVal[i] += OOFP*invdr3*x*srcDen[j*4+3]*dot_sum;
- }
- }
- return;
- }
- void laplaceGradSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal)
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double OOFP = 1.0/(4.0*const_pi<double>());
- __m128d tempx; __m128d tempy; __m128d tempz;
- double aux_arr[3*SIMD_LEN+1];
- double *tempvalx, *tempvaly, *tempvalz;
- // if aux_arr is misaligned
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) tempvalx = aux_arr + 1;
- else tempvalx = aux_arr;
- if (size_t(tempvalx)%IDEAL_ALIGNMENT) abort();
- tempvaly=tempvalx+SIMD_LEN;
- tempvalz=tempvaly+SIMD_LEN;
- /*! One over four pi */
- __m128d oofp = _mm_set1_pd (OOFP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- tempx = _mm_setzero_pd();
- tempy = _mm_setzero_pd();
- tempz = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sden = _mm_set_pd (srcDen[j+1], srcDen[j]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- __m128d S2;
- __m128d S3;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d reqzero = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (reqzero, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- S2 = _mm_mul_pd (S, S);
- S3 = _mm_mul_pd (S2, S);
- __m128d S3_sden=_mm_mul_pd(S3, sden);
- tempx = _mm_add_pd(_mm_mul_pd(S3_sden,dX),tempx);
- tempy = _mm_add_pd(_mm_mul_pd(S3_sden,dY),tempy);
- tempz = _mm_add_pd(_mm_mul_pd(S3_sden,dZ),tempz);
- }
- tempx = _mm_mul_pd (tempx, oofp);
- tempy = _mm_mul_pd (tempy, oofp);
- tempz = _mm_mul_pd (tempz, oofp);
- _mm_store_pd(tempvalx, tempx);
- _mm_store_pd(tempvaly, tempy);
- _mm_store_pd(tempvalz, tempz);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i*3 ] += tempvalx[k];
- trgVal[i*3+1] += tempvaly[k];
- trgVal[i*3+2] += tempvalz[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double invdr2=invdr*invdr;
- double invdr3=invdr2*invdr;
- trgVal[i*3 ] += OOFP*invdr3*x*srcDen[j];
- trgVal[i*3+1] += OOFP*invdr3*y*srcDen[j];
- trgVal[i*3+2] += OOFP*invdr3*z*srcDen[j];
- }
- }
- return;
- }
- #undef SIMD_LEN
- #define X(s,k) (s)[(k)*COORD_DIM]
- #define Y(s,k) (s)[(k)*COORD_DIM+1]
- #define Z(s,k) (s)[(k)*COORD_DIM+2]
- void laplaceSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- // TODO
- }
- void laplaceSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- double* buff=NULL;
- if(mem_mgr) buff=(double*)mem_mgr->malloc((ns+1+nt)*3*sizeof(double));
- else buff=(double*)malloc((ns+1+nt)*3*sizeof(double));
- double* buff_=buff;
- pvfmm::Vector<double> xs(ns+1,buff_,false); buff_+=ns+1;
- pvfmm::Vector<double> ys(ns+1,buff_,false); buff_+=ns+1;
- pvfmm::Vector<double> zs(ns+1,buff_,false); buff_+=ns+1;
- pvfmm::Vector<double> xt(nt ,buff_,false); buff_+=nt ;
- pvfmm::Vector<double> yt(nt ,buff_,false); buff_+=nt ;
- pvfmm::Vector<double> zt(nt ,buff_,false); buff_+=nt ;
- //std::vector<double> xs(ns+1);
- //std::vector<double> ys(ns+1);
- //std::vector<double> zs(ns+1);
- //std::vector<double> xt(nt );
- //std::vector<double> yt(nt );
- //std::vector<double> zt(nt );
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- laplaceSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
- if(mem_mgr) mem_mgr->free(buff);
- else free(buff);
- return;
- }
- void laplaceDblSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- // TODO
- }
- void laplaceDblSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- std::vector<double> xs(ns+1); std::vector<double> xt(nt);
- std::vector<double> ys(ns+1); std::vector<double> yt(nt);
- std::vector<double> zs(ns+1); std::vector<double> zt(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- laplaceDblSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
- return;
- }
- void laplaceGradSSEShuffle(const int ns, const int nt, float const src[], float const trg[], float const den[], float pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- // TODO
- }
- void laplaceGradSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- int tid=omp_get_thread_num();
- static std::vector<std::vector<double> > xs_(100); static std::vector<std::vector<double> > xt_(100);
- static std::vector<std::vector<double> > ys_(100); static std::vector<std::vector<double> > yt_(100);
- static std::vector<std::vector<double> > zs_(100); static std::vector<std::vector<double> > zt_(100);
- std::vector<double>& xs=xs_[tid]; std::vector<double>& xt=xt_[tid];
- std::vector<double>& ys=ys_[tid]; std::vector<double>& yt=yt_[tid];
- std::vector<double>& zs=zs_[tid]; std::vector<double>& zt=zt_[tid];
- xs.resize(ns+1); xt.resize(nt);
- ys.resize(ns+1); yt.resize(nt);
- zs.resize(ns+1); zt.resize(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- laplaceGradSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
- return;
- }
- #undef X
- #undef Y
- #undef Z
- #undef IDEAL_ALIGNMENT
- #undef DECL_SIMD_ALIGNED
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(12*dof));
- if(dof==1){
- laplaceSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
- return;
- }
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(19*dof));
- if(dof==1){
- laplaceDblSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
- return;
- }
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(10+12*dof));
- if(dof==1){
- laplaceGradSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src, k_out, mem_mgr);
- return;
- }
- }
- #endif
- #endif
- ////////////////////////////////////////////////////////////////////////////////
- //////// STOKES KERNEL ////////
- ////////////////////////////////////////////////////////////////////////////////
- /**
- * \brief Green's function for the Stokes's equation. Kernel tensor
- * dimension = 3x3.
- */
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(28*dof));
- #endif
- const T mu=1.0;
- const T OOEPMU = 1.0/(8.0*const_pi<T>()*mu);
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[3]={0,0,0};
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R!=0){
- T invR2=1.0/R;
- T invR=sqrt(invR2);
- T v_src[3]={v_src_[(s*dof+i)*3 ],
- v_src_[(s*dof+i)*3+1],
- v_src_[(s*dof+i)*3+2]};
- T inner_prod=(v_src[0]*dR[0] +
- v_src[1]*dR[1] +
- v_src[2]*dR[2])* invR2;
- p[0] += (v_src[0] + dR[0]*inner_prod)*invR;
- p[1] += (v_src[1] + dR[1]*inner_prod)*invR;
- p[2] += (v_src[2] + dR[2]*inner_prod)*invR;
- }
- }
- k_out[(t*dof+i)*3+0] += p[0]*OOEPMU;
- k_out[(t*dof+i)*3+1] += p[1]*OOEPMU;
- k_out[(t*dof+i)*3+2] += p[2]*OOEPMU;
- }
- }
- }
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(47*dof));
- #endif
- const T mu=1.0;
- const T OOEPMU = -1.0/(8.0*const_pi<T>()*mu);
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[3]={0,0,0};
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R!=0){
- T invR2=1.0/R;
- T invR=sqrt(invR2);
- T invR3=invR2*invR;
- T* f=&v_src[(s*dof+i)*6+0];
- T* n=&v_src[(s*dof+i)*6+3];
- T r_dot_n=(n[0]*dR[0]+n[1]*dR[1]+n[2]*dR[2]);
- T r_dot_f=(f[0]*dR[0]+f[1]*dR[1]+f[2]*dR[2]);
- T n_dot_f=(f[0]* n[0]+f[1]* n[1]+f[2]* n[2]);
- p[0] += dR[0]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
- p[1] += dR[1]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
- p[2] += dR[2]*(n_dot_f - 3*r_dot_n*r_dot_f*invR2)*invR3;
- }
- }
- k_out[(t*dof+i)*3+0] += p[0]*OOEPMU;
- k_out[(t*dof+i)*3+1] += p[1]*OOEPMU;
- k_out[(t*dof+i)*3+2] += p[2]*OOEPMU;
- }
- }
- }
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(17*dof));
- #endif
- const T OOFP = 1.0/(4.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p=0;
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R!=0){
- T invR2=1.0/R;
- T invR=sqrt(invR2);
- T invR3=invR2*invR;
- T v_src[3]={v_src_[(s*dof+i)*3 ],
- v_src_[(s*dof+i)*3+1],
- v_src_[(s*dof+i)*3+2]};
- T inner_prod=(v_src[0]*dR[0] +
- v_src[1]*dR[1] +
- v_src[2]*dR[2])* invR3;
- p += inner_prod;
- }
- }
- k_out[t*dof+i] += p*OOFP;
- }
- }
- }
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(45*dof));
- #endif
- const T TOFP = -3.0/(4.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[9]={0,0,0,
- 0,0,0,
- 0,0,0};
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R!=0){
- T invR2=1.0/R;
- T invR=sqrt(invR2);
- T invR3=invR2*invR;
- T invR5=invR3*invR2;
- T v_src[3]={v_src_[(s*dof+i)*3 ],
- v_src_[(s*dof+i)*3+1],
- v_src_[(s*dof+i)*3+2]};
- T inner_prod=(v_src[0]*dR[0] +
- v_src[1]*dR[1] +
- v_src[2]*dR[2])* invR5;
- p[0] += inner_prod*dR[0]*dR[0]; p[1] += inner_prod*dR[1]*dR[0]; p[2] += inner_prod*dR[2]*dR[0];
- p[3] += inner_prod*dR[0]*dR[1]; p[4] += inner_prod*dR[1]*dR[1]; p[5] += inner_prod*dR[2]*dR[1];
- p[6] += inner_prod*dR[0]*dR[2]; p[7] += inner_prod*dR[1]*dR[2]; p[8] += inner_prod*dR[2]*dR[2];
- }
- }
- k_out[(t*dof+i)*9+0] += p[0]*TOFP;
- k_out[(t*dof+i)*9+1] += p[1]*TOFP;
- k_out[(t*dof+i)*9+2] += p[2]*TOFP;
- k_out[(t*dof+i)*9+3] += p[3]*TOFP;
- k_out[(t*dof+i)*9+4] += p[4]*TOFP;
- k_out[(t*dof+i)*9+5] += p[5]*TOFP;
- k_out[(t*dof+i)*9+6] += p[6]*TOFP;
- k_out[(t*dof+i)*9+7] += p[7]*TOFP;
- k_out[(t*dof+i)*9+8] += p[8]*TOFP;
- }
- }
- }
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(89*dof));
- #endif
- const T mu=1.0;
- const T OOEPMU = 1.0/(8.0*const_pi<T>()*mu);
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[9]={0,0,0,
- 0,0,0,
- 0,0,0};
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R!=0){
- T invR2=1.0/R;
- T invR=sqrt(invR2);
- T invR3=invR2*invR;
- T v_src[3]={v_src_[(s*dof+i)*3 ],
- v_src_[(s*dof+i)*3+1],
- v_src_[(s*dof+i)*3+2]};
- T inner_prod=(v_src[0]*dR[0] +
- v_src[1]*dR[1] +
- v_src[2]*dR[2]);
- p[0] += ( inner_prod*(1-3*dR[0]*dR[0]*invR2))*invR3; //6
- p[1] += (dR[1]*v_src[0]-v_src[1]*dR[0]+inner_prod*( -3*dR[1]*dR[0]*invR2))*invR3; //9
- p[2] += (dR[2]*v_src[0]-v_src[2]*dR[0]+inner_prod*( -3*dR[2]*dR[0]*invR2))*invR3;
- p[3] += (dR[0]*v_src[1]-v_src[0]*dR[1]+inner_prod*( -3*dR[0]*dR[1]*invR2))*invR3;
- p[4] += ( inner_prod*(1-3*dR[1]*dR[1]*invR2))*invR3;
- p[5] += (dR[2]*v_src[1]-v_src[2]*dR[1]+inner_prod*( -3*dR[2]*dR[1]*invR2))*invR3;
- p[6] += (dR[0]*v_src[2]-v_src[0]*dR[2]+inner_prod*( -3*dR[0]*dR[2]*invR2))*invR3;
- p[7] += (dR[1]*v_src[2]-v_src[1]*dR[2]+inner_prod*( -3*dR[1]*dR[2]*invR2))*invR3;
- p[8] += ( inner_prod*(1-3*dR[2]*dR[2]*invR2))*invR3;
- }
- }
- k_out[(t*dof+i)*9+0] += p[0]*OOEPMU;
- k_out[(t*dof+i)*9+1] += p[1]*OOEPMU;
- k_out[(t*dof+i)*9+2] += p[2]*OOEPMU;
- k_out[(t*dof+i)*9+3] += p[3]*OOEPMU;
- k_out[(t*dof+i)*9+4] += p[4]*OOEPMU;
- k_out[(t*dof+i)*9+5] += p[5]*OOEPMU;
- k_out[(t*dof+i)*9+6] += p[6]*OOEPMU;
- k_out[(t*dof+i)*9+7] += p[7]*OOEPMU;
- k_out[(t*dof+i)*9+8] += p[8]*OOEPMU;
- }
- }
- }
- #ifndef __MIC__
- #ifdef USE_SSE
- namespace
- {
- #define IDEAL_ALIGNMENT 16
- #define SIMD_LEN (int)(IDEAL_ALIGNMENT / sizeof(double))
- #define DECL_SIMD_ALIGNED __declspec(align(IDEAL_ALIGNMENT))
- void stokesDirectVecSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal,
- const double cof )
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double mu = cof;
- double OOEP = 1.0/(8.0*const_pi<double>());
- __m128d tempx;
- __m128d tempy;
- __m128d tempz;
- double oomeu = 1/mu;
- double aux_arr[3*SIMD_LEN+1];
- double *tempvalx;
- double *tempvaly;
- double *tempvalz;
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
- {
- tempvalx = aux_arr + 1;
- if (size_t(tempvalx)%IDEAL_ALIGNMENT)
- abort();
- }
- else
- tempvalx = aux_arr;
- tempvaly=tempvalx+SIMD_LEN;
- tempvalz=tempvaly+SIMD_LEN;
- /*! One over eight pi */
- __m128d ooep = _mm_set1_pd (OOEP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- __m128d oomu = _mm_set1_pd (1/mu);
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- tempx = _mm_setzero_pd();
- tempy = _mm_setzero_pd();
- tempz = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
- __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
- __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d temp = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (temp, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- __m128d dotx = _mm_mul_pd (dX, sdenx);
- __m128d doty = _mm_mul_pd (dY, sdeny);
- __m128d dotz = _mm_mul_pd (dZ, sdenz);
- __m128d dot_sum = _mm_add_pd (dotx, doty);
- dot_sum = _mm_add_pd (dot_sum, dotz);
- dot_sum = _mm_mul_pd (dot_sum, S);
- dot_sum = _mm_mul_pd (dot_sum, S);
- dotx = _mm_mul_pd (dot_sum, dX);
- doty = _mm_mul_pd (dot_sum, dY);
- dotz = _mm_mul_pd (dot_sum, dZ);
- sdenx = _mm_add_pd (sdenx, dotx);
- sdeny = _mm_add_pd (sdeny, doty);
- sdenz = _mm_add_pd (sdenz, dotz);
- sdenx = _mm_mul_pd (sdenx, S);
- sdeny = _mm_mul_pd (sdeny, S);
- sdenz = _mm_mul_pd (sdenz, S);
- tempx = _mm_add_pd (sdenx, tempx);
- tempy = _mm_add_pd (sdeny, tempy);
- tempz = _mm_add_pd (sdenz, tempz);
- }
- tempx = _mm_mul_pd (tempx, ooep);
- tempy = _mm_mul_pd (tempy, ooep);
- tempz = _mm_mul_pd (tempz, ooep);
- tempx = _mm_mul_pd (tempx, oomu);
- tempy = _mm_mul_pd (tempy, oomu);
- tempz = _mm_mul_pd (tempz, oomu);
- _mm_store_pd(tempvalx, tempx);
- _mm_store_pd(tempvaly, tempy);
- _mm_store_pd(tempvalz, tempz);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i*3] += tempvalx[k];
- trgVal[i*3+1] += tempvaly[k];
- trgVal[i*3+2] += tempvalz[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr * invdr;
- double denx = srcDen[j*3] + dot*x;
- double deny = srcDen[j*3+1] + dot*y;
- double denz = srcDen[j*3+2] + dot*z;
- trgVal[i*3] += denx*invdr*OOEP*oomeu;
- trgVal[i*3+1] += deny*invdr*OOEP*oomeu;
- trgVal[i*3+2] += denz*invdr*OOEP*oomeu;
- }
- }
- return;
- }
- void stokesPressureSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal)
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double OOFP = 1.0/(4.0*const_pi<double>());
- __m128d temp_press;
- double aux_arr[SIMD_LEN+1];
- double *tempval_press;
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
- {
- tempval_press = aux_arr + 1;
- if (size_t(tempval_press)%IDEAL_ALIGNMENT)
- abort();
- }
- else
- tempval_press = aux_arr;
- /*! One over eight pi */
- __m128d oofp = _mm_set1_pd (OOFP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- temp_press = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
- __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
- __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d temp = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (temp, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- __m128d dotx = _mm_mul_pd (dX, sdenx);
- __m128d doty = _mm_mul_pd (dY, sdeny);
- __m128d dotz = _mm_mul_pd (dZ, sdenz);
- __m128d dot_sum = _mm_add_pd (dotx, doty);
- dot_sum = _mm_add_pd (dot_sum, dotz);
- dot_sum = _mm_mul_pd (dot_sum, S);
- dot_sum = _mm_mul_pd (dot_sum, S);
- dot_sum = _mm_mul_pd (dot_sum, S);
- temp_press = _mm_add_pd (dot_sum, temp_press);
- }
- temp_press = _mm_mul_pd (temp_press, oofp);
- _mm_store_pd(tempval_press, temp_press);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i] += tempval_press[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr * invdr * invdr;
- trgVal[i] += dot*OOFP;
- }
- }
- return;
- }
- void stokesStressSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal)
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double TOFP = -3.0/(4.0*const_pi<double>());
- __m128d tempxx; __m128d tempxy; __m128d tempxz;
- __m128d tempyx; __m128d tempyy; __m128d tempyz;
- __m128d tempzx; __m128d tempzy; __m128d tempzz;
- double aux_arr[9*SIMD_LEN+1];
- double *tempvalxx, *tempvalxy, *tempvalxz;
- double *tempvalyx, *tempvalyy, *tempvalyz;
- double *tempvalzx, *tempvalzy, *tempvalzz;
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
- {
- tempvalxx = aux_arr + 1;
- if (size_t(tempvalxx)%IDEAL_ALIGNMENT)
- abort();
- }
- else
- tempvalxx = aux_arr;
- tempvalxy=tempvalxx+SIMD_LEN;
- tempvalxz=tempvalxy+SIMD_LEN;
- tempvalyx=tempvalxz+SIMD_LEN;
- tempvalyy=tempvalyx+SIMD_LEN;
- tempvalyz=tempvalyy+SIMD_LEN;
- tempvalzx=tempvalyz+SIMD_LEN;
- tempvalzy=tempvalzx+SIMD_LEN;
- tempvalzz=tempvalzy+SIMD_LEN;
- /*! One over eight pi */
- __m128d tofp = _mm_set1_pd (TOFP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d zero = _mm_setzero_pd ();
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- tempxx = _mm_setzero_pd(); tempxy = _mm_setzero_pd(); tempxz = _mm_setzero_pd();
- tempyx = _mm_setzero_pd(); tempyy = _mm_setzero_pd(); tempyz = _mm_setzero_pd();
- tempzx = _mm_setzero_pd(); tempzy = _mm_setzero_pd(); tempzz = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
- __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
- __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- __m128d S2;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d temp = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (temp, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- S2 = _mm_mul_pd (S, S);
- __m128d dotx = _mm_mul_pd (dX, sdenx);
- __m128d doty = _mm_mul_pd (dY, sdeny);
- __m128d dotz = _mm_mul_pd (dZ, sdenz);
- __m128d dot_sum = _mm_add_pd (dotx, doty);
- dot_sum = _mm_add_pd (dot_sum, dotz);
- dot_sum = _mm_mul_pd (dot_sum, S);
- dot_sum = _mm_mul_pd (dot_sum, S2);
- dot_sum = _mm_mul_pd (dot_sum, S2);
- dotx = _mm_mul_pd (dot_sum, dX);
- doty = _mm_mul_pd (dot_sum, dY);
- dotz = _mm_mul_pd (dot_sum, dZ);
- tempxx = _mm_add_pd (_mm_mul_pd(dotx,dX), tempxx);
- tempxy = _mm_add_pd (_mm_mul_pd(dotx,dY), tempxy);
- tempxz = _mm_add_pd (_mm_mul_pd(dotx,dZ), tempxz);
- tempyx = _mm_add_pd (_mm_mul_pd(doty,dX), tempyx);
- tempyy = _mm_add_pd (_mm_mul_pd(doty,dY), tempyy);
- tempyz = _mm_add_pd (_mm_mul_pd(doty,dZ), tempyz);
- tempzx = _mm_add_pd (_mm_mul_pd(dotz,dX), tempzx);
- tempzy = _mm_add_pd (_mm_mul_pd(dotz,dY), tempzy);
- tempzz = _mm_add_pd (_mm_mul_pd(dotz,dZ), tempzz);
- }
- tempxx = _mm_mul_pd (tempxx, tofp);
- tempxy = _mm_mul_pd (tempxy, tofp);
- tempxz = _mm_mul_pd (tempxz, tofp);
- tempyx = _mm_mul_pd (tempyx, tofp);
- tempyy = _mm_mul_pd (tempyy, tofp);
- tempyz = _mm_mul_pd (tempyz, tofp);
- tempzx = _mm_mul_pd (tempzx, tofp);
- tempzy = _mm_mul_pd (tempzy, tofp);
- tempzz = _mm_mul_pd (tempzz, tofp);
- _mm_store_pd(tempvalxx, tempxx); _mm_store_pd(tempvalxy, tempxy); _mm_store_pd(tempvalxz, tempxz);
- _mm_store_pd(tempvalyx, tempyx); _mm_store_pd(tempvalyy, tempyy); _mm_store_pd(tempvalyz, tempyz);
- _mm_store_pd(tempvalzx, tempzx); _mm_store_pd(tempvalzy, tempzy); _mm_store_pd(tempvalzz, tempzz);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i*9 ] += tempvalxx[k];
- trgVal[i*9+1] += tempvalxy[k];
- trgVal[i*9+2] += tempvalxz[k];
- trgVal[i*9+3] += tempvalyx[k];
- trgVal[i*9+4] += tempvalyy[k];
- trgVal[i*9+5] += tempvalyz[k];
- trgVal[i*9+6] += tempvalzx[k];
- trgVal[i*9+7] += tempvalzy[k];
- trgVal[i*9+8] += tempvalzz[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double invdr2=invdr*invdr;
- double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]) * invdr2 * invdr2 * invdr;
- double denx = dot*x;
- double deny = dot*y;
- double denz = dot*z;
- trgVal[i*9 ] += denx*x*TOFP;
- trgVal[i*9+1] += denx*y*TOFP;
- trgVal[i*9+2] += denx*z*TOFP;
- trgVal[i*9+3] += deny*x*TOFP;
- trgVal[i*9+4] += deny*y*TOFP;
- trgVal[i*9+5] += deny*z*TOFP;
- trgVal[i*9+6] += denz*x*TOFP;
- trgVal[i*9+7] += denz*y*TOFP;
- trgVal[i*9+8] += denz*z*TOFP;
- }
- }
- return;
- }
- void stokesGradSSE(
- const int ns,
- const int nt,
- const double *sx,
- const double *sy,
- const double *sz,
- const double *tx,
- const double *ty,
- const double *tz,
- const double *srcDen,
- double *trgVal,
- const double cof )
- {
- if ( size_t(sx)%IDEAL_ALIGNMENT || size_t(sy)%IDEAL_ALIGNMENT || size_t(sz)%IDEAL_ALIGNMENT )
- abort();
- double mu = cof;
- double OOEP = 1.0/(8.0*const_pi<double>());
- __m128d tempxx; __m128d tempxy; __m128d tempxz;
- __m128d tempyx; __m128d tempyy; __m128d tempyz;
- __m128d tempzx; __m128d tempzy; __m128d tempzz;
- double oomeu = 1/mu;
- double aux_arr[9*SIMD_LEN+1];
- double *tempvalxx, *tempvalxy, *tempvalxz;
- double *tempvalyx, *tempvalyy, *tempvalyz;
- double *tempvalzx, *tempvalzy, *tempvalzz;
- if (size_t(aux_arr)%IDEAL_ALIGNMENT) // if aux_arr is misaligned
- {
- tempvalxx = aux_arr + 1;
- if (size_t(tempvalxx)%IDEAL_ALIGNMENT)
- abort();
- }
- else
- tempvalxx = aux_arr;
- tempvalxy=tempvalxx+SIMD_LEN;
- tempvalxz=tempvalxy+SIMD_LEN;
- tempvalyx=tempvalxz+SIMD_LEN;
- tempvalyy=tempvalyx+SIMD_LEN;
- tempvalyz=tempvalyy+SIMD_LEN;
- tempvalzx=tempvalyz+SIMD_LEN;
- tempvalzy=tempvalzx+SIMD_LEN;
- tempvalzz=tempvalzy+SIMD_LEN;
- /*! One over eight pi */
- __m128d ooep = _mm_set1_pd (OOEP);
- __m128d half = _mm_set1_pd (0.5);
- __m128d opf = _mm_set1_pd (1.5);
- __m128d three = _mm_set1_pd (3.0);
- __m128d zero = _mm_setzero_pd ();
- __m128d oomu = _mm_set1_pd (1/mu);
- __m128d ooepmu = _mm_mul_pd(ooep,oomu);
- // loop over sources
- int i = 0;
- for (; i < nt; i++) {
- tempxx = _mm_setzero_pd(); tempxy = _mm_setzero_pd(); tempxz = _mm_setzero_pd();
- tempyx = _mm_setzero_pd(); tempyy = _mm_setzero_pd(); tempyz = _mm_setzero_pd();
- tempzx = _mm_setzero_pd(); tempzy = _mm_setzero_pd(); tempzz = _mm_setzero_pd();
- __m128d txi = _mm_load1_pd (&tx[i]);
- __m128d tyi = _mm_load1_pd (&ty[i]);
- __m128d tzi = _mm_load1_pd (&tz[i]);
- int j = 0;
- // Load and calculate in groups of SIMD_LEN
- for (; j + SIMD_LEN <= ns; j+=SIMD_LEN) {
- __m128d sxj = _mm_load_pd (&sx[j]);
- __m128d syj = _mm_load_pd (&sy[j]);
- __m128d szj = _mm_load_pd (&sz[j]);
- __m128d sdenx = _mm_set_pd (srcDen[(j+1)*3], srcDen[j*3]);
- __m128d sdeny = _mm_set_pd (srcDen[(j+1)*3+1], srcDen[j*3+1]);
- __m128d sdenz = _mm_set_pd (srcDen[(j+1)*3+2], srcDen[j*3+2]);
- __m128d dX, dY, dZ;
- __m128d dR2;
- __m128d S;
- __m128d S2;
- __m128d S3;
- dX = _mm_sub_pd(txi , sxj);
- dY = _mm_sub_pd(tyi , syj);
- dZ = _mm_sub_pd(tzi , szj);
- sxj = _mm_mul_pd(dX, dX);
- syj = _mm_mul_pd(dY, dY);
- szj = _mm_mul_pd(dZ, dZ);
- dR2 = _mm_add_pd(sxj, syj);
- dR2 = _mm_add_pd(szj, dR2);
- __m128d temp = _mm_cmpeq_pd (dR2, zero);
- __m128d xhalf = _mm_mul_pd (half, dR2);
- __m128 dR2_s = _mm_cvtpd_ps(dR2);
- __m128 S_s = _mm_rsqrt_ps(dR2_s);
- __m128d S_d = _mm_cvtps_pd(S_s);
- // To handle the condition when src and trg coincide
- S_d = _mm_andnot_pd (temp, S_d);
- S = _mm_mul_pd (S_d, S_d);
- S = _mm_mul_pd (S, xhalf);
- S = _mm_sub_pd (opf, S);
- S = _mm_mul_pd (S, S_d);
- S2 = _mm_mul_pd (S, S);
- S3 = _mm_mul_pd (S2, S);
- __m128d dotx = _mm_mul_pd (dX, sdenx);
- __m128d doty = _mm_mul_pd (dY, sdeny);
- __m128d dotz = _mm_mul_pd (dZ, sdenz);
- __m128d dot_sum = _mm_add_pd (dotx, doty);
- dot_sum = _mm_add_pd (dot_sum, dotz);
- dot_sum = _mm_mul_pd (dot_sum, S2);
- 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);
- 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);
- 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);
- 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);
- 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);
- 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);
- 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);
- 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);
- 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);
- }
- tempxx = _mm_mul_pd (tempxx, ooepmu);
- tempxy = _mm_mul_pd (tempxy, ooepmu);
- tempxz = _mm_mul_pd (tempxz, ooepmu);
- tempyx = _mm_mul_pd (tempyx, ooepmu);
- tempyy = _mm_mul_pd (tempyy, ooepmu);
- tempyz = _mm_mul_pd (tempyz, ooepmu);
- tempzx = _mm_mul_pd (tempzx, ooepmu);
- tempzy = _mm_mul_pd (tempzy, ooepmu);
- tempzz = _mm_mul_pd (tempzz, ooepmu);
- _mm_store_pd(tempvalxx, tempxx); _mm_store_pd(tempvalxy, tempxy); _mm_store_pd(tempvalxz, tempxz);
- _mm_store_pd(tempvalyx, tempyx); _mm_store_pd(tempvalyy, tempyy); _mm_store_pd(tempvalyz, tempyz);
- _mm_store_pd(tempvalzx, tempzx); _mm_store_pd(tempvalzy, tempzy); _mm_store_pd(tempvalzz, tempzz);
- for (int k = 0; k < SIMD_LEN; k++) {
- trgVal[i*9 ] += tempvalxx[k];
- trgVal[i*9+1] += tempvalxy[k];
- trgVal[i*9+2] += tempvalxz[k];
- trgVal[i*9+3] += tempvalyx[k];
- trgVal[i*9+4] += tempvalyy[k];
- trgVal[i*9+5] += tempvalyz[k];
- trgVal[i*9+6] += tempvalzx[k];
- trgVal[i*9+7] += tempvalzy[k];
- trgVal[i*9+8] += tempvalzz[k];
- }
- for (; j < ns; j++) {
- double x = tx[i] - sx[j];
- double y = ty[i] - sy[j];
- double z = tz[i] - sz[j];
- double r2 = x*x + y*y + z*z;
- double r = sqrt(r2);
- double invdr;
- if (r == 0)
- invdr = 0;
- else
- invdr = 1/r;
- double invdr2=invdr*invdr;
- double invdr3=invdr2*invdr;
- double dot = (x*srcDen[j*3] + y*srcDen[j*3+1] + z*srcDen[j*3+2]);
- trgVal[i*9 ] += OOEP*oomeu*invdr3*( x*srcDen[j*3 ] - srcDen[j*3 ]*x + dot*(1-3*x*x*invdr2) );
- trgVal[i*9+1] += OOEP*oomeu*invdr3*( y*srcDen[j*3 ] - srcDen[j*3+1]*x + dot*(0-3*y*x*invdr2) );
- trgVal[i*9+2] += OOEP*oomeu*invdr3*( z*srcDen[j*3 ] - srcDen[j*3+2]*x + dot*(0-3*z*x*invdr2) );
- trgVal[i*9+3] += OOEP*oomeu*invdr3*( x*srcDen[j*3+1] - srcDen[j*3 ]*y + dot*(0-3*x*y*invdr2) );
- trgVal[i*9+4] += OOEP*oomeu*invdr3*( y*srcDen[j*3+1] - srcDen[j*3+1]*y + dot*(1-3*y*y*invdr2) );
- trgVal[i*9+5] += OOEP*oomeu*invdr3*( z*srcDen[j*3+1] - srcDen[j*3+2]*y + dot*(0-3*z*y*invdr2) );
- trgVal[i*9+6] += OOEP*oomeu*invdr3*( x*srcDen[j*3+2] - srcDen[j*3 ]*z + dot*(0-3*x*z*invdr2) );
- trgVal[i*9+7] += OOEP*oomeu*invdr3*( y*srcDen[j*3+2] - srcDen[j*3+1]*z + dot*(0-3*y*z*invdr2) );
- trgVal[i*9+8] += OOEP*oomeu*invdr3*( z*srcDen[j*3+2] - srcDen[j*3+2]*z + dot*(1-3*z*z*invdr2) );
- }
- }
- return;
- }
- #undef SIMD_LEN
- #define X(s,k) (s)[(k)*COORD_DIM]
- #define Y(s,k) (s)[(k)*COORD_DIM+1]
- #define Z(s,k) (s)[(k)*COORD_DIM+2]
- 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)
- {
- std::vector<double> xs(ns+1); std::vector<double> xt(nt);
- std::vector<double> ys(ns+1); std::vector<double> yt(nt);
- std::vector<double> zs(ns+1); std::vector<double> zt(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- stokesDirectVecSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot,kernel_coef);
- return;
- }
- void stokesPressureSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- std::vector<double> xs(ns+1); std::vector<double> xt(nt);
- std::vector<double> ys(ns+1); std::vector<double> yt(nt);
- std::vector<double> zs(ns+1); std::vector<double> zt(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- stokesPressureSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
- return;
- }
- void stokesStressSSEShuffle(const int ns, const int nt, double const src[], double const trg[], double const den[], double pot[], mem::MemoryManager* mem_mgr=NULL)
- {
- std::vector<double> xs(ns+1); std::vector<double> xt(nt);
- std::vector<double> ys(ns+1); std::vector<double> yt(nt);
- std::vector<double> zs(ns+1); std::vector<double> zt(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- stokesStressSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot);
- return;
- }
- 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)
- {
- std::vector<double> xs(ns+1); std::vector<double> xt(nt);
- std::vector<double> ys(ns+1); std::vector<double> yt(nt);
- std::vector<double> zs(ns+1); std::vector<double> zt(nt);
- int x_shift = size_t(&xs[0]) % IDEAL_ALIGNMENT ? 1:0;
- int y_shift = size_t(&ys[0]) % IDEAL_ALIGNMENT ? 1:0;
- int z_shift = size_t(&zs[0]) % IDEAL_ALIGNMENT ? 1:0;
- //1. reshuffle memory
- for (int k =0;k<ns;k++){
- xs[k+x_shift]=X(src,k);
- ys[k+y_shift]=Y(src,k);
- zs[k+z_shift]=Z(src,k);
- }
- for (int k=0;k<nt;k++){
- xt[k]=X(trg,k);
- yt[k]=Y(trg,k);
- zt[k]=Z(trg,k);
- }
- //2. perform caclulation
- stokesGradSSE(ns,nt,&xs[x_shift],&ys[y_shift],&zs[z_shift],&xt[0],&yt[0],&zt[0],den,pot,kernel_coef);
- return;
- }
- #undef X
- #undef Y
- #undef Z
- #undef IDEAL_ALIGNMENT
- #undef DECL_SIMD_ALIGNED
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(28*dof));
- const T mu=1.0;
- stokesDirectSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mu, mem_mgr);
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(17*dof));
- stokesPressureSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mem_mgr);
- return;
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(45*dof));
- stokesStressSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mem_mgr);
- }
- template <>
- 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){
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(89*dof));
- const T mu=1.0;
- stokesGradSSEShuffle(src_cnt, trg_cnt, r_src, r_trg, v_src_, k_out, mu, mem_mgr);
- }
- #endif
- #endif
- ////////////////////////////////////////////////////////////////////////////////
- //////// BIOT-SAVART KERNEL ////////
- ////////////////////////////////////////////////////////////////////////////////
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(26*dof));
- #endif
- const T OOFP = -1.0/(4.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[3]={0,0,0};
- for(int s=0;s<src_cnt;s++){
- T dR[3]={r_trg[3*t ]-r_src[3*s ],
- r_trg[3*t+1]-r_src[3*s+1],
- r_trg[3*t+2]-r_src[3*s+2]};
- T R2 = (dR[0]*dR[0]+dR[1]*dR[1]+dR[2]*dR[2]);
- if (R2!=0){
- T invR2=1.0/R2;
- T invR=sqrt(invR2);
- T invR3=invR*invR2;
- T v_src[3]={v_src_[(s*dof+i)*3 ],
- v_src_[(s*dof+i)*3+1],
- v_src_[(s*dof+i)*3+2]};
- p[0] -= (v_src[1]*dR[2]-v_src[2]*dR[1])*invR3;
- p[1] -= (v_src[2]*dR[0]-v_src[0]*dR[2])*invR3;
- p[2] -= (v_src[0]*dR[1]-v_src[1]*dR[0])*invR3;
- }
- }
- k_out[(t*dof+i)*3+0] += p[0]*OOFP;
- k_out[(t*dof+i)*3+1] += p[1]*OOFP;
- k_out[(t*dof+i)*3+2] += p[2]*OOFP;
- }
- }
- }
- ////////////////////////////////////////////////////////////////////////////////
- //////// HELMHOLTZ KERNEL ////////
- ////////////////////////////////////////////////////////////////////////////////
- /**
- * \brief Green's function for the Helmholtz's equation. Kernel tensor
- * dimension = 2x2.
- */
- template <class T>
- 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){
- #ifndef __MIC__
- Profile::Add_FLOP((long long)trg_cnt*(long long)src_cnt*(24*dof));
- #endif
- const T mu = (20.0*const_pi<T>());
- for(int t=0;t<trg_cnt;t++){
- for(int i=0;i<dof;i++){
- T p[2]={0,0};
- for(int s=0;s<src_cnt;s++){
- T dX_reg=r_trg[3*t ]-r_src[3*s ];
- T dY_reg=r_trg[3*t+1]-r_src[3*s+1];
- T dZ_reg=r_trg[3*t+2]-r_src[3*s+2];
- T R = (dX_reg*dX_reg+dY_reg*dY_reg+dZ_reg*dZ_reg);
- if (R!=0){
- R = sqrt(R);
- T invR=1.0/R;
- T G[2]={cos(mu*R)*invR, sin(mu*R)*invR};
- p[0] += v_src[(s*dof+i)*2+0]*G[0] - v_src[(s*dof+i)*2+1]*G[1];
- p[1] += v_src[(s*dof+i)*2+0]*G[1] + v_src[(s*dof+i)*2+1]*G[0];
- }
- }
- k_out[(t*dof+i)*2+0] += p[0];
- k_out[(t*dof+i)*2+1] += p[1];
- }
- }
- }
- template <class T>
- 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){
- //TODO Implement this.
- }
- }//end namespace
|