1
0
Fork 0
mirror of https://github.com/cosmo-sims/MUSIC.git synced 2024-09-19 17:03:46 +02:00

new version number 1.5

fixed a compile error in tetmesh output
made k-space sampling the default
This commit is contained in:
Oliver Hahn 2013-11-14 11:04:08 +01:00
parent 144efb78ea
commit caec0b9d76
4 changed files with 19 additions and 19 deletions

View file

@ -393,23 +393,23 @@ namespace convolution{
cparam_.lz = dx * cparam_.nz;
cparam_.pcf = pcf_;
fftw_real *rkernel = reinterpret_cast<fftw_real*>( &kdata_[0] );
fftw_real *rkernel = reinterpret_cast<fftw_real*>( &kdata_[0] );
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_complex *kkernel = reinterpret_cast<fftwf_complex*> (&rkernel[0]);
fftwf_complex *kkernel = reinterpret_cast<fftwf_complex*> (&rkernel[0]);
fftwf_plan plan = fftwf_plan_dft_r2c_3d( cparam_.nx, cparam_.ny, cparam_.nz, rkernel, kkernel, FFTW_ESTIMATE);
fftwf_execute(plan);
fftwf_destroy_plan(plan);
#else
fftw_complex *kkernel = reinterpret_cast<fftw_complex*> (&rkernel[0]);
fftw_complex *kkernel = reinterpret_cast<fftw_complex*> (&rkernel[0]);
fftw_plan plan = fftw_plan_dft_r2c_3d( cparam_.nx, cparam_.ny, cparam_.nz, rkernel, kkernel, FFTW_ESTIMATE);
fftw_execute(plan);
fftw_destroy_plan(plan);
#endif
#else
rfftwnd_plan plan = rfftw3d_create_plan( cparam_.nx, cparam_.ny, cparam_.nz,
FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE);
rfftwnd_plan plan = rfftw3d_create_plan( cparam_.nx, cparam_.ny, cparam_.nz,
FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE);
#ifndef SINGLETHREAD_FFTW
@ -562,9 +562,7 @@ namespace convolution{
int ref_fac = (deconv&&kspacepoisson)? 2 : 0;
const int ql = -ref_fac/2+1, qr = ql+ref_fac;
const double rf8 = pow(ref_fac,3);
const double dx05 = 0.5*dx, dx025 = 0.25*dx;
std::cerr << ">>>>>>>>>>>> " << ref_fac << " <<<<<<<<<<<<<<<<" << std::endl;
const double dx05 = 0.5*dx, dx025 = 0.25*dx;
#endif
if( bperiodic )
@ -1134,8 +1132,7 @@ namespace convolution{
}
#endif // OLD_KERNEL_SAMPLING
#endif // #OLD_KERNEL_SAMPLING
sprintf(cachefname,"temp_kernel_level%03d.tmp",ilevel);
LOGUSER("Storing kernel in temp file \'%s\'.",cachefname);
fp = fopen(cachefname,"w+");

View file

@ -160,9 +160,9 @@ void fft_interpolate( m1& V, m2& v, bool from_basegrid=false )
{
int ii(i),jj(j),kk(k);
if( i> nxc/2 ) ii += nxf/2;
if( j> nyc/2 ) jj += nyf/2;
if( k> nzc/2 ) kk += nzf/2;
if( i> (int)nxc/2 ) ii += (int)nxf/2;
if( j> (int)nyc/2 ) jj += (int)nyf/2;
if( k> (int)nzc/2 ) kk += (int)nzf/2;
size_t qc,qf;
@ -347,7 +347,7 @@ void GenerateDensityUnigrid( config_file& cf, transfer_function *ptf, tf_type ty
levelmin = cf.getValueSafe<unsigned>("setup","levelmin_TF",levelminPoisson);
levelmax = cf.getValue<unsigned>("setup","levelmax");
bool kspace = cf.getValueSafe<unsigned>("setup","kspace_TF",false);
bool kspace = cf.getValue<bool>("setup","kspace_TF");
unsigned nbase = 1<<levelmin;
@ -438,9 +438,9 @@ void GenerateDensityHierarchy( config_file& cf, transfer_function *ptf, tf_type
levelminPoisson = cf.getValue<unsigned>("setup","levelmin");
levelmin = cf.getValueSafe<unsigned>("setup","levelmin_TF",levelminPoisson);
levelmax = cf.getValue<unsigned>("setup","levelmax");
kspaceTF = cf.getValueSafe<bool>("setup", "kspace_TF", false);
kspaceTF = cf.getValue<bool>("setup", "kspace_TF");
blend_sharpness = cf.getValueSafe<double>("setup","kspace_filter",0.333);
blend_sharpness = cf.getValueSafe<double>("setup","kspace_filter", blend_sharpness);
unsigned nbase = 1<<levelmin;

View file

@ -37,7 +37,7 @@
#include "transfer_function.hh"
#define THE_CODE_NAME "music!"
#define THE_CODE_VERSION "1.4b"
#define THE_CODE_VERSION "1.5"
namespace music
@ -359,7 +359,10 @@ int main (int argc, const char * argv[])
cf.insertValue("setup","levelmin_TF",cf.getValue<std::string>("setup","levelmin"));
}
if( cf.getValueSafe<bool>( "setup", "kspace_TF", true ) )
LOGINFO("Using k-space sampled transfer functions...");
else
LOGINFO("Using real space sampled transfer functions...");
//------------------------------------------------------------------------------
//... initialize multithread FFTW

View file

@ -1110,7 +1110,7 @@ public:
}
if( gh.is_in_mask(ilevel,i,j,k) && !gh.is_refined(ilevel,i,j,k) )
if( gh.is_in_mask(ilevel,ix,iy,iz) && !gh.is_refined(ilevel,ix,iy,iz) )
{
dorefine = false;
break;