CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
ModelFunction_TwoPointCorrelation_projected.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 
47 std::vector<double> cbl::modelling::twopt::wp_from_xi_approx (FunctionVectorVectorPtrVectorRef func, const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
48 {
49  // structure contaning the required input data
50  shared_ptr<STR_data_HOD> pp = static_pointer_cast<STR_data_HOD>(inputs);
51 
52  vector<double> wp(rp.size(), 0);
53 
54 #pragma omp parallel num_threads(omp_get_max_threads())
55  {
56 
57 #pragma omp for schedule(static, 2)
58  for (size_t i=0; i<wp.size(); i++) {
59 
60  auto integrand = [&] (double rad) { return func({rad}, inputs, parameter)[0]/sqrt(rad*rad-rp[i]*rp[i]); };
61 
62  double r_out = sqrt(pow(rp[i], 2)+pow(pp->r_max_int, 2));
63  wp[i] = 2.*wrapper::gsl::GSL_integrate_qag(integrand, rp[i], r_out);
64 
65  }
66 
67  }
68 
69  return wp;
70 }
71 
72 
73 // ============================================================================================
74 
75 
76 std::vector<double> cbl::modelling::twopt::wp_1halo_approx (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
77 {
78  return wp_from_xi_approx(modelling::twopt::xi_1halo, rp, inputs, parameter);
79 }
80 
81 
82 // ============================================================================================
83 
84 
85 std::vector<double> cbl::modelling::twopt::wp_2halo_approx (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
86 {
87  return wp_from_xi_approx(modelling::twopt::xi_2halo, rp, inputs, parameter);
88 }
89 
90 
91 // ============================================================================================
92 
93 
94 std::vector<double> cbl::modelling::twopt::wp_HOD_approx (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
95 {
96  return wp_from_xi_approx(modelling::twopt::xi_HOD, rp, inputs, parameter);
97 }
98 
99 
100 // ============================================================================================
101 
102 
103 std::vector<double> cbl::modelling::twopt::wp_from_xi (FunctionDoubleDoubleDoublePtrVectorRef func, const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
104 {
105  // structure contaning the required input data
106  shared_ptr<STR_data_HOD> pp = static_pointer_cast<STR_data_HOD>(inputs);
107 
108  vector<double> wp(rp.size(), 0);
109 
110 #pragma omp parallel num_threads(omp_get_max_threads())
111  {
112 
113 #pragma omp for schedule(static, 2)
114  for (size_t i=0; i<wp.size(); i++) {
115 
116  auto integrand = [&] (double pi) { return func(rp[i], pi, inputs, parameter); };
117 
118  wp[i] = 2.*wrapper::gsl::GSL_integrate_qag(integrand, 0., pp->pi_max);
119 
120  }
121 
122  }
123 
124  return wp;
125 }
126 
127 
128 // ============================================================================================
129 
130 
131 std::vector<double> cbl::modelling::twopt::wp_1halo (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
132 {
133  return wp_from_xi(modelling::twopt::xi_1halo_zspace, rp, inputs, parameter);
134 }
135 
136 
137 // ============================================================================================
138 
139 
140 std::vector<double> cbl::modelling::twopt::wp_2halo (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
141 {
142  return wp_from_xi(modelling::twopt::xi_2halo_zspace, rp, inputs, parameter);
143 }
144 
145 
146 // ============================================================================================
147 
148 
149 std::vector<double> cbl::modelling::twopt::wp_HOD (const std::vector<double> rp, const std::shared_ptr<void> inputs, std::vector<double> &parameter)
150 {
151  return wp_from_xi(modelling::twopt::xi_HOD_zspace, rp, inputs, parameter);
152 }
Functions to model the projected two-point correlation function.
static const double pi
: the ratio of a circle's circumference to its diameter
Definition: Constants.h:199
std::vector< double > wp_from_xi_approx(FunctionVectorVectorPtrVectorRef func, const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
function used to compute the projected two-point correlation function
double xi_HOD_zspace(const double rp, const double pi, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
HOD model of the redshift-space monopole of the two-point correlation function.
std::vector< double > wp_HOD_approx(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
HOD model of the projected two-point correlation function.
std::vector< double > xi_HOD(const std::vector< double > rad, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
HOD model of the monopole of the two-point correlation function.
std::vector< double > wp_2halo(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 2-halo term of the projected two-point correlation function
std::vector< double > wp_1halo(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 1-halo term of the projected two-point correlation function
double xi_2halo_zspace(const double rp, const double pi, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 2-halo term of the redshift-space monopole of the two-point correlation function
std::vector< double > wp_1halo_approx(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 1-halo term of the projected two-point correlation function
double xi_1halo_zspace(const double rp, const double pi, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 1-halo term of the redshift-space monopole of the two-point correlation function
std::vector< double > wp_from_xi(FunctionDoubleDoubleDoublePtrVectorRef func, const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
function used to compute the projected two-point correlation function
std::vector< double > xi_2halo(const std::vector< double > rad, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 2-halo term of the monopole of the two-point correlation function
std::vector< double > wp_HOD(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
HOD model of the projected two-point correlation function.
std::vector< double > wp_2halo_approx(const std::vector< double > rp, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 2-halo term of the projected two-point correlation function
std::vector< double > xi_1halo(const std::vector< double > rad, const std::shared_ptr< void > inputs, std::vector< double > &parameter)
model for the 1-halo term of the monopole of the two-point correlation function
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::function< std::vector< double >std::vector< double >, std::shared_ptr< void >, std::vector< double > &)> FunctionVectorVectorPtrVectorRef
typedef of a function returning a vector with a vector, a pointer and a vector reference in input
Definition: Kernel.h:705
std::function< double(double, double, std::shared_ptr< void >, std::vector< double > &)> FunctionDoubleDoubleDoublePtrVectorRef
typedef of a function returning a double with two double, a pointer and a vector reference in input
Definition: Kernel.h:699
double wp(const double rp, const std::vector< double > rr, const std::vector< double > xi, const double r_max=100.)
the projected two-point correlation function
Definition: FuncXi.cpp:192