mirror of
https://github.com/cosmo-sims/MUSIC.git
synced 2024-09-19 17:03:46 +02:00
* Added support for velocity potentials for baryons and CDM
* Added new transfer function type 'music' which reads in a file with density and velocity transfer functions for baryons and CDM. * Powerspectrum evolution in two-fluid sim appears to scale as predicted now at the per cent level * Tool to compute velocity transfer functions is not included as of now. Let's see.. * Velocity potentials are not yet supported with 2LPT.
This commit is contained in:
parent
37c0d4be36
commit
f3364e9877
8 changed files with 478 additions and 120 deletions
|
@ -978,15 +978,19 @@ namespace convolution{
|
||||||
/*************************************************************************************/
|
/*************************************************************************************/
|
||||||
/*************************************************************************************/
|
/*************************************************************************************/
|
||||||
/******* perform deconvolution *******************************************************/
|
/******* perform deconvolution *******************************************************/
|
||||||
|
bool deconv = pcf_->getValueSafe<bool>("setup","deconvolve",true);
|
||||||
if( pcf_->getValueSafe<bool>("setup","deconvolve",true) )
|
bool deconv_baryons = pcf_->getValueSafe<bool>("setup","deconvolve_baryons",false);
|
||||||
|
bool bsmooth_baryons = false;//type==baryon;// && !deconv_baryons;
|
||||||
|
bool baryons = type==baryon||type==vbaryon;
|
||||||
|
//if( deconv )
|
||||||
|
if( deconv && !(type==baryon&&!deconv_baryons) )
|
||||||
{
|
{
|
||||||
|
|
||||||
LOGUSER("Deconvolving fine kernel...");
|
LOGUSER("Deconvolving fine kernel...");
|
||||||
std::cout << " - Deconvoling density kernel...\n";
|
std::cout << " - Deconvoling density kernel...\n";
|
||||||
|
|
||||||
bool kspacepoisson = (pcf_->getValueSafe<bool>("poisson","fft_fine",true)|
|
bool kspacepoisson = ((pcf_->getValueSafe<bool>("poisson","fft_fine",true)|
|
||||||
pcf_->getValueSafe<bool>("poisson","kspace",false));
|
pcf_->getValueSafe<bool>("poisson","kspace",false)))&!baryons ;
|
||||||
|
|
||||||
double fftnorm = 1.0/((size_t)nx*(size_t)ny*(size_t)nz);
|
double fftnorm = 1.0/((size_t)nx*(size_t)ny*(size_t)nz);
|
||||||
double k0 = rkernel[0];
|
double k0 = rkernel[0];
|
||||||
|
@ -994,7 +998,8 @@ namespace convolution{
|
||||||
fftw_complex *kkernel = reinterpret_cast<fftw_complex*>( &rkernel[0] );
|
fftw_complex *kkernel = reinterpret_cast<fftw_complex*>( &rkernel[0] );
|
||||||
|
|
||||||
//... subtract white noise component before deconvolution
|
//... subtract white noise component before deconvolution
|
||||||
rkernel[0] = 0.0;
|
if(!bsmooth_baryons)
|
||||||
|
rkernel[0] = 0.0;
|
||||||
|
|
||||||
#ifdef FFTW3
|
#ifdef FFTW3
|
||||||
fftw_plan
|
fftw_plan
|
||||||
|
@ -1036,7 +1041,7 @@ namespace convolution{
|
||||||
double kkmax = kmax;
|
double kkmax = kmax;
|
||||||
size_t q = ((size_t)i*ny+(size_t)j)*(nz/2+1)+(size_t)k;
|
size_t q = ((size_t)i*ny+(size_t)j)*(nz/2+1)+(size_t)k;
|
||||||
|
|
||||||
if( true )//!cparam_.smooth )
|
if( !bsmooth_baryons )
|
||||||
{
|
{
|
||||||
if( kspacepoisson )
|
if( kspacepoisson )
|
||||||
{
|
{
|
||||||
|
@ -1070,9 +1075,14 @@ namespace convolution{
|
||||||
kkernel[q].im *= ipix;
|
kkernel[q].im *= ipix;
|
||||||
#endif
|
#endif
|
||||||
}
|
}
|
||||||
}else{
|
}
|
||||||
|
#if 1
|
||||||
|
else{
|
||||||
//... if smooth==true, convolve with
|
//... if smooth==true, convolve with
|
||||||
//... NGP kernel to get CIC smoothness
|
//... NGP kernel to get CIC smoothness
|
||||||
|
|
||||||
|
//kkmax *= 2.0;
|
||||||
|
|
||||||
double ipix = 1.0;
|
double ipix = 1.0;
|
||||||
if( i > 0 )
|
if( i > 0 )
|
||||||
ipix /= sin(kx*2.0*kkmax)/(kx*2.0*kkmax);
|
ipix /= sin(kx*2.0*kkmax)/(kx*2.0*kkmax);
|
||||||
|
@ -1089,7 +1099,7 @@ namespace convolution{
|
||||||
kkernel[q].im /= ipix;
|
kkernel[q].im /= ipix;
|
||||||
#endif
|
#endif
|
||||||
}
|
}
|
||||||
|
#endif
|
||||||
//... store k-space average
|
//... store k-space average
|
||||||
#ifdef FFTW3
|
#ifdef FFTW3
|
||||||
if( k==0 || k==nz/2 )
|
if( k==0 || k==nz/2 )
|
||||||
|
@ -1118,7 +1128,8 @@ namespace convolution{
|
||||||
dk = k0-ksum/kcount;
|
dk = k0-ksum/kcount;
|
||||||
|
|
||||||
//... set white noise component to zero if smoothing is enabled
|
//... set white noise component to zero if smoothing is enabled
|
||||||
if( false )//cparam_.smooth )
|
//if( false )//cparam_.smooth )
|
||||||
|
if( bsmooth_baryons )
|
||||||
dk = 0.0;
|
dk = 0.0;
|
||||||
|
|
||||||
//... enforce the r=0 component by adjusting the k-space mean
|
//... enforce the r=0 component by adjusting the k-space mean
|
||||||
|
|
|
@ -1,60 +1,59 @@
|
||||||
[setup]
|
[setup]
|
||||||
boxlength = 100
|
boxlength = 100
|
||||||
zstart = 50
|
zstart = 50
|
||||||
levelmin = 7
|
levelmin = 7
|
||||||
levelmin_TF = 7
|
levelmin_TF = 7
|
||||||
levelmax = 8
|
levelmax = 8
|
||||||
padding = 8
|
padding = 8
|
||||||
overlap = 4
|
overlap = 4
|
||||||
ref_offset = 0.4, 0.4, 0.4
|
ref_offset = 0.4, 0.4, 0.4
|
||||||
ref_extent = 0.2, 0.2, 0.2
|
ref_extent = 0.2, 0.2, 0.2
|
||||||
align_top = yes
|
align_top = yes
|
||||||
baryons = no
|
baryons = no
|
||||||
use_2LPT = no
|
use_2LPT = no
|
||||||
use_LLA = no
|
use_LLA = no
|
||||||
periodic_TF = yes
|
periodic_TF = yes
|
||||||
|
|
||||||
[cosmology]
|
[cosmology]
|
||||||
Omega_m = 0.276
|
Omega_m = 0.276
|
||||||
Omega_L = 0.724
|
Omega_L = 0.724
|
||||||
Omega_b = 0.045
|
Omega_b = 0.045
|
||||||
H0 = 70.3
|
H0 = 70.3
|
||||||
sigma_8 = 0.811
|
sigma_8 = 0.811
|
||||||
nspec = 0.961
|
nspec = 0.961
|
||||||
transfer = eisenstein ##this cannot be changed yet!
|
transfer = eisenstein
|
||||||
|
|
||||||
[random]
|
[random]
|
||||||
seed[7] = 12345
|
seed[7] = 12345
|
||||||
seed[8] = 23456
|
seed[8] = 23456
|
||||||
disk_cached = no
|
|
||||||
|
|
||||||
[output]
|
[output]
|
||||||
##generic FROLIC data format (used for testing)
|
##generic FROLIC data format (used for testing)
|
||||||
#format = generic
|
#format = generic
|
||||||
#filename = debug.hdf5
|
#filename = debug.hdf5
|
||||||
|
|
||||||
##ENZO - also outputs the settings for the parameter file
|
##ENZO - also outputs the settings for the parameter file
|
||||||
#format = enzo
|
#format = enzo
|
||||||
#filename = ic.enzo
|
#filename = ic.enzo
|
||||||
|
|
||||||
##Gadget-2 (type=1: high-res particles, type=5: rest)
|
##Gadget-2 (type=1: high-res particles, type=5: rest)
|
||||||
##no gas possible at the moment
|
##no gas possible at the moment
|
||||||
format = gadget2
|
format = gadget2
|
||||||
filename = ics_gadget.dat
|
filename = ics_gadget.dat
|
||||||
shift_back = yes
|
shift_back = yes
|
||||||
|
|
||||||
##Grafic2 compatible format for use with RAMSES
|
##Grafic2 compatible format for use with RAMSES
|
||||||
##option 'ramses_nml'=yes writes out a startup nml file
|
##option 'ramses_nml'=yes writes out a startup nml file
|
||||||
#format = grafic2
|
#format = grafic2
|
||||||
#filename = ics_ramses
|
#filename = ics_ramses
|
||||||
#ramses_nml = yes
|
#ramses_nml = yes
|
||||||
|
|
||||||
[poisson]
|
[poisson]
|
||||||
fft_fine = no
|
fft_fine = no
|
||||||
accuracy = 1e-5
|
accuracy = 1e-5
|
||||||
pre_smooth = 3
|
pre_smooth = 3
|
||||||
post_smooth = 3
|
post_smooth = 3
|
||||||
smoother = gs
|
smoother = gs
|
||||||
laplace_order = 6
|
laplace_order = 6
|
||||||
grad_order = 6
|
grad_order = 6
|
||||||
|
|
||||||
|
|
131
main.cc
131
main.cc
|
@ -438,6 +438,10 @@ int main (int argc, const char * argv[])
|
||||||
LOGUSER("Computing dark matter displacements...");
|
LOGUSER("Computing dark matter displacements...");
|
||||||
|
|
||||||
grid_hierarchy f( nbnd ), u(nbnd);
|
grid_hierarchy f( nbnd ), u(nbnd);
|
||||||
|
tf_type my_tf_type = cdm;
|
||||||
|
if( !do_baryons || !the_transfer_function_plugin->tf_is_distinct() )
|
||||||
|
my_tf_type = total;
|
||||||
|
|
||||||
|
|
||||||
GenerateDensityHierarchy( cf, the_transfer_function_plugin, cdm , rh_TF, rand, f, true, false );
|
GenerateDensityHierarchy( cf, the_transfer_function_plugin, cdm , rh_TF, rand, f, true, false );
|
||||||
coarsen_density(rh_Poisson, f);
|
coarsen_density(rh_Poisson, f);
|
||||||
|
@ -507,45 +511,122 @@ int main (int argc, const char * argv[])
|
||||||
LOGUSER("Computing velocitites...");
|
LOGUSER("Computing velocitites...");
|
||||||
|
|
||||||
//... velocities
|
//... velocities
|
||||||
if( do_baryons )
|
if( !the_transfer_function_plugin->tf_has_velocities() || !do_baryons )
|
||||||
|
{
|
||||||
|
if( do_baryons )
|
||||||
|
{
|
||||||
|
GenerateDensityHierarchy( cf, the_transfer_function_plugin, total , rh_TF, rand, f, true, false );
|
||||||
|
coarsen_density(rh_Poisson, f);
|
||||||
|
normalize_density(f);
|
||||||
|
|
||||||
|
u = f; u.zero();
|
||||||
|
|
||||||
|
err = the_poisson_solver->solve(f, u);
|
||||||
|
|
||||||
|
if(!bdefd)
|
||||||
|
f.deallocate();
|
||||||
|
}
|
||||||
|
grid_hierarchy data_forIO(u);
|
||||||
|
for( int icoord = 0; icoord < 3; ++icoord )
|
||||||
|
{
|
||||||
|
//... displacement
|
||||||
|
if(bdefd)
|
||||||
|
{
|
||||||
|
data_forIO.zero();
|
||||||
|
*data_forIO.get_grid(data_forIO.levelmax()) = *f.get_grid(f.levelmax());
|
||||||
|
poisson_hybrid(*data_forIO.get_grid(data_forIO.levelmax()), icoord, grad_order, data_forIO.levelmin()==data_forIO.levelmax());
|
||||||
|
*data_forIO.get_grid(data_forIO.levelmax()) /= 1<<f.levelmax();
|
||||||
|
the_poisson_solver->gradient_add(icoord, u, data_forIO );
|
||||||
|
}
|
||||||
|
else
|
||||||
|
the_poisson_solver->gradient(icoord, u, data_forIO );
|
||||||
|
|
||||||
|
//... multiply to get velocity
|
||||||
|
data_forIO *= cosmo.vfact;
|
||||||
|
|
||||||
|
if(do_CVM)
|
||||||
|
subtract_finest_mean(data_forIO);
|
||||||
|
|
||||||
|
the_output_plugin->write_dm_velocity(icoord, data_forIO);
|
||||||
|
if( do_baryons )
|
||||||
|
the_output_plugin->write_gas_velocity(icoord, data_forIO);
|
||||||
|
}
|
||||||
|
|
||||||
|
}
|
||||||
|
else
|
||||||
{
|
{
|
||||||
GenerateDensityHierarchy( cf, the_transfer_function_plugin, total , rh_TF, rand, f, true, false );
|
//... we do baryons and have velocity transfer functions
|
||||||
|
//... do DM first
|
||||||
|
GenerateDensityHierarchy( cf, the_transfer_function_plugin, vcdm , rh_TF, rand, f, true, false );
|
||||||
|
//GenerateDensityHierarchy( cf, the_transfer_function_plugin, vcdm , rh_TF, rand, f, false, false );
|
||||||
coarsen_density(rh_Poisson, f);
|
coarsen_density(rh_Poisson, f);
|
||||||
normalize_density(f);
|
normalize_density(f);
|
||||||
|
|
||||||
u = f; u.zero();
|
u = f; u.zero();
|
||||||
|
|
||||||
err = the_poisson_solver->solve(f, u);
|
err = the_poisson_solver->solve(f, u);
|
||||||
|
|
||||||
if(!bdefd)
|
if(!bdefd)
|
||||||
f.deallocate();
|
f.deallocate();
|
||||||
}
|
|
||||||
grid_hierarchy data_forIO(u);
|
grid_hierarchy data_forIO(u);
|
||||||
for( int icoord = 0; icoord < 3; ++icoord )
|
for( int icoord = 0; icoord < 3; ++icoord )
|
||||||
{
|
|
||||||
//... displacement
|
|
||||||
if(bdefd)
|
|
||||||
{
|
{
|
||||||
data_forIO.zero();
|
//... displacement
|
||||||
*data_forIO.get_grid(data_forIO.levelmax()) = *f.get_grid(f.levelmax());
|
if(bdefd)
|
||||||
poisson_hybrid(*data_forIO.get_grid(data_forIO.levelmax()), icoord, grad_order, data_forIO.levelmin()==data_forIO.levelmax());
|
{
|
||||||
*data_forIO.get_grid(data_forIO.levelmax()) /= 1<<f.levelmax();
|
data_forIO.zero();
|
||||||
the_poisson_solver->gradient_add(icoord, u, data_forIO );
|
*data_forIO.get_grid(data_forIO.levelmax()) = *f.get_grid(f.levelmax());
|
||||||
|
poisson_hybrid(*data_forIO.get_grid(data_forIO.levelmax()), icoord, grad_order, data_forIO.levelmin()==data_forIO.levelmax());
|
||||||
|
*data_forIO.get_grid(data_forIO.levelmax()) /= 1<<f.levelmax();
|
||||||
|
the_poisson_solver->gradient_add(icoord, u, data_forIO );
|
||||||
|
}
|
||||||
|
else
|
||||||
|
the_poisson_solver->gradient(icoord, u, data_forIO );
|
||||||
|
|
||||||
|
//... multiply to get velocity
|
||||||
|
data_forIO *= cosmo.vfact;
|
||||||
|
|
||||||
|
//... we have two velocity contributions, can't do averaging at the moment
|
||||||
|
//if(do_CVM)
|
||||||
|
// subtract_finest_mean(data_forIO);
|
||||||
|
|
||||||
|
the_output_plugin->write_dm_velocity(icoord, data_forIO);
|
||||||
}
|
}
|
||||||
else
|
data_forIO.deallocate();
|
||||||
|
|
||||||
|
//... do baryons
|
||||||
|
GenerateDensityHierarchy( cf, the_transfer_function_plugin, vbaryon , rh_TF, rand, f, false, false );
|
||||||
|
coarsen_density(rh_Poisson, f);
|
||||||
|
normalize_density(f);
|
||||||
|
|
||||||
|
u = f; u.zero();
|
||||||
|
|
||||||
|
err = the_poisson_solver->solve(f, u);
|
||||||
|
|
||||||
|
//if(!bdefd)
|
||||||
|
f.deallocate();
|
||||||
|
|
||||||
|
data_forIO = u;
|
||||||
|
for( int icoord = 0; icoord < 3; ++icoord )
|
||||||
|
{
|
||||||
|
//... displacement
|
||||||
the_poisson_solver->gradient(icoord, u, data_forIO );
|
the_poisson_solver->gradient(icoord, u, data_forIO );
|
||||||
|
|
||||||
//... multiply to get velocity
|
//... multiply to get velocity
|
||||||
data_forIO *= cosmo.vfact;
|
data_forIO *= cosmo.vfact;
|
||||||
|
|
||||||
if(do_CVM)
|
//... we have two velocity contributions, can't do averaging at the moment
|
||||||
subtract_finest_mean(data_forIO);
|
//if(do_CVM)
|
||||||
|
// subtract_finest_mean(data_forIO);
|
||||||
the_output_plugin->write_dm_velocity(icoord, data_forIO);
|
|
||||||
if( do_baryons )
|
|
||||||
the_output_plugin->write_gas_velocity(icoord, data_forIO);
|
the_output_plugin->write_gas_velocity(icoord, data_forIO);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
/*********************************************************************************************/
|
||||||
|
/*********************************************************************************************/
|
||||||
|
/*** 2LPT ************************************************************************************/
|
||||||
|
/*********************************************************************************************/
|
||||||
}else {
|
}else {
|
||||||
//.. use 2LPT ...
|
//.. use 2LPT ...
|
||||||
LOGUSER("Entering 2LPT branch");
|
LOGUSER("Entering 2LPT branch");
|
||||||
|
|
|
@ -63,7 +63,7 @@ public:
|
||||||
m_Gamma = FreeGamma;
|
m_Gamma = FreeGamma;
|
||||||
|
|
||||||
tf_distinct_ = false;
|
tf_distinct_ = false;
|
||||||
|
tf_withvel_ = false;
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
|
@ -135,6 +135,7 @@ public:
|
||||||
gsl_spline_init (spline_baryon, &m_tab_k[0], &m_tab_Tk_baryon[0], m_tab_k.size() );
|
gsl_spline_init (spline_baryon, &m_tab_k[0], &m_tab_Tk_baryon[0], m_tab_k.size() );
|
||||||
|
|
||||||
tf_distinct_ = true;
|
tf_distinct_ = true;
|
||||||
|
tf_withvel_ = false;
|
||||||
}
|
}
|
||||||
|
|
||||||
~transfer_CAMB_plugin()
|
~transfer_CAMB_plugin()
|
||||||
|
|
|
@ -186,7 +186,10 @@ public:
|
||||||
TFset_parameters( (cosmo_.Omega_m)*cosmo_.H0*cosmo_.H0*(0.01*0.01),
|
TFset_parameters( (cosmo_.Omega_m)*cosmo_.H0*cosmo_.H0*(0.01*0.01),
|
||||||
cosmo_.Omega_b/(cosmo_.Omega_m-cosmo_.Omega_b),//-aCosm.Omega_b),
|
cosmo_.Omega_b/(cosmo_.Omega_m-cosmo_.Omega_b),//-aCosm.Omega_b),
|
||||||
Tcmb);
|
Tcmb);
|
||||||
|
|
||||||
|
std::cerr << "CHECK!!\n";
|
||||||
tf_distinct_ = false;
|
tf_distinct_ = false;
|
||||||
|
tf_withvel_ = false;
|
||||||
}
|
}
|
||||||
|
|
||||||
//! Computes the transfer function for k in Mpc/h by calling TFfit_onek
|
//! Computes the transfer function for k in Mpc/h by calling TFfit_onek
|
||||||
|
|
296
plugins/transfer_music.cc
Normal file
296
plugins/transfer_music.cc
Normal file
|
@ -0,0 +1,296 @@
|
||||||
|
/*
|
||||||
|
|
||||||
|
transfer_camb.cc - This file is part of MUSIC -
|
||||||
|
a code to generate multi-scale initial conditions
|
||||||
|
for cosmological simulations
|
||||||
|
|
||||||
|
Copyright (C) 2010 Oliver Hahn
|
||||||
|
|
||||||
|
This program is free software: you can redistribute it and/or modify
|
||||||
|
it under the terms of the GNU General Public License as published by
|
||||||
|
the Free Software Foundation, either version 3 of the License, or
|
||||||
|
(at your option) any later version.
|
||||||
|
|
||||||
|
This program is distributed in the hope that it will be useful,
|
||||||
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||||
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||||
|
GNU General Public License for more details.
|
||||||
|
|
||||||
|
You should have received a copy of the GNU General Public License
|
||||||
|
along with this program. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
|
||||||
|
*/
|
||||||
|
|
||||||
|
#include "transfer_function.hh"
|
||||||
|
|
||||||
|
class transfer_MUSIC_plugin : public transfer_function_plugin
|
||||||
|
{
|
||||||
|
|
||||||
|
private:
|
||||||
|
//Cosmology m_Cosmology;
|
||||||
|
|
||||||
|
std::string m_filename_Pk, m_filename_Tk;
|
||||||
|
std::vector<double> m_tab_k, m_tab_Tk_tot, m_tab_Tk_cdm, m_tab_Tk_baryon, m_tab_Tvk_cdm, m_tab_Tvk_baryon;
|
||||||
|
|
||||||
|
Spline_interp *m_psinterp;
|
||||||
|
gsl_interp_accel *acc_dtot, *acc_dcdm, *acc_dbaryon, *acc_vcdm, *acc_vbaryon;
|
||||||
|
gsl_spline *spline_dtot, *spline_dcdm, *spline_dbaryon, *spline_vcdm, *spline_vbaryon;
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
void read_table( void ){
|
||||||
|
#ifdef WITH_MPI
|
||||||
|
if( MPI::COMM_WORLD.Get_rank() == 0 ){
|
||||||
|
#endif
|
||||||
|
std::cerr
|
||||||
|
<< " - reading tabulated transfer function data from file \n"
|
||||||
|
<< " \'" << m_filename_Tk << "\'\n";
|
||||||
|
|
||||||
|
std::string line;
|
||||||
|
std::ifstream ifs( m_filename_Tk.c_str() );
|
||||||
|
|
||||||
|
if(! ifs.good() )
|
||||||
|
throw std::runtime_error("Could not find transfer function file \'"+m_filename_Tk+"\'");
|
||||||
|
|
||||||
|
m_tab_k.clear();
|
||||||
|
m_tab_Tk_tot.clear();
|
||||||
|
m_tab_Tk_cdm.clear();
|
||||||
|
m_tab_Tk_baryon.clear();
|
||||||
|
m_tab_Tvk_cdm.clear();
|
||||||
|
m_tab_Tvk_baryon.clear();
|
||||||
|
|
||||||
|
while( !ifs.eof() ){
|
||||||
|
getline(ifs,line);
|
||||||
|
|
||||||
|
if(ifs.eof()) break;
|
||||||
|
|
||||||
|
std::stringstream ss(line);
|
||||||
|
|
||||||
|
double k, Tkc, Tkb, Tktot, Tkvc, Tkvb;
|
||||||
|
ss >> k;
|
||||||
|
ss >> Tktot;
|
||||||
|
ss >> Tkc;
|
||||||
|
ss >> Tkb;
|
||||||
|
ss >> Tkvc;
|
||||||
|
ss >> Tkvb;
|
||||||
|
|
||||||
|
|
||||||
|
//std::cerr << k << " " << Tktot << " " << Tkc << std::endl;
|
||||||
|
|
||||||
|
m_tab_k.push_back( log10(k) );
|
||||||
|
|
||||||
|
m_tab_Tk_tot.push_back( log10(Tktot) );
|
||||||
|
m_tab_Tk_baryon.push_back( log10(Tkb) );
|
||||||
|
m_tab_Tk_cdm.push_back( log10(Tkc) );
|
||||||
|
m_tab_Tvk_cdm.push_back( log10(Tkvc) );
|
||||||
|
m_tab_Tvk_baryon.push_back( log10(Tkvb) );
|
||||||
|
|
||||||
|
}
|
||||||
|
|
||||||
|
ifs.close();
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
#ifdef WITH_MPI
|
||||||
|
}
|
||||||
|
|
||||||
|
unsigned n=m_tab_k.size();
|
||||||
|
MPI::COMM_WORLD.Bcast( &n, 1, MPI_UNSIGNED, 0 );
|
||||||
|
|
||||||
|
if( MPI::COMM_WORLD.Get_rank() > 0 ){
|
||||||
|
m_tab_k.assign(n,0);
|
||||||
|
m_tab_Tk_tot.assign(n,0);
|
||||||
|
m_tab_Tk_cdm.assign(n,0);
|
||||||
|
m_tab_Tk_baryon.assign(n,0);
|
||||||
|
m_tab_Tvk_cdm.assign(n,0);
|
||||||
|
m_tab_Tvk_baryon.assign(n,0);
|
||||||
|
}
|
||||||
|
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_k[0], n, MPI_DOUBLE, 0 );
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_Tk_tot[0], n, MPI_DOUBLE, 0 );
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_Tk_cdm[0], n, MPI_DOUBLE, 0 );
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_Tk_baryon[0], n, MPI_DOUBLE, 0 );
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_Tvk_cdm[0], n, MPI_DOUBLE, 0 );
|
||||||
|
MPI::COMM_WORLD.Bcast( &m_tab_Tvk_baryon[0], n, MPI_DOUBLE, 0 );
|
||||||
|
#endif
|
||||||
|
|
||||||
|
}
|
||||||
|
|
||||||
|
public:
|
||||||
|
transfer_MUSIC_plugin( config_file& cf )
|
||||||
|
: transfer_function_plugin( cf )
|
||||||
|
{
|
||||||
|
m_filename_Tk = pcf_->getValue<std::string>("cosmology","transfer_file");
|
||||||
|
|
||||||
|
read_table( );
|
||||||
|
|
||||||
|
acc_dtot = gsl_interp_accel_alloc();
|
||||||
|
acc_dcdm = gsl_interp_accel_alloc();
|
||||||
|
acc_dbaryon = gsl_interp_accel_alloc();
|
||||||
|
|
||||||
|
|
||||||
|
spline_dtot = gsl_spline_alloc( gsl_interp_akima, m_tab_k.size() );
|
||||||
|
spline_dcdm = gsl_spline_alloc( gsl_interp_akima, m_tab_k.size() );
|
||||||
|
spline_dbaryon = gsl_spline_alloc( gsl_interp_akima, m_tab_k.size() );
|
||||||
|
spline_vcdm = gsl_spline_alloc( gsl_interp_akima, m_tab_k.size() );
|
||||||
|
spline_vbaryon = gsl_spline_alloc( gsl_interp_akima, m_tab_k.size() );
|
||||||
|
|
||||||
|
gsl_spline_init (spline_dtot, &m_tab_k[0], &m_tab_Tk_tot[0], m_tab_k.size() );
|
||||||
|
gsl_spline_init (spline_dcdm, &m_tab_k[0], &m_tab_Tk_cdm[0], m_tab_k.size() );
|
||||||
|
gsl_spline_init (spline_dbaryon, &m_tab_k[0], &m_tab_Tk_baryon[0], m_tab_k.size() );
|
||||||
|
gsl_spline_init (spline_vcdm, &m_tab_k[0], &m_tab_Tvk_cdm[0], m_tab_k.size() );
|
||||||
|
gsl_spline_init (spline_vbaryon, &m_tab_k[0], &m_tab_Tvk_baryon[0], m_tab_k.size() );
|
||||||
|
|
||||||
|
tf_distinct_ = true;
|
||||||
|
tf_withvel_ = true;
|
||||||
|
}
|
||||||
|
|
||||||
|
~transfer_MUSIC_plugin()
|
||||||
|
{
|
||||||
|
gsl_spline_free (spline_dtot);
|
||||||
|
gsl_spline_free (spline_dcdm);
|
||||||
|
gsl_spline_free (spline_dbaryon);
|
||||||
|
gsl_spline_free (spline_vcdm);
|
||||||
|
gsl_spline_free (spline_vbaryon);
|
||||||
|
|
||||||
|
gsl_interp_accel_free (acc_dtot);
|
||||||
|
gsl_interp_accel_free (acc_dcdm);
|
||||||
|
gsl_interp_accel_free (acc_dbaryon);
|
||||||
|
gsl_interp_accel_free (acc_vcdm);
|
||||||
|
gsl_interp_accel_free (acc_vbaryon);
|
||||||
|
}
|
||||||
|
|
||||||
|
inline double extrap_left( double k, const tf_type& type )
|
||||||
|
{
|
||||||
|
if( k<1e-8 )
|
||||||
|
return 1.0;
|
||||||
|
|
||||||
|
double v1(1.0), v2(1.0);
|
||||||
|
switch( type )
|
||||||
|
{
|
||||||
|
case cdm:
|
||||||
|
v1 = m_tab_Tk_cdm[0];
|
||||||
|
v2 = m_tab_Tk_cdm[1];
|
||||||
|
break;
|
||||||
|
case baryon:
|
||||||
|
v1 = m_tab_Tk_baryon[0];
|
||||||
|
v2 = m_tab_Tk_baryon[1];
|
||||||
|
break;
|
||||||
|
case vcdm:
|
||||||
|
v1 = m_tab_Tvk_cdm[0];
|
||||||
|
v2 = m_tab_Tvk_cdm[1];
|
||||||
|
break;
|
||||||
|
case vbaryon:
|
||||||
|
v1 = m_tab_Tvk_baryon[0];
|
||||||
|
v2 = m_tab_Tvk_baryon[1];
|
||||||
|
break;
|
||||||
|
case total:
|
||||||
|
v1 = m_tab_Tk_tot[0];
|
||||||
|
v2 = m_tab_Tk_tot[1];
|
||||||
|
break;
|
||||||
|
|
||||||
|
default:
|
||||||
|
throw std::runtime_error("Invalid type requested in transfer function evaluation");
|
||||||
|
}
|
||||||
|
|
||||||
|
double lk = log10(k);
|
||||||
|
double dk = m_tab_k[1]-m_tab_k[0];
|
||||||
|
double delk = lk-m_tab_k[0];
|
||||||
|
|
||||||
|
//double xi = (v2-v1)/dk;
|
||||||
|
return pow(10.0,(v2-v1)/dk*(delk)+v1);
|
||||||
|
}
|
||||||
|
|
||||||
|
inline double extrap_right( double k, const tf_type& type )
|
||||||
|
{
|
||||||
|
double v1(1.0), v2(1.0);
|
||||||
|
|
||||||
|
int n=m_tab_k.size()-1, n1=n-1;
|
||||||
|
switch( type )
|
||||||
|
{
|
||||||
|
case cdm:
|
||||||
|
v1 = m_tab_Tk_cdm[n1];
|
||||||
|
v2 = m_tab_Tk_cdm[n];
|
||||||
|
break;
|
||||||
|
case baryon:
|
||||||
|
v1 = m_tab_Tk_baryon[n1];
|
||||||
|
v2 = m_tab_Tk_baryon[n];
|
||||||
|
break;
|
||||||
|
case vcdm:
|
||||||
|
v1 = m_tab_Tvk_cdm[n1];
|
||||||
|
v2 = m_tab_Tvk_cdm[n];
|
||||||
|
break;
|
||||||
|
case vbaryon:
|
||||||
|
v1 = m_tab_Tvk_baryon[n1];
|
||||||
|
v2 = m_tab_Tvk_baryon[n];
|
||||||
|
break;
|
||||||
|
case total:
|
||||||
|
v1 = m_tab_Tk_tot[n1];
|
||||||
|
v2 = m_tab_Tk_tot[n];
|
||||||
|
break;
|
||||||
|
|
||||||
|
default:
|
||||||
|
throw std::runtime_error("Invalid type requested in transfer function evaluation");
|
||||||
|
}
|
||||||
|
|
||||||
|
double lk = log10(k);
|
||||||
|
double dk = m_tab_k[n]-m_tab_k[n1];
|
||||||
|
double delk = lk-m_tab_k[n];
|
||||||
|
|
||||||
|
//double xi = (v2-v1)/dk;
|
||||||
|
return pow(10.0,(v2-v1)/dk*(delk)+v2);
|
||||||
|
}
|
||||||
|
|
||||||
|
inline double compute( double k, tf_type type ){
|
||||||
|
|
||||||
|
double lk = log10(k);
|
||||||
|
|
||||||
|
//if( lk<m_tab_k[1])
|
||||||
|
// return 1.0;
|
||||||
|
|
||||||
|
//if( lk>m_tab_k[m_tab_k.size()-2] );
|
||||||
|
// return m_tab_Tk_cdm[m_tab_k.size()-2]/k/k;
|
||||||
|
|
||||||
|
if( k<get_kmin() )
|
||||||
|
return extrap_left(k, type );
|
||||||
|
|
||||||
|
if( k>get_kmax() )
|
||||||
|
return extrap_right(k,type );
|
||||||
|
|
||||||
|
|
||||||
|
switch( type )
|
||||||
|
{
|
||||||
|
case cdm:
|
||||||
|
return pow(10.0, gsl_spline_eval (spline_dcdm, lk, acc_dcdm) );
|
||||||
|
case baryon:
|
||||||
|
return pow(10.0, gsl_spline_eval (spline_dbaryon, lk, acc_dbaryon) );
|
||||||
|
case vcdm:
|
||||||
|
return pow(10.0, gsl_spline_eval (spline_vcdm, lk, acc_vcdm) );
|
||||||
|
case vbaryon:
|
||||||
|
return pow(10.0, gsl_spline_eval (spline_vbaryon, lk, acc_vbaryon) );
|
||||||
|
case total:
|
||||||
|
return pow(10.0, gsl_spline_eval (spline_dtot, lk, acc_dtot) );
|
||||||
|
|
||||||
|
default:
|
||||||
|
throw std::runtime_error("Invalid type requested in transfer function evaluation");
|
||||||
|
}
|
||||||
|
|
||||||
|
return 1.0;
|
||||||
|
}
|
||||||
|
|
||||||
|
inline double get_kmin( void ){
|
||||||
|
return pow(10.0,m_tab_k[0]);
|
||||||
|
}
|
||||||
|
|
||||||
|
inline double get_kmax( void ){
|
||||||
|
return pow(10.0,m_tab_k[m_tab_k.size()-1]);
|
||||||
|
}
|
||||||
|
|
||||||
|
};
|
||||||
|
|
||||||
|
namespace{
|
||||||
|
transfer_function_plugin_creator_concrete< transfer_MUSIC_plugin > creator("music");
|
||||||
|
}
|
||||||
|
|
||||||
|
|
|
@ -46,7 +46,7 @@
|
||||||
#include "config_file.hh"
|
#include "config_file.hh"
|
||||||
|
|
||||||
enum tf_type{
|
enum tf_type{
|
||||||
total, cdm, baryon
|
total, cdm, baryon, vcdm, vbaryon
|
||||||
};
|
};
|
||||||
|
|
||||||
//! Abstract base class for transfer functions
|
//! Abstract base class for transfer functions
|
||||||
|
@ -58,7 +58,8 @@ class transfer_function_plugin{
|
||||||
public:
|
public:
|
||||||
Cosmology cosmo_; //!< cosmological parameter, read from config_file
|
Cosmology cosmo_; //!< cosmological parameter, read from config_file
|
||||||
config_file *pcf_; //!< pointer to config_file from which to read parameters
|
config_file *pcf_; //!< pointer to config_file from which to read parameters
|
||||||
bool tf_distinct_; //!< bool if transfer function is distinct for baryons and DM
|
bool tf_distinct_; //!< bool if density transfer function is distinct for baryons and DM
|
||||||
|
bool tf_withvel_; //!< bool if also have velocity transfer functions
|
||||||
|
|
||||||
public:
|
public:
|
||||||
|
|
||||||
|
@ -89,9 +90,13 @@ public:
|
||||||
//! return minimum wave number allowed
|
//! return minimum wave number allowed
|
||||||
virtual double get_kmin( void ) = 0;
|
virtual double get_kmin( void ) = 0;
|
||||||
|
|
||||||
//! return if transfer function is distinct for baryons and DM
|
//! return if density transfer function is distinct for baryons and DM
|
||||||
bool tf_is_distinct( void )
|
bool tf_is_distinct( void )
|
||||||
{ return tf_distinct_; }
|
{ return tf_distinct_; }
|
||||||
|
|
||||||
|
//! return if we also have velocity transfer functions
|
||||||
|
bool tf_has_velocities( void )
|
||||||
|
{ return tf_withvel_; }
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|
||||||
|
@ -278,7 +283,7 @@ protected:
|
||||||
in[i][1] = 0.0;
|
in[i][1] = 0.0;
|
||||||
|
|
||||||
sum_in += in[i][0];
|
sum_in += in[i][0];
|
||||||
ofsk << std::setw(16) << k <<std::setw(16) << in[i][0] << std::endl;
|
ofsk << std::setw(16) << k <<std::setw(16) << in[i][0] << std::setw(16) << T << std::endl;
|
||||||
#else
|
#else
|
||||||
|
|
||||||
#ifdef XI_SAMPLE
|
#ifdef XI_SAMPLE
|
||||||
|
@ -289,7 +294,7 @@ protected:
|
||||||
in[i].im = 0.0;
|
in[i].im = 0.0;
|
||||||
|
|
||||||
sum_in += in[i].re;
|
sum_in += in[i].re;
|
||||||
ofsk << std::setw(16) << k <<std::setw(16) << in[i].re << std::endl;
|
ofsk << std::setw(16) << k <<std::setw(16) << in[i].re << std::setw(16) << T << std::endl;
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
}
|
}
|
||||||
|
@ -509,43 +514,6 @@ public:
|
||||||
Tr0_ = res;
|
Tr0_ = res;
|
||||||
}
|
}
|
||||||
|
|
||||||
#if 0
|
|
||||||
{
|
|
||||||
double sum = 0.0;
|
|
||||||
unsigned long long count = 0;
|
|
||||||
|
|
||||||
int nf2 = nfull/2;
|
|
||||||
double dk = 2.0*M_PI/boxlength;
|
|
||||||
|
|
||||||
#pragma omp parallel for reduction(+:sum,count)
|
|
||||||
for( int i=0; i<=nf2; ++i )
|
|
||||||
for( int j=0; j<=nf2; ++j )
|
|
||||||
for( int k=0; k<=nf2; ++k )
|
|
||||||
{
|
|
||||||
int ii(i),jj(j),kk(k);
|
|
||||||
double rr = (double)ii*(double)ii+(double)jj*(double)jj+(double)kk*(double)kk;
|
|
||||||
double kp = sqrt(rr)*dk;
|
|
||||||
|
|
||||||
if( i==0||j==0||k==0||i==nf2||j==nf2||k==nf2 )
|
|
||||||
{
|
|
||||||
sum += pow(kp,0.5*nspec_)*ptf_->compute( kp, type_ );
|
|
||||||
count++;
|
|
||||||
}
|
|
||||||
else
|
|
||||||
{
|
|
||||||
sum += 2.0*pow(kp,0.5*nspec_)*ptf_->compute( kp, type_ );
|
|
||||||
count+=2;
|
|
||||||
}
|
|
||||||
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
Tr0_ = pow(2.0*M_PI/boxlength*nfull,3.0)*dplus*sqrt(pnorm)*sum/(double)count;
|
|
||||||
|
|
||||||
}
|
|
||||||
#endif
|
|
||||||
|
|
||||||
|
|
||||||
/*****************************************************************/
|
/*****************************************************************/
|
||||||
//... store as table for spline interpolation
|
//... store as table for spline interpolation
|
||||||
|
|
||||||
|
@ -653,7 +621,7 @@ public:
|
||||||
|
|
||||||
gsl_integration_workspace_free (wx);
|
gsl_integration_workspace_free (wx);
|
||||||
|
|
||||||
return resx;
|
return resx*cos(0.5*M_PI*kx/kmax);
|
||||||
}
|
}
|
||||||
|
|
||||||
static double call_y( double ky, void *arg )
|
static double call_y( double ky, void *arg )
|
||||||
|
@ -674,24 +642,23 @@ public:
|
||||||
|
|
||||||
gsl_integration_workspace_free (wy);
|
gsl_integration_workspace_free (wy);
|
||||||
|
|
||||||
return resy;
|
return resy*cos(0.5*M_PI*ky/kmax);
|
||||||
}
|
}
|
||||||
|
|
||||||
static double call_z( double kz, void *arg )
|
static double call_z( double kz, void *arg )
|
||||||
{
|
{
|
||||||
double *a = (double*)arg;
|
double *a = (double*)arg;
|
||||||
//double kmin = a[3], kmax = a[4];
|
//double kmin = a[3], kmax = a[4];
|
||||||
double kx = a[0], ky = a[1];
|
double kx = a[0], ky = a[1], kmax = a[4];
|
||||||
|
|
||||||
double kk = sqrt(kx*kx+ky*ky+kz*kz);
|
double kk = sqrt(kx*kx+ky*ky+kz*kz);
|
||||||
double T = ptf_->compute( kk, type_ );
|
double T = ptf_->compute( kk, type_ );
|
||||||
|
|
||||||
|
|
||||||
return pow(kk,0.5*nspec_)*T;
|
return pow(kk,0.5*nspec_)*T*cos(0.5*M_PI*kz/kmax);
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
~TransferFunction_real()
|
~TransferFunction_real()
|
||||||
{ }
|
{ }
|
||||||
|
|
||||||
|
|
Loading…
Reference in a new issue