CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
Modelling_TwoPointCorrelation2D_cartesian.cpp
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 
38 
39 using namespace std;
40 
41 using namespace cbl;
42 
43 
44 // ============================================================================================
45 
46 
48 {
49  cout << endl; coutCBL << "Setting up the fiducial dark matter two-point correlation function model" << endl;
50 
51  const vector<double> rad = logarithmic_bin_vector(m_data_model->step, max(m_data_model->r_min, 1.e-4), min(m_data_model->r_max, 100.));
52  vector<double> xiDM(m_data_model->step);
53 
54  if (m_data_model->sigmaNL==0) {
55  for (size_t i=0; i<(size_t)m_data_model->step; i++)
56  xiDM[i] = m_data_model->cosmology->xi_matter(rad[i], m_data_model->method_Pk, m_data_model->NL, m_data_model->redshift, true, m_data_model->output_root, m_data_model->norm, m_data_model->k_min, m_data_model->k_max, m_data_model->aa, m_data_model->GSL, m_data_model->prec, m_data_model->file_par);
57  if (!m_data_model->store_output)
58  m_data_model->cosmology->remove_output_Pk_tables(m_data_model->method_Pk, m_data_model->NL, m_data_model->redshift, m_data_model->output_root);
59  }
60 
61  else {
62  const vector<double> kk = logarithmic_bin_vector(m_data_model->step, max(m_data_model->k_min, 1.e-4), min(m_data_model->k_max, 500.));
63  vector<double> Pk = m_data_model->cosmology->Pk_matter_DeWiggled("CAMB", "EisensteinHu", kk, m_data_model->redshift, m_data_model->sigmaNL, 4, 10, 0.05, m_data_model->store_output, m_data_model->output_root, m_data_model->norm, m_data_model->prec);
64  xiDM = Xi0(rad, kk, Pk);
65  }
66 
67  m_data_model->func_xi = make_shared<glob::FuncGrid>(glob::FuncGrid(rad, xiDM, "Spline"));
68 
69  vector<double> xiDM_(m_data_model->step), xiDM__(m_data_model->step);
70 
71  auto integrand_ = [&] (const double rr) {
72  return m_data_model->func_xi->operator()(rr)*rr*rr;
73  };
74 
75  auto integrand__ = [&] (const double rr) {
76  return m_data_model->func_xi->operator()(rr)*rr*rr*rr*rr;
77  };
78 
79  for (size_t i=0; i<(size_t)m_data_model->step; i++) {
80  xiDM_[i] = 3.*wrapper::gsl::GSL_integrate_qag(integrand_, 0., rad[i])*pow(rad[i], -3);
81  xiDM__[i] = 5.*wrapper::gsl::GSL_integrate_qag(integrand__, 0., rad[i])*pow(rad[i], -5);
82  }
83 
84  m_data_model->func_xi_ = make_shared<glob::FuncGrid>(glob::FuncGrid(rad, xiDM_, "Spline"));
85  m_data_model->func_xi__ = make_shared<glob::FuncGrid>(glob::FuncGrid(rad, xiDM__, "Spline"));
86 }
87 
88 
89 // ============================================================================================
90 
91 
93 {
94  // compute the fiducial dark matter two-point correlation function and associated quantities
95  set_fiducial_xiDM();
96 
97  // set the model parameters
98  const int nparameters = 5;
99 
100  vector<statistics::ParameterType> parameterType(nparameters, statistics::ParameterType::_Base_);
101 
102  vector<string> parameterName(nparameters);
103 
104  parameterName[0] = "f*sigma8";
105  parameterName[1] = "b*sigma8";
106  parameterName[2] = "sigmav";
107  parameterName[3] = "alpha_perp";
108  parameterName[4] = "alpha_par";
109 
110  vector<statistics::PriorDistribution> priors = {fsigma8_prior, bsigma8_prior, sigmav_prior, alpha_perp_prior, alpha_par_prior};
111 
112  // set the priors
113  m_set_prior(priors);
114 
115  // construct the model
116  m_model = make_shared<statistics::Model2D>(statistics::Model2D(&xi2D_dispersion, nparameters, parameterType, parameterName, m_data_model));
117 
118 }
#define coutCBL
CBL print message.
Definition: Kernel.h:734
The class Modelling_TwoPointCorrelation2D_cartesian.
The class FuncGrid.
Definition: FuncGrid.h:55
void set_fiducial_xiDM()
set the fiducial model for dark matter two point correlation function
void set_model_dispersion(const statistics::PriorDistribution fsigma8_prior={}, const statistics::PriorDistribution bsigma8_prior={}, const statistics::PriorDistribution sigmav_prior={}, const statistics::PriorDistribution alpha_perp_prior={cbl::glob::DistributionType::_Constant_, 1.}, const statistics::PriorDistribution alpha_par_prior={cbl::glob::DistributionType::_Constant_, 1.})
set the dispersion model to fit the 2D two-point correlation function, in Cartesian coordinates
The class Model2D.
Definition: Model2D.h:60
The class PriorDistribution.
std::vector< std::vector< double > > xi2D_dispersion(const std::vector< double > rp, const std::vector< double > pi, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 2D two-point correlation function, in Cartesian coordinates
double GSL_integrate_qag(gsl_function Func, const double a, const double b, const double rel_err=1.e-3, const double abs_err=0, const int limit_size=1000, const int rule=6)
integral, computed using the GSL qag method
Definition: GSLwrapper.cpp:180
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38
std::vector< T > logarithmic_bin_vector(const size_t nn, const T min, const T max)
fill a std::vector with logarithmically spaced values
Definition: Kernel.h:1621
std::vector< double > Xi0(const std::vector< double > r, const std::vector< double > kk, const std::vector< double > Pk0, const double k_cut=0.7, const double cut_pow=2, const int IntegrationMethod=1)
function to obtain the two point correlation funciton monopole