CLASS MANUAL
transfer.h
Go to the documentation of this file.
1 
3 #ifndef __TRANSFER__
4 #define __TRANSFER__
5 
6 #include "fourier.h"
7 #include "hyperspherical.h"
8 #include "errno.h"
9 
10 /* macro: test if index_tt is in the range between index and index+num, while the flag is true */
11 #define _index_tt_in_range_(index,num,flag) (flag == _TRUE_) && (index_tt >= index) && (index_tt < index+num)
12 /* macro: test if index_tt corresponds to an integrated nCl/sCl contribution */
13 #define _integrated_ncl_ (_index_tt_in_range_(ptr->index_tt_lensing, ppt->selection_num, ppt->has_cl_lensing_potential)) || \
14  (_index_tt_in_range_(ptr->index_tt_nc_lens, ppt->selection_num, ppt->has_nc_lens)) || \
15  (_index_tt_in_range_(ptr->index_tt_nc_g4, ppt->selection_num, ppt->has_nc_gr)) || \
16  (_index_tt_in_range_(ptr->index_tt_nc_g5, ppt->selection_num, ppt->has_nc_gr))
17 /* macro: test if index_tt corresponds to an non-integrated nCl/sCl contribution */
18 #define _nonintegrated_ncl_ (_index_tt_in_range_(ptr->index_tt_density, ppt->selection_num, ppt->has_nc_density)) || \
19  (_index_tt_in_range_(ptr->index_tt_rsd, ppt->selection_num, ppt->has_nc_rsd)) || \
20  (_index_tt_in_range_(ptr->index_tt_d0, ppt->selection_num, ppt->has_nc_rsd)) || \
21  (_index_tt_in_range_(ptr->index_tt_d1, ppt->selection_num, ppt->has_nc_rsd)) || \
22  (_index_tt_in_range_(ptr->index_tt_nc_g1, ppt->selection_num, ppt->has_nc_gr)) || \
23  (_index_tt_in_range_(ptr->index_tt_nc_g2, ppt->selection_num, ppt->has_nc_gr)) || \
24  (_index_tt_in_range_(ptr->index_tt_nc_g3, ppt->selection_num, ppt->has_nc_gr))
25 /* macro: bin number associated to particular redshift bin and selection function for non-integrated contributions*/
26 #define _get_bin_nonintegrated_ncl_(index_tt) \
27  if (_index_tt_in_range_(ptr->index_tt_density, ppt->selection_num, ppt->has_nc_density)) \
28  bin = index_tt - ptr->index_tt_density; \
29  if (_index_tt_in_range_(ptr->index_tt_rsd, ppt->selection_num, ppt->has_nc_rsd)) \
30  bin = index_tt - ptr->index_tt_rsd; \
31  if (_index_tt_in_range_(ptr->index_tt_d0, ppt->selection_num, ppt->has_nc_rsd)) \
32  bin = index_tt - ptr->index_tt_d0; \
33  if (_index_tt_in_range_(ptr->index_tt_d1, ppt->selection_num, ppt->has_nc_rsd)) \
34  bin = index_tt - ptr->index_tt_d1; \
35  if (_index_tt_in_range_(ptr->index_tt_nc_g1, ppt->selection_num, ppt->has_nc_gr)) \
36  bin = index_tt - ptr->index_tt_nc_g1; \
37  if (_index_tt_in_range_(ptr->index_tt_nc_g2, ppt->selection_num, ppt->has_nc_gr)) \
38  bin = index_tt - ptr->index_tt_nc_g2; \
39  if (_index_tt_in_range_(ptr->index_tt_nc_g3, ppt->selection_num, ppt->has_nc_gr)) \
40  bin = index_tt - ptr->index_tt_nc_g3;
41 /* macro: bin number associated to particular redshift bin and selection function for integrated contributions*/
42 #define _get_bin_integrated_ncl_(index_tt) \
43  if (_index_tt_in_range_(ptr->index_tt_lensing, ppt->selection_num, ppt->has_cl_lensing_potential)) \
44  bin = index_tt - ptr->index_tt_lensing; \
45  if (_index_tt_in_range_(ptr->index_tt_nc_lens, ppt->selection_num, ppt->has_nc_lens)) \
46  bin = index_tt - ptr->index_tt_nc_lens; \
47  if (_index_tt_in_range_(ptr->index_tt_nc_g4, ppt->selection_num, ppt->has_nc_gr)) \
48  bin = index_tt - ptr->index_tt_nc_g4; \
49  if (_index_tt_in_range_(ptr->index_tt_nc_g5, ppt->selection_num, ppt->has_nc_gr)) \
50  bin = index_tt - ptr->index_tt_nc_g5;
74 struct transfer {
75 
81 
82  double lcmb_rescale;
85  double lcmb_tilt;
88  double lcmb_pivot;
94  short has_nz_file;
96  FileName nz_file_name;
97  int nz_size;
98  double * nz_z;
99  double * nz_nz;
100  double * nz_ddnz;
104  FileName nz_evo_file_name;
106  double * nz_evo_z;
107  double * nz_evo_nz;
108  double * nz_evo_dlog_nz;
109  double * nz_evo_dd_dlog_nz;
112 
116 
117  short has_cls;
120 
124 
125  int md_size;
146  int * tt_size;
149 
153 
154  int ** l_size_tt;
156  int * l_size;
160  int * l;
162  //int * l_size_bessel; /**< for each wavenumber, maximum value of l at which bessel functions must be evaluated */
163 
167 
171 
172  size_t q_size;
174  double * q;
176  double ** k;
182  size_t q_size_limber;
184  double * q_limber;
186  double ** k_limber;
189 
193 
194  double ** transfer;
196  double ** transfer_limber;
199 
203 
206  ErrorMsg error_message;
208  short is_allocated;
211 };
212 
221 
225 
226  HyperInterpStruct HIS;
230  HyperInterpStruct * pBIS;
232  int l_size;
235 
239 
240  int tau_size;
246  double * sources;
251  double * tau0_minus_tau;
252  double * w_trapz;
253  double * chi;
257  double * cscKgen;
258  double * cotKgen;
261 
265 
266  double K;
267  int sgnK;
270 
273 };
274 
282 typedef enum {SCALAR_TEMPERATURE_0,
283  SCALAR_TEMPERATURE_1,
284  SCALAR_TEMPERATURE_2,
285  SCALAR_POLARISATION_E,
286  VECTOR_TEMPERATURE_1,
287  VECTOR_TEMPERATURE_2,
288  VECTOR_POLARISATION_E,
289  VECTOR_POLARISATION_B,
290  TENSOR_TEMPERATURE_2,
291  TENSOR_POLARISATION_E,
292  TENSOR_POLARISATION_B,
293  NC_RSD} radial_function_type;
294 
295 enum Hermite_Interpolation_Order {HERMITE3, HERMITE4, HERMITE6};
296 
297 /*************************************************************************************************************/
298 /* @cond INCLUDE_WITH_DOXYGEN */
299 /*
300  * Boilerplate for C++
301  */
302 #ifdef __cplusplus
303 extern "C" {
304 #endif
305 
307  struct transfer * ptr,
308  int index_md,
309  int index_ic,
310  int index_type,
311  int index_l,
312  double q,
313  double * ptransfer_local
314  );
315 
316  int transfer_init(
317  struct precision * ppr,
318  struct background * pba,
319  struct thermodynamics * pth,
320  struct perturbations * ppt,
321  struct fourier * pfo,
322  struct transfer * ptr
323  );
324 
325  int transfer_free(
326  struct transfer * ptr
327  );
328 
329  int transfer_indices(
330  struct precision * ppr,
331  struct perturbations * ppt,
332  struct transfer * ptr,
333  double q_period,
334  double K,
335  int sgnK
336  );
337 
338  int transfer_perturbation_copy_sources_and_nl_corrections(
339  struct perturbations * ppt,
340  struct fourier * pfo,
341  struct transfer * ptr,
342  double *** sources
343  );
344 
345  int transfer_perturbation_source_spline(
346  struct perturbations * ppt,
347  struct transfer * ptr,
348  double *** sources,
349  double *** sources_spline
350  );
351 
352  int transfer_perturbation_sources_free(
353  struct perturbations * ppt,
354  struct fourier * pfo,
355  struct transfer * ptr,
356  double *** sources
357  );
358 
359  int transfer_perturbation_sources_spline_free(
360  struct perturbations * ppt,
361  struct transfer * ptr,
362  double *** sources_spline
363  );
364 
366  struct precision * ppr,
367  struct perturbations * ppt,
368  struct transfer * ptr
369  );
370 
372  struct precision * ppr,
373  struct perturbations * ppt,
374  struct transfer * ptr,
375  double q_period,
376  double K,
377  int sgnK
378  );
379 
381  struct precision * ppr,
382  struct perturbations * ppt,
383  struct transfer * ptr,
384  double K,
385  int sgnK
386  );
387 
389  struct perturbations * ppt,
390  struct transfer * ptr,
391  double K
392  );
393 
395  struct perturbations * ppt,
396  struct transfer * ptr,
397  int ** tp_of_tt
398  );
399 
400  int transfer_free_source_correspondence(
401  struct transfer * ptr,
402  int ** tp_of_tt
403  );
404 
405  int transfer_source_tau_size_max(
406  struct precision * ppr,
407  struct background * pba,
408  struct perturbations * ppt,
409  struct transfer * ptr,
410  double tau_rec,
411  double tau0,
412  int * tau_size_max
413  );
414 
416  struct precision * ppr,
417  struct background * pba,
418  struct perturbations * ppt,
419  struct transfer * ptr,
420  double tau_rec,
421  double tau0,
422  int index_md,
423  int index_tt,
424  int * tau_size
425  );
426 
428  struct precision * ppr,
429  struct background * pba,
430  struct perturbations * ppt,
431  struct transfer * ptr,
432  int ** tp_of_tt,
433  int index_q,
434  int tau_size_max,
435  double tau_rec,
436  double *** sources,
437  double *** sources_spline,
438  double * window,
439  struct transfer_workspace * ptw,
440  short use_full_limber
441  );
442 
443  int transfer_radial_coordinates(
444  struct transfer * ptr,
445  struct transfer_workspace * ptw,
446  int index_md,
447  int index_q
448  );
449 
451  struct perturbations * ppt,
452  struct transfer * ptr,
453  double k,
454  int index_md,
455  int index_ic,
456  int index_type,
457  double * sources,
458  double * source_spline,
459  double * interpolated_sources
460  );
461 
462  int transfer_sources(
463  struct precision * ppr,
464  struct background * pba,
465  struct perturbations * ppt,
466  struct transfer * ptr,
467  double * interpolated_sources,
468  double tau_rec,
469  double k,
470  int index_md,
471  int index_tt,
472  double * sources,
473  double * window,
474  int tau_size_max,
475  double * tau0_minus_tau,
476  double * delta_tau,
477  int * tau_size_out
478  );
479 
481  struct precision * ppr,
482  struct perturbations * ppt,
483  struct transfer * ptr,
484  int bin,
485  double z,
486  double * selection);
487 
489  struct transfer * ptr,
490  double z,
491  double * dNdz,
492  double * dln_dNdz_dz);
493 
495  struct precision * ppr,
496  struct background * pba,
497  struct perturbations * ppt,
498  struct transfer * ptr,
499  int bin,
500  double * tau0_minus_tau,
501  int tau_size);
502 
504  struct precision * ppr,
505  struct background * pba,
506  struct perturbations * ppt,
507  struct transfer * ptr,
508  int bin,
509  double tau0,
510  double * tau0_minus_tau,
511  int tau_size);
512 
514  struct precision * ppr,
515  struct background * pba,
516  struct perturbations * ppt,
517  struct transfer * ptr,
518  int bin,
519  double * tau0_minus_tau,
520  int tau_size,
521  int index_md,
522  double tau0,
523  double * interpolated_sources,
524  double * sources);
525 
527  struct precision * ppr,
528  struct background * pba,
529  struct perturbations * ppt,
530  struct transfer * ptr,
531  int bin,
532  double * tau_min,
533  double * tau_mean,
534  double * tau_max);
535 
537  struct precision * ppr,
538  struct background * pba,
539  struct perturbations * ppt,
540  struct transfer * ptr,
541  double * selection,
542  double * tau0_minus_tau,
543  double * delta_tau,
544  int tau_size,
545  double * pvecback,
546  double tau0,
547  int bin);
548 
550  struct transfer_workspace * ptw,
551  struct precision * ppr,
552  struct perturbations * ppt,
553  struct transfer * ptr,
554  int index_q,
555  int index_md,
556  int index_ic,
557  int index_tt,
558  int index_l,
559  double l,
560  double q_max_bessel,
561  radial_function_type radial_type,
562  short use_full_limber
563  );
564 
565  int transfer_use_limber(
566  struct precision * ppr,
567  struct perturbations * ppt,
568  struct transfer * ptr,
569  double q_max_bessel,
570  int index_md,
571  int index_tt,
572  double q,
573  double l,
574  short * use_limber
575  );
576 
577  int transfer_integrate(
578  struct perturbations * ppt,
579  struct transfer * ptr,
580  struct transfer_workspace *ptw,
581  int index_q,
582  int index_md,
583  int index_tt,
584  double l,
585  int index_l,
586  double q,
587  radial_function_type radial_type,
588  double * trsf
589  );
590 
591  int transfer_limber(
592  struct transfer * ptr,
593  struct transfer_workspace * ptw,
594  int index_md,
595  int index_q,
596  double l,
597  double q,
598  radial_function_type radial_type,
599  double * trsf
600  );
601 
603  struct transfer * ptr,
604  double * tau0_minus_tau,
605  double * sources,
606  int tau_size,
607  double tau0_minus_tau_limber,
608  double * S
609  );
610 
611  int transfer_limber2(
612  int tau_size,
613  struct transfer * ptr,
614  int index_md,
615  int index_q,
616  double l,
617  double q,
618  double * tau0_minus_tau,
619  double * sources,
620  radial_function_type radial_type,
621  double * trsf
622  );
623 
624  int transfer_can_be_neglected(
625  struct precision * ppr,
626  struct perturbations * ppt,
627  struct transfer * ptr,
628  int index_md,
629  int index_ic,
630  int index_tt,
631  double ra_rec,
632  double q,
633  double l,
634  short * neglect
635  );
636 
637  int transfer_late_source_can_be_neglected(
638  struct precision * ppr,
639  struct perturbations * ppt,
640  struct transfer * ptr,
641  int index_md,
642  int index_tt,
643  double l,
644  short * neglect);
645 
646  int transfer_select_radial_function(
647  struct perturbations * ppt,
648  struct transfer * ptr,
649  int index_md,
650  int index_tt,
651  radial_function_type *radial_type
652  );
653 
654  int transfer_radial_function(
655  struct transfer_workspace * ptw,
656  struct perturbations * ppt,
657  struct transfer * ptr,
658  double k,
659  int index_q,
660  int index_l,
661  int x_size,
662  double * radial_function,
663  radial_function_type radial_type
664  );
665 
666  int transfer_init_HIS_from_bessel(
667  struct transfer * ptr,
668  HyperInterpStruct *pHIS
669  );
670 
671  int transfer_global_selection_read(
672  struct transfer * ptr
673  );
674 
675  int transfer_workspace_init(
676  struct transfer * ptr,
677  struct precision * ppr,
678  struct transfer_workspace *ptw,
679  int perturbations_tau_size,
680  int tau_size_max,
681  double K,
682  int sgnK,
683  double tau0_minus_tau_cut,
684  HyperInterpStruct * pBIS
685  );
686 
687  int transfer_workspace_free(
688  struct transfer * ptr,
689  struct transfer_workspace *ptw
690  );
691 
692  int transfer_update_HIS(
693  struct precision * ppr,
694  struct transfer * ptr,
695  struct transfer_workspace * ptw,
696  int index_q,
697  double tau0
698  );
699 
700  int transfer_get_lmax(int (*get_xmin_generic)(int sgnK,
701  int l,
702  double nu,
703  double xtol,
704  double phiminabs,
705  double *x_nonzero,
706  int *fevals),
707  int sgnK,
708  double nu,
709  int *lvec,
710  int lsize,
711  double phiminabs,
712  double xmax,
713  double xtol,
714  int *index_l_left,
715  int *index_l_right,
716  ErrorMsg error_message);
717 
719  struct precision * ppr,
720  struct background * pba,
721  struct perturbations * ppt,
722  struct transfer * ptr,
723  double tau_rec,
724  int tau_size_max,
725  double ** window
726  );
727 
728  int transfer_f_evo(
729  struct background* pba,
730  struct transfer * ptr,
731  double* pvecback,
732  int last_index,
733  double cotKgen,
734  double* f_evo
735  );
736 
737 #ifdef __cplusplus
738 }
739 #endif
740 
741 #endif
742 /* @endcond */
Definition: common.h:406
Definition: fourier.h:33
#define _SELECTION_NUM_MAX_
Definition: perturbations.h:70
Definition: perturbations.h:98
Definition: background.h:44
Definition: thermodynamics.h:59
int transfer_compute_for_each_q(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, int **tp_of_tt, int index_q, int tau_size_max, double tau_rec, double ***pert_sources, double ***pert_sources_spline, double *window, struct transfer_workspace *ptw, short use_full_limber)
Definition: transfer.c:1810
int transfer_limber_interpolate(struct transfer *ptr, double *tau0_minus_tau, double *sources, int tau_size, double tau0_minus_tau_limber, double *S)
Definition: transfer.c:3728
int transfer_functions_at_q(struct transfer *ptr, int index_md, int index_ic, int index_tt, int index_l, double q, double *transfer_function)
Definition: transfer.c:62
int transfer_get_source_correspondence(struct perturbations *ppt, struct transfer *ptr, int **tp_of_tt)
Definition: transfer.c:1453
int transfer_compute_for_each_l(struct transfer_workspace *ptw, struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, int index_q, int index_md, int index_ic, int index_tt, int index_l, double l, double q_max_bessel, radial_function_type radial_type, short use_full_limber)
Definition: transfer.c:3189
int transfer_get_k_list(struct perturbations *ppt, struct transfer *ptr, double K)
Definition: transfer.c:1344
int transfer_indices(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, double q_period, double K, int sgnK)
Definition: transfer.c:485
int transfer_source_tau_size(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, double tau_rec, double tau0, int index_md, int index_tt, int *tau_size)
Definition: transfer.c:1637
int transfer_integrate(struct perturbations *ppt, struct transfer *ptr, struct transfer_workspace *ptw, int index_q, int index_md, int index_tt, double l, int index_l, double k, radial_function_type radial_type, double *trsf)
Definition: transfer.c:3390
int transfer_selection_sampling(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, int bin, double *tau0_minus_tau, int tau_size)
Definition: transfer.c:2796
int transfer_get_q_list(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, double q_period, double K, int sgnK)
Definition: transfer.c:1032
int transfer_interpolate_sources(struct perturbations *ppt, struct transfer *ptr, double k, int index_md, int index_ic, int index_type, double *pert_source, double *pert_source_spline, double *interpolated_sources)
Definition: transfer.c:2191
int transfer_dNdz_analytic(struct transfer *ptr, double z, double *dNdz, double *dln_dNdz_dz)
Definition: transfer.c:2752
int transfer_limber2(int tau_size, struct transfer *ptr, int index_md, int index_k, double l, double k, double *tau0_minus_tau, double *sources, radial_function_type radial_type, double *trsf)
Definition: transfer.c:3815
int transfer_selection_times(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, int bin, double *tau_min, double *tau_mean, double *tau_max)
Definition: transfer.c:2988
int transfer_init(struct precision *ppr, struct background *pba, struct thermodynamics *pth, struct perturbations *ppt, struct fourier *pfo, struct transfer *ptr)
Definition: transfer.c:116
int transfer_limber(struct transfer *ptr, struct transfer_workspace *ptw, int index_md, int index_q, double l, double q, radial_function_type radial_type, double *trsf)
Definition: transfer.c:3571
int transfer_get_q_limber_list(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, double K, int sgnK)
Definition: transfer.c:1261
int transfer_lensing_sampling(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, int bin, double tau0, double *tau0_minus_tau, int tau_size)
Definition: transfer.c:2864
int transfer_source_resample(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, int bin, double *tau0_minus_tau, int tau_size, int index_md, double tau0, double *interpolated_sources, double *sources)
Definition: transfer.c:2921
int transfer_precompute_selection(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, double tau_rec, int tau_size_max, double **window)
Definition: transfer.c:4793
int transfer_selection_function(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, int bin, double z, double *selection)
Definition: transfer.c:2607
int transfer_free(struct transfer *ptr)
Definition: transfer.c:417
int transfer_get_l_list(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr)
Definition: transfer.c:833
int transfer_sources(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, double *interpolated_sources, double tau_rec, double k, int index_md, int index_tt, double *sources, double *window, int tau_size_max, double *tau0_minus_tau, double *w_trapz, int *tau_size_out)
Definition: transfer.c:2276
int transfer_selection_compute(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, double *selection, double *tau0_minus_tau, double *w_trapz, int tau_size, double *pvecback, double tau0, int bin)
Definition: transfer.c:3068
double angular_rescaling
Definition: transfer.h:164
double ** k_limber
Definition: transfer.h:186
int * l_size
Definition: transfer.h:156
double * cscKgen
Definition: transfer.h:257
double selection_magnification_bias[_SELECTION_NUM_MAX_]
Definition: transfer.h:92
double tau0_minus_tau_cut
Definition: transfer.h:271
double * nz_evo_nz
Definition: transfer.h:107
int index_tt_d1
Definition: transfer.h:138
double * nz_evo_z
Definition: transfer.h:106
double K
Definition: transfer.h:266
double * nz_z
Definition: transfer.h:98
short has_nz_file
Definition: transfer.h:94
int index_tt_nc_g3
Definition: transfer.h:142
double * chi
Definition: transfer.h:253
double * w_trapz
Definition: transfer.h:252
double ** transfer_limber
Definition: transfer.h:196
int index_tt_lcmb
Definition: transfer.h:132
short transfer_verbose
Definition: transfer.h:204
int nz_size
Definition: transfer.h:97
short has_nz_evo_analytic
Definition: transfer.h:103
double * nz_evo_dd_dlog_nz
Definition: transfer.h:109
int * l
Definition: transfer.h:160
int nz_evo_size
Definition: transfer.h:105
short has_nz_analytic
Definition: transfer.h:95
int index_tt_nc_g1
Definition: transfer.h:140
int index_tt_e
Definition: transfer.h:130
int tau_size
Definition: transfer.h:240
double lcmb_tilt
Definition: transfer.h:85
double * nz_evo_dlog_nz
Definition: transfer.h:108
double ** k
Definition: transfer.h:176
double selection_bias[_SELECTION_NUM_MAX_]
Definition: transfer.h:91
double * sources
Definition: transfer.h:246
FileName nz_evo_file_name
Definition: transfer.h:104
short do_lcmb_full_limber
Definition: transfer.h:180
int l_size
Definition: transfer.h:232
int ** l_size_tt
Definition: transfer.h:154
double * nz_nz
Definition: transfer.h:99
short is_allocated
Definition: transfer.h:208
int index_tt_nc_lens
Definition: transfer.h:139
int index_q_flat_approximation
Definition: transfer.h:178
int index_tt_nc_g4
Definition: transfer.h:143
size_t q_size_limber
Definition: transfer.h:182
HyperInterpStruct HIS
Definition: transfer.h:226
int HIS_allocated
Definition: transfer.h:228
ErrorMsg error_message
Definition: transfer.h:206
int index_tt_t0
Definition: transfer.h:127
int index_tt_density
Definition: transfer.h:133
double * q
Definition: transfer.h:174
HyperInterpStruct * pBIS
Definition: transfer.h:230
double * tau0_minus_tau
Definition: transfer.h:251
int index_tt_b
Definition: transfer.h:131
int md_size
Definition: transfer.h:125
int index_tt_lensing
Definition: transfer.h:134
short neglect_late_source
Definition: transfer.h:272
double * cotKgen
Definition: transfer.h:258
double * q_limber
Definition: transfer.h:184
int * tt_size
Definition: transfer.h:146
int index_tt_nc_g2
Definition: transfer.h:141
int index_tt_t1
Definition: transfer.h:128
size_t q_size
Definition: transfer.h:172
int index_tt_d0
Definition: transfer.h:137
double * interpolated_sources
Definition: transfer.h:242
FileName nz_file_name
Definition: transfer.h:96
double lcmb_rescale
Definition: transfer.h:82
int index_tt_nc_g5
Definition: transfer.h:144
short has_cls
Definition: transfer.h:117
int tau_size_max
Definition: transfer.h:241
short has_nz_evo_file
Definition: transfer.h:102
double ** transfer
Definition: transfer.h:194
int index_tt_t2
Definition: transfer.h:129
radial_function_type
Definition: transfer.h:282
int sgnK
Definition: transfer.h:267
int index_tt_rsd
Definition: transfer.h:136
double lcmb_pivot
Definition: transfer.h:88
int l_size_max
Definition: transfer.h:158
double * nz_ddnz
Definition: transfer.h:100
Definition: transfer.h:74
Definition: transfer.h:220