CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
TwoPointCorrelation2D_cartesian.cpp
Go to the documentation of this file.
1 /********************************************************************
2  * Copyright (C) 2015 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 
39 
40 using namespace std;
41 
42 using namespace cbl;
43 using namespace catalogue;
44 using namespace chainmesh;
45 using namespace data;
46 using namespace pairs;
47 using namespace measure;
48 using namespace twopt;
49 
50 
51 // ============================================================================================
52 
53 
54 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::set_parameters (const BinType binType_rp, const double rpMin, const double rpMax, const int nbins_rp, const double shift_rp, const BinType binType_pi, const double piMin, const double piMax, const int nbins_pi, const double shift_pi, const CoordinateUnits angularUnits, std::function<double(double)> angularWeight, const bool compute_extra_info)
55 {
56  if (!compute_extra_info) {
57  if (binType_rp==BinType::_logarithmic_) {
58  if (binType_pi==BinType::_logarithmic_) {
59  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
60  }
61  else {
62  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
63  }
64  }
65  else {
66  if (binType_pi==BinType::_logarithmic_) {
67  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
68  }
69  else {
70  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
71  }
72  }
73  }
74  else {
75  if (binType_rp==BinType::_logarithmic_) {
76  if (binType_pi==BinType::_logarithmic_) {
77  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_extra_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
78  }
79  else {
80  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_extra_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
81  }
82  }
83  else {
84  if (binType_pi==BinType::_logarithmic_) {
85  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_extra_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
86  }
87  else {
88  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_extra_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits, angularWeight));
89  }
90  }
91  }
92 
93  if (binType_rp==BinType::_logarithmic_) {
94  if (binType_pi==BinType::_logarithmic_) {
95  m_rr = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
96  m_dr = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
97  }
98  else {
99  m_rr = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
100  m_dr = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
101  }
102  }
103  else {
104  if (binType_pi==BinType::_logarithmic_) {
105  m_rr = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
106  m_dr = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
107  }
108  else {
109  m_rr = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
110  m_dr = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, nbins_rp, shift_rp, piMin, piMax, nbins_pi, shift_pi, angularUnits));
111  }
112  }
113 }
114 
115 
116 // ============================================================================================
117 
118 
119 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::set_parameters (const BinType binType_rp, const double rpMin, const double rpMax, const double binSize_rp, const double shift_rp, const BinType binType_pi, const double piMin, const double piMax, const double binSize_pi, const double shift_pi, const CoordinateUnits angularUnits, std::function<double(double)> angularWeight, const bool compute_extra_info)
120 {
121  if (!compute_extra_info) {
122  if (binType_rp==BinType::_logarithmic_) {
123  if (binType_pi==BinType::_logarithmic_) {
124  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
125  }
126  else {
127  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
128  }
129  }
130  else {
131  if (binType_pi==BinType::_logarithmic_) {
132  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
133  }
134  else {
135  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
136  }
137  }
138  }
139  else {
140  if (binType_rp==BinType::_logarithmic_) {
141  if (binType_pi==BinType::_logarithmic_) {
142  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_extra_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
143  }
144  else {
145  m_dd = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_extra_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
146  }
147  }
148  else {
149  if (binType_pi==BinType::_logarithmic_) {
150  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_extra_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
151  }
152  else {
153  m_dd = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_extra_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits, angularWeight));
154  }
155  }
156  }
157 
158 
159  if (binType_rp==BinType::_logarithmic_) {
160  if (binType_pi==BinType::_logarithmic_) {
161  m_rr = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
162  m_dr = move(Pair::Create(PairType::_comovingCartesian_loglog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
163  }
164  else {
165  m_rr = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
166  m_dr = move(Pair::Create(PairType::_comovingCartesian_loglin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
167  }
168  }
169  else {
170  if (binType_pi==BinType::_logarithmic_) {
171  m_rr = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
172  m_dr = move(Pair::Create(PairType::_comovingCartesian_linlog_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
173  }
174  else {
175  m_rr = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
176  m_dr = move(Pair::Create(PairType::_comovingCartesian_linlin_, PairInfo::_standard_, rpMin, rpMax, binSize_rp, shift_rp, piMin, piMax, binSize_pi, shift_pi, angularUnits));
177  }
178  }
179 }
180 
181 
182 // ============================================================================================
183 
184 
185 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::read (const std::string dir, const std::string file)
186 {
187  m_dataset->read(dir+file);
188 }
189 
190 
191 // ============================================================================================
192 
193 
194 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::write (const std::string dir, const std::string file, const bool full, const int rank) const
195 {
196  vector<double> xx = m_dataset->xx(), yy = m_dataset->yy();
197 
198  checkDim(xx, m_dd->nbins_D1(), "rp"); checkDim(yy, m_dd->nbins_D2(), "pi");
199 
200  string header = "[1] perpendicular separation at the bin centre # [2] parallel separation at the bin centre # [3] 2D two-point correlation function # [4] error";
201  if (m_compute_extra_info) header += " # [5] mean perpendicular separation # [6] standard deviation of the distribution of perpendicular separations # [7] mean parallel separation # [8] standard deviation of the distribution of parallel separations # [9] mean redshift # [10] standard deviation of the redshift distribution";
202 
203  m_dataset->write(dir, file, header, full, 5, rank);
204 }
205 
206 
207 // ============================================================================================
208 
209 
210 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::measure (const ErrorType errorType, const std::string dir_output_pairs, const std::vector<std::string> dir_input_pairs, const std::string dir_output_resample, const int nMocks, const bool count_dd, const bool count_rr, const bool count_dr, const bool tcount, const Estimator estimator, const double fact, const int seed)
211 {
212  switch (errorType) {
213  case (ErrorType::_Poisson_) :
214  measurePoisson(dir_output_pairs, dir_input_pairs, count_dd, count_rr, count_dr, tcount, estimator, fact);
215  break;
216  case (ErrorType::_Jackknife_) :
217  measureJackknife(dir_output_pairs, dir_input_pairs, dir_output_resample, count_dd, count_rr, count_dr, tcount, estimator, fact);
218  break;
219  case (ErrorType::_Bootstrap_) :
220  measureBootstrap(nMocks, dir_output_pairs, dir_input_pairs, dir_output_resample, count_dd, count_rr, count_dr, tcount, estimator, fact, seed);
221  break;
222  default:
223  ErrorCBL("unknown type of error!", "measure", "TwoPointCorrelation2D_cartesian.cpp");
224  }
225 }
226 
227 
228 // ============================================================================================
229 
230 
231 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::measurePoisson (const std::string dir_output_pairs, const std::vector<std::string> dir_input_pairs, const bool count_dd, const bool count_rr, const bool count_dr, const bool tcount, const Estimator estimator, const double fact)
232 {
233  // ----------- count the data-data, random-random and data-random pairs, or read them from file -----------
234 
235  count_allPairs(m_twoPType, dir_output_pairs, dir_input_pairs, count_dd, count_rr, count_dr, tcount, estimator, fact);
236 
237 
238  // ----------- compute the monopole of the two-point correlation function -----------
239 
240  if (estimator==Estimator::_natural_)
241  m_dataset = correlation_NaturalEstimator(m_dd, m_rr);
242  else if (estimator==Estimator::_LandySzalay_) {
243  m_dataset = correlation_LandySzalayEstimator(m_dd, m_rr, m_dr);
244  }
245  else
246  ErrorCBL("the chosen estimator is not implemented!", "measurePoisson", "TwoPointCorrelation2D_cartesian.cpp");
247 
248 }
249 
250 
251 // ============================================================================================
252 
253 
254 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::measureJackknife (const std::string dir_output_pairs, const std::vector<std::string> dir_input_pairs, const std::string dir_output_resample, const bool count_dd, const bool count_rr, const bool count_dr, const bool tcount, const Estimator estimator, const double fact)
255 {
256  if (dir_output_resample!=par::defaultString && dir_output_resample!="") {
257  string mkdir = "mkdir -p "+dir_output_resample;
258  if (system(mkdir.c_str())) {}
259  }
260 
261  vector<long> region_list = m_data->region_list();
262  size_t nRegions = region_list.size();
263 
264  vector< vector<double > > xi_SubSamples, covariance;
265 
266  vector<shared_ptr<Pair> > dd_regions, rr_regions, dr_regions;
267  count_allPairs_region(dd_regions, rr_regions, dr_regions, m_twoPType, dir_output_pairs, dir_input_pairs, count_dd, count_rr, count_dr, tcount, estimator, fact);
268 
269  vector<shared_ptr<Data> > data_SS = (estimator==Estimator::_natural_) ? XiJackknife(dd_regions, rr_regions) : XiJackknife(dd_regions, rr_regions, dr_regions);
270 
271  for (size_t i=0; i<nRegions; i++) {
272 
273  if (dir_output_resample!=par::defaultString && dir_output_resample!="") {
274  string file = "xi_Jackknife_"+conv(i, par::fINT)+".dat";
275  string header = "[1] perpendicular separation at the bin centre # [2] parallel separation at the bin centre # [3] 2D two-point correlation function # [4] error";
276  if (m_compute_extra_info) header += " # [5] mean perpendicular separation # [6] standard deviation of the distribution of perpendicular separations # [7] mean parallel separation # [8] standard deviation of the distribution of parallel separations # [9] mean redshift # [10] standard deviation of the redshift distribution";
277  data_SS[i]->write(dir_output_resample, file, header, false, 10, 0);
278  }
279 
280  xi_SubSamples.push_back(data_SS[i]->data());
281  }
282 
283  covariance_matrix(xi_SubSamples, covariance, 1);
284 
285  if (estimator==Estimator::_natural_)
286  m_dataset = correlation_NaturalEstimator(m_dd, m_rr);
287  else if (estimator==Estimator::_LandySzalay_)
288  m_dataset = correlation_LandySzalayEstimator(m_dd, m_rr, m_dr);
289  else
290  ErrorCBL("the chosen estimator is not implemented!", "measureJackknife", "TwoPointCorrelation2D_cartesian.cpp");
291 
292  m_dataset->set_covariance(covariance);
293 }
294 
295 // ============================================================================================
296 
297 
298 void cbl::measure::twopt::TwoPointCorrelation2D_cartesian::measureBootstrap (const int nMocks, const std::string dir_output_pairs, const std::vector<std::string> dir_input_pairs, const std::string dir_output_resample, const bool count_dd, const bool count_rr, const bool count_dr, const bool tcount, const Estimator estimator, const double fact, const int seed)
299 {
300  if (nMocks<1)
301  ErrorCBL("the number of mocks must be >0", "measureBootstrap", "TwoPointCorrelation2D_cartesian.cpp");
302 
303  if (dir_output_resample!=par::defaultString && dir_output_resample!="") {
304  string mkdir = "mkdir -p "+dir_output_resample;
305  if (system(mkdir.c_str())) {}
306  }
307 
308  vector< vector<double > > xi_SubSamples, covariance;
309 
310  vector<shared_ptr<Pair> > dd_regions, rr_regions, dr_regions;
311  count_allPairs_region(dd_regions, rr_regions, dr_regions, m_twoPType, dir_output_pairs, dir_input_pairs, count_dd, count_rr, count_dr, tcount, estimator, fact);
312 
313  vector<shared_ptr<Data> > data_SS = (estimator==Estimator::_natural_) ? XiBootstrap(nMocks, dd_regions, rr_regions, seed) : XiBootstrap(nMocks, dd_regions, rr_regions, dr_regions, seed);
314 
315  for (int i=0; i<nMocks; i++) {
316 
317  if (dir_output_resample!=par::defaultString && dir_output_resample!="") {
318  string file = "xi_Bootstrap_"+conv(i, par::fINT)+".dat";
319  string header = "[1] perpendicular separation at the bin centre # [2] parallel separation at the bin centre # [3] 2D two-point correlation function # [4] error";
320  if (m_compute_extra_info) header += " # [5] mean perpendicular separation # [6] standard deviation of the distribution of perpendicular separations # [7] mean parallel separation # [8] standard deviation of the distribution of parallel separations # [9] mean redshift # [10] standard deviation of the redshift distribution";
321  data_SS[i]->write(dir_output_resample, file, header, false, 10, 0);
322  }
323 
324  xi_SubSamples.push_back(data_SS[i]->data());
325  }
326 
327  covariance_matrix(xi_SubSamples, covariance, 1);
328 
329  if (estimator==Estimator::_natural_)
330  m_dataset = correlation_NaturalEstimator(m_dd, m_rr);
331  else if (estimator==Estimator::_LandySzalay_)
332  m_dataset = correlation_LandySzalayEstimator(m_dd, m_rr, m_dr);
333  else
334  ErrorCBL("the chosen estimator is not implemented!", "measureBootstrap", "TwoPointCorrelation2D_cartesian.cpp");
335 
336  m_dataset->set_covariance(covariance);
337 }
The class TwoPointCorrelation2D_cartesian.
void measureJackknife(const std::string dir_output_pairs=par::defaultString, const std::vector< std::string > dir_input_pairs={}, const std::string dir_output_resample=par::defaultString, const bool count_dd=true, const bool count_rr=true, const bool count_dr=true, const bool tcount=true, const Estimator estimator=Estimator::_LandySzalay_, const double fact=0.1) override
measure the 2D two-point correlation function in Cartesian coordinates, estimating the covariance wit...
void measurePoisson(const std::string dir_output_pairs=par::defaultString, const std::vector< std::string > dir_input_pairs={}, const bool count_dd=true, const bool count_rr=true, const bool count_dr=true, const bool tcount=true, const Estimator estimator=Estimator::_LandySzalay_, const double fact=0.1) override
measure the 2D two-point correlation function in Cartesian coordinates, with Poisson errors
void measure(const ErrorType errorType=ErrorType::_Poisson_, const std::string dir_output_pairs=par::defaultString, const std::vector< std::string > dir_input_pairs={}, const std::string dir_output_resample=par::defaultString, const int nMocks=0, const bool count_dd=true, const bool count_rr=true, const bool count_dr=true, const bool tcount=true, const Estimator estimator=Estimator::_LandySzalay_, const double fact=0.1, const int seed=3213) override
measure the 2D two-point correlation function in Cartesian coordinates
void set_parameters(const BinType binType_rp, const double rpMin, const double rpMax, const int nbins_rp, const double shift_rp, const BinType binType_pi, const double piMin, const double piMax, const int nbins_pi, const double shift_pi, const CoordinateUnits angularUnits=CoordinateUnits::_radians_, std::function< double(double)> angularWeight=nullptr, const bool compute_extra_info=false)
set the binning parameters
void write(const std::string dir, const std::string file, const bool full, const int rank=0) const override
write the 2D two-point correlation function
void read(const std::string dir, const std::string file) override
read the 2D two-point correlation function
void measureBootstrap(const int nMocks, const std::string dir_output_pairs=par::defaultString, const std::vector< std::string > dir_input_pairs={}, const std::string dir_output_resample=par::defaultString, const bool count_dd=true, const bool count_rr=true, const bool count_dr=true, const bool tcount=true, const Estimator estimator=Estimator::_LandySzalay_, const double fact=0.1, const int seed=3213) override
measure the 2D two-point correlation function in Cartesian coordinates, estimating the covariance wit...
static const char fINT[]
conversion symbol for: int -> std::string
Definition: Constants.h:121
static const std::string defaultString
default std::string value
Definition: Constants.h:336
Estimator
the two-point correlation estimator
ErrorType
the two-point correlation function error type
Definition: Measure.h:59
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38
std::string conv(const T val, const char *fact)
convert a number to a std::string
Definition: Kernel.h:903
void checkDim(const std::vector< T > vect, const int val, const std::string vector, bool equal=true)
check if the dimension of a std::vector is equal/lower than an input value
Definition: Kernel.h:1532
int ErrorCBL(const std::string msg, const std::string functionCBL, const std::string fileCBL, const cbl::glob::ExitCode exitCode=cbl::glob::ExitCode::_error_)
throw an exception: it is used for handling exceptions inside the CosmoBolognaLib
Definition: Kernel.h:780
void covariance_matrix(const std::vector< std::vector< double >> mat, std::vector< std::vector< double >> &cov, const bool JK=false)
compute the covariance matrix from an input dataset
Definition: Func.cpp:761
CoordinateUnits
the coordinate units
Definition: Kernel.h:562
BinType
the binning type
Definition: Kernel.h:505