CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
ModelFunction_NumberCounts.h
Go to the documentation of this file.
1 /********************************************************************
2  * Copyright (C) 2016 by Federico Marulli and Alfonso Veropalumbo *
3  * federico.marulli3@unibo.it *
4  * *
5  * This program is free software; you can redistribute it and/or *
6  * modify it under the terms of the GNU General Public License as *
7  * published by the Free Software Foundation; either version 2 of *
8  * the License, or (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public *
16  * License along with this program; if not, write to the Free *
17  * Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  ********************************************************************/
20 
35 #ifndef __MODFUNCNC__
36 #define __MODFUNCNC__
37 
38 #include "Cosmology.h"
39 #include "NumberCounts.h"
40 
41 
42 // ============================================================================
43 
44 
45 namespace cbl {
46 
47  namespace modelling {
48 
49  namespace numbercounts {
50 
59 
61  bool isSnapshot;
62 
64  std::vector<double> edges_x;
65 
67  std::vector<double> edges_y;
68 
70  std::shared_ptr<cosmology::Cosmology> cosmology;
71 
73  std::vector<cosmology::CosmologicalParameter> Cpar;
74 
76  std::function<double(const double, const double, const std::shared_ptr<void>)> fz;
77 
79  double (*response_fact)(const double, const double, const double, const double, const std::string, const double, const std::string, std::shared_ptr<void>);
80 
82  std::function<double(const double, const double)> z_error;
83 
85  std::function<double(const double, const double)> proxy_error;
86 
88  double redshift;
89 
91  std::string method_Pk;
92 
94  bool NL;
95 
97  double k_min;
98 
100  double k_max;
101 
103  int step;
104 
106  std::vector<double> kk;
107 
110 
112  std::string output_root;
113 
115  int norm;
116 
118  std::string file_par;
119 
121  double prec;
122 
124  double Delta;
125 
128 
130  std::string model_MF;
131 
133  std::string model_bias;
134 
136  double z_min;
137 
139  double z_max;
140 
142  int z_step;
143 
145  std::vector<double> z_vector;
146 
148  double Mass_min;
149 
151  double Mass_max;
152 
155 
157  std::vector<double> Mass_vector;
158 
160  double area_rad;
161 
163  double Volume;
164 
166  bool use_SF;
167 
170 
172  std::shared_ptr<glob::FuncGrid2D> interp_SelectionFunction;
173 
175  std::vector<double> SF_weights;
176 
178  double z_pivot;
179 
181  double proxy_pivot;
182 
184  double mass_pivot;
185 
187  double log_base;
188 
192  STR_NC_data_model () = default;
193  };
194 
195 
196 
205 
207  std::vector<double> radii;
208 
210  std::shared_ptr<cosmology::Cosmology> cosmology;
211 
213  std::vector<cosmology::CosmologicalParameter> Cpar;
214 
216  double redshift;
217 
219  std::string model_SF;
220 
222  double b_eff;
223 
225  double b_slope;
226 
228  double b_offset;
229 
231  double deltav_NL;
232 
234  double delta_c;
235 
237  std::string method_Pk;
238 
241  double k_Pk_ratio;
242 
247 
250  std::string output_root;
251 
253  std::string interpType;
254 
257  double k_max;
258 
260  std::string input_file;
261 
265 
269  STR_NCSF_data_model () = default;
270  };
271 
272 
282  double Filter_sigmaR (const double kk, const double radius);
283 
294  double Filter_dsigmaR (const double kk, const double radius);
295 
313  void sigmaM_dlnsigmaM (double &sigmaM, double &dlnsigmaM, const double mass, const cbl::glob::FuncGrid interp_Pk, const double kmax, const double rho);
314 
336  void sigmaM_dlnsigmaM (std::vector<double> &sigmaM, std::vector<double> &dlnsigmaM, const std::vector<double> mass, const std::vector<double> kk, const std::vector<double> Pk, const std::string interpType, const double kmax, const double rho);
337 
357  std::vector<cbl::glob::FuncGrid> sigmaM_dlnsigmaM (const std::vector<double> mass, cosmology::Cosmology cosmology, const std::vector<double> kk, const std::vector<double> Pk, const std::string interpType, const double kmax);
358 
397  double mass_function (const double mass, cosmology::Cosmology cosmology, const double redshift, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const cbl::glob::FuncGrid interp_Pk, const double kmax);
398 
441  std::vector<double> mass_function (const std::vector<double> mass, cosmology::Cosmology cosmology, const double redshift, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const std::vector<double> kk, const std::vector<double> Pk, const std::string interpType, const double kmax);
442 
487  std::vector<std::vector<double>> mass_function (const std::vector<double> redshift, const std::vector<double> mass, cosmology::Cosmology cosmology, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const std::vector<double> kk, const std::vector<double> Pk, const std::string interpType, const double kmax);
488 
555  std::vector<double> size_function (cosmology::Cosmology cosmology, const std::vector<double> radii, const double redshift, const std::string model, const double b_eff, double slope=0.854, double offset=0.420, const double deltav_NL=-0.795, const double del_c=1.69, const std::string method_Pk="Eisensteinhu", const bool store_output=true, const std::string output_root="test", const std::string interpType="Linear", const double k_max=100., const std::string input_file=par::defaultString, const bool is_parameter_file=true);
556 
605  double number_counts (const double redshift_min, const double redshift_max, const double Mass_min, const double Mass_max, cosmology::Cosmology cosmology, const double Area, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const glob::FuncGrid interp_sigmaM, const glob::FuncGrid interp_DlnsigmaM);
606 
709  double counts_proxy (const double alpha, const double beta, const double gamma, const double scatter0, const double scatterM, const double scatterM_exp, const double scatterz, const double scatterz_exp, const double z_bias, const double proxy_bias, const double z_err, const double proxy_err, const double Plambda_a, const double Plambda_b, const double Plambda_c, std::function<double(const double, const double, const std::shared_ptr<void>)> fz, std::function<double(const double, const double)> z_error, std::function<double(const double, const double)> proxy_error, double (*response_fact)(const double, const double, const double, const double, const std::string, const double, const std::string, std::shared_ptr<void>), const double redshift_min, const double redshift_max, const double proxy_min, const double proxy_max, cbl::cosmology::Cosmology cosmology, const double Area, const std::string model_MF, const std::string model_bias, const bool store_output, const double Delta, const bool isDelta_critical, const cbl::glob::FuncGrid interp_sigmaM, const cbl::glob::FuncGrid interp_DlnsigmaM, const cbl::glob::FuncGrid interp_DN, const double proxy_pivot, const double z_pivot, const double mass_pivot, const double log_base, const double weight);
710 
711  }
712  }
713 }
714 
715 #endif
The class Cosmology.
The class NumberCounts.
The class Cosmology.
Definition: Cosmology.h:277
The class FuncGrid.
Definition: FuncGrid.h:55
static const std::string defaultString
default std::string value
Definition: Constants.h:336
static const double alpha
: the fine-structure constant
Definition: Constants.h:233
void sigmaM_dlnsigmaM(double &sigmaM, double &dlnsigmaM, const double mass, const cbl::glob::FuncGrid interp_Pk, const double kmax, const double rho)
compute
double counts_proxy(const double alpha, const double beta, const double gamma, const double scatter0, const double scatterM, const double scatterM_exp, const double scatterz, const double scatterz_exp, const double z_bias, const double proxy_bias, const double z_err, const double proxy_err, const double Plambda_a, const double Plambda_b, const double Plambda_c, std::function< double(const double, const double, const std::shared_ptr< void >)> fz, std::function< double(const double, const double)> z_error, std::function< double(const double, const double)> proxy_error, double(*response_fact)(const double, const double, const double, const double, const std::string, const double, const std::string, std::shared_ptr< void >), const double redshift_min, const double redshift_max, const double proxy_min, const double proxy_max, cbl::cosmology::Cosmology cosmology, const double Area, const std::string model_MF, const std::string model_bias, const bool store_output, const double Delta, const bool isDelta_critical, const cbl::glob::FuncGrid interp_sigmaM, const cbl::glob::FuncGrid interp_DlnsigmaM, const cbl::glob::FuncGrid interp_DN, const double proxy_pivot, const double z_pivot, const double mass_pivot, const double log_base, const double weight)
compute the number counts as function of mass proxy and redshift
double mass_function(const double mass, cosmology::Cosmology cosmology, const double redshift, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const cbl::glob::FuncGrid interp_Pk, const double kmax)
compute the mass function
std::vector< double > size_function(cosmology::Cosmology cosmology, const std::vector< double > radii, const double redshift, const std::string model, const double b_eff, double slope=0.854, double offset=0.420, const double deltav_NL=-0.795, const double del_c=1.69, const std::string method_Pk="Eisensteinhu", const bool store_output=true, const std::string output_root="test", const std::string interpType="Linear", const double k_max=100., const std::string input_file=par::defaultString, const bool is_parameter_file=true)
the void size function
double number_counts(const double redshift_min, const double redshift_max, const double Mass_min, const double Mass_max, cosmology::Cosmology cosmology, const double Area, const std::string model_MF, const bool store_output, const double Delta, const bool isDelta_critical, const glob::FuncGrid interp_sigmaM, const glob::FuncGrid interp_DlnsigmaM)
compute the number counts as function of mass and redshift
double Filter_sigmaR(const double kk, const double radius)
the filter to compute
double Filter_dsigmaR(const double kk, const double radius)
the filter to compute
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38
std::vector< cosmology::CosmologicalParameter > Cpar
cosmological parameters
std::shared_ptr< cosmology::Cosmology > cosmology
fiducial cosmology
STR_NCSF_data_model()=default
default constructor
double delta_c
critical value of the linear density field
std::string input_file
either the parameter file or the power spectrum file
std::string interpType
interpType method to interpolate the power spectrum
double b_offset
second coefficent to convert the effective bias
std::string method_Pk
method to compute the dark matter power spectrum
double b_slope
first coefficent to convert the effective bias
double z_pivot
redshift pivot in the mass-observable scaling relation
std::shared_ptr< cosmology::Cosmology > cosmology
fiducial cosmology
double proxy_pivot
mass proxy pivot in the mass-observable scaling relation
std::vector< double > kk
vector of wave vector modules
std::vector< cosmology::CosmologicalParameter > Cpar
cosmological parameters
double k_max
maximum wave vector module up to which the power spectrum is computed
std::vector< double > SF_weights
number counts weights derived from the selection function
std::string method_Pk
method to compute the dark matter power spectrum
double(* response_fact)(const double, const double, const double, const double, const std::string, const double, const std::string, std::shared_ptr< void >)
pointer to the multiplicative term in the response function
std::string output_root
output root of the parameter file used to compute the dark matter power spectrum
bool use_SF
false → don't use the selection function; true → use the selection function
bool store_output
true the output files created by the Boltzmann solver are stored; false the output files are remove...
STR_NC_data_model()=default
default constructor
std::function< double(const double, const double, const std::shared_ptr< void >)> fz
the redshift evolution function in the scaling relation
std::function< double(const double, const double)> proxy_error
function returning the absolute mass proxy error
bool isSnapshot
false → data not from a simulation snapshot; true → data from a simulation snapshot
std::function< double(const double, const double)> z_error
function returning the absolute redshift error
std::vector< double > z_vector
vector of redshifts
double log_base
the base of the logarithm used in the mass-mass proxy scaling relation
double mass_pivot
mass pivot in the mass-observable scaling relation
int norm
0 → don't normalize the power spectrum; 1 → normalize the power spectrum
std::string model_bias
author(s) who proposed the bias function
int step
number of steps used to compute the binned dark matter correlation function
std::vector< double > edges_x
the x variable bin edges
bool NL
false → linear power spectrum; true → non-linear power spectrum
double k_min
minimum wave vector module up to which the power spectrum is computed
int Mass_step
number of mass steps used to compute the binned mass function
std::vector< double > edges_y
the y variable bin edges
bool is_sigma8_free
true → sigma8 is a free parameter; false → sigma8 can be considered a derived parameter
std::string model_MF
author(s) who proposed the mass function
std::shared_ptr< glob::FuncGrid2D > interp_SelectionFunction
function to interpolate the selection function in mass and redshift
int z_step
number of redshift steps used to compute the binned mass function