#pragma once #include #include // for unlink #include #include #include #include namespace particle{ //! implement Marcos et al. PLT calculation inline void test_plt( void ){ csoca::ilog << "-------------------------------------------------------------------------------" << std::endl; csoca::ilog << "Testing PLT implementation..." << std::endl; real_t boxlen = 1.0; size_t ngrid = 64; size_t npgrid = 1; size_t dpg = ngrid/npgrid; size_t nump = npgrid*npgrid*npgrid; real_t pweight = 1.0/real_t(nump); real_t eta = 2.0 * boxlen/ngrid; const real_t alpha = 1.0/std::sqrt(2)/eta; const real_t alpha2 = alpha*alpha; const real_t alpha3 = alpha2*alpha; const real_t sqrtpi = std::sqrt(M_PI); const real_t pi3halfs = std::pow(M_PI,1.5); const real_t dV( std::pow( boxlen/ngrid, 3 ) ); Grid_FFT rho({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); std::vector< vec3 > gpos ; auto kronecker = []( int i, int j ) -> real_t { return (i==j)? 1.0 : 0.0; }; auto greensftide_sr = [&]( int mu, int nu, const vec3& vR, const vec3& vP ) -> real_t { auto d = vR-vP; d.x = (d.x>0.5)? d.x-1.0 : (d.x<-0.5)? d.x+1.0 : d.x; d.y = (d.y>0.5)? d.y-1.0 : (d.y<-0.5)? d.y+1.0 : d.y; d.z = (d.z>0.5)? d.z-1.0 : (d.z<-0.5)? d.z+1.0 : d.z; auto r = d.norm(); if( r< 1e-14 ) return 0.0; real_t val = 0.0; val -= d[mu]*d[nu]/(r*r) * alpha3/pi3halfs * std::exp(-alpha*alpha*r*r); val += 1.0/(4.0*M_PI)*(kronecker(mu,nu)/std::pow(r,3) - 3.0 * (d[mu]*d[nu])/std::pow(r,5)) * (std::erfc(alpha*r) + 2.0*alpha/sqrtpi*std::exp(-alpha*alpha*r*r)*r); return pweight * val; }; gpos.reserve(nump); // sc for( size_t i=0; i ccomplex_t { real_t kmod = k.norm(); return -x * std::exp(-0.5*eta*eta*kmod*kmod) / (kmod*kmod); }); rho.zero_DC_mode(); auto evaluate_D = [&]( int mu, int nu, const vec3& v ) -> real_t{ real_t sr = 0.0; for( auto& p : gpos ){ sr += greensftide_sr( mu, nu, v, p); } if( v.norm()<1e-14 ) return 0.0; return sr; }; // std::random_device rd; //Will be used to obtain a seed for the random number engine // std::mt19937 gen(rd()); //Standard mersenne_twister_engine seeded with rd() // std::uniform_real_distribution<> dis(-0.25,0.25); Grid_FFT D_xx({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); Grid_FFT D_xy({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); Grid_FFT D_xz({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); Grid_FFT D_yy({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); Grid_FFT D_yz({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); Grid_FFT D_zz({ngrid, ngrid, ngrid}, {boxlen, boxlen, boxlen}); #pragma omp parallel for for( size_t i=0; i p; p.x = real_t(i)/ngrid; for( size_t j=0; j D; vec3 eval, evec1, evec2, evec3; for( size_t j=0; j kv = D_xx.get_k(i,j,k); D = { std::real(D_xx.kelem(i,j,k) - kv[0]*kv[0] * rho.kelem(i,j,k) ), std::real(D_xy.kelem(i,j,k) - kv[0]*kv[1] * rho.kelem(i,j,k) ), std::real(D_xz.kelem(i,j,k) - kv[0]*kv[2] * rho.kelem(i,j,k) ), std::real(D_yy.kelem(i,j,k) - kv[1]*kv[1] * rho.kelem(i,j,k) ), std::real(D_yz.kelem(i,j,k) - kv[1]*kv[2] * rho.kelem(i,j,k) ), std::real(D_zz.kelem(i,j,k) - kv[2]*kv[2] * rho.kelem(i,j,k) ) }; D.eigen(eval, evec1, evec2, evec3); ofs << std::setw(16) << kv.norm() / kNyquist << std::setw(16) << eval[0] *nfac + 1.0/3.0 << std::setw(16) << eval[1] *nfac + 1.0/3.0 << std::setw(16) << eval[2] *nfac + 1.0/3.0 << std::setw(16) << kv[0] << std::setw(16) << kv[1] << std::setw(16) << kv[2] << std::endl; } } } // std::string filename("plt_test.hdf5"); // unlink(filename.c_str()); // #if defined(USE_MPI) // MPI_Barrier(MPI_COMM_WORLD); // #endif // rho.Write_to_HDF5(filename, "rho"); } }