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 <sys/shm.h>
9 #include <sys/stat.h>
10 #include "errno.h"
11 
12 /* macro: test if index_tt is in the range between index and index+num, while the flag is true */
13 #define _index_tt_in_range_(index,num,flag) (flag == _TRUE_) && (index_tt >= index) && (index_tt < index+num)
14 /* macro: test if index_tt corresponds to an integrated nCl/sCl contribution */
15 #define _integrated_ncl_ (_index_tt_in_range_(ptr->index_tt_lensing, ppt->selection_num, ppt->has_cl_lensing_potential)) || \
16  (_index_tt_in_range_(ptr->index_tt_nc_lens, ppt->selection_num, ppt->has_nc_lens)) || \
17  (_index_tt_in_range_(ptr->index_tt_nc_g4, ppt->selection_num, ppt->has_nc_gr)) || \
18  (_index_tt_in_range_(ptr->index_tt_nc_g5, ppt->selection_num, ppt->has_nc_gr))
19 /* macro: test if index_tt corresponds to an non-integrated nCl/sCl contribution */
20 #define _nonintegrated_ncl_ (_index_tt_in_range_(ptr->index_tt_density, ppt->selection_num, ppt->has_nc_density)) || \
21  (_index_tt_in_range_(ptr->index_tt_rsd, ppt->selection_num, ppt->has_nc_rsd)) || \
22  (_index_tt_in_range_(ptr->index_tt_d0, ppt->selection_num, ppt->has_nc_rsd)) || \
23  (_index_tt_in_range_(ptr->index_tt_d1, ppt->selection_num, ppt->has_nc_rsd)) || \
24  (_index_tt_in_range_(ptr->index_tt_nc_g1, ppt->selection_num, ppt->has_nc_gr)) || \
25  (_index_tt_in_range_(ptr->index_tt_nc_g2, ppt->selection_num, ppt->has_nc_gr)) || \
26  (_index_tt_in_range_(ptr->index_tt_nc_g3, ppt->selection_num, ppt->has_nc_gr))
27 /* macro: bin number associated to particular redshift bin and selection function for non-integrated contributions*/
28 #define _get_bin_nonintegrated_ncl_(index_tt) \
29  if (_index_tt_in_range_(ptr->index_tt_density, ppt->selection_num, ppt->has_nc_density)) \
30  bin = index_tt - ptr->index_tt_density; \
31  if (_index_tt_in_range_(ptr->index_tt_rsd, ppt->selection_num, ppt->has_nc_rsd)) \
32  bin = index_tt - ptr->index_tt_rsd; \
33  if (_index_tt_in_range_(ptr->index_tt_d0, ppt->selection_num, ppt->has_nc_rsd)) \
34  bin = index_tt - ptr->index_tt_d0; \
35  if (_index_tt_in_range_(ptr->index_tt_d1, ppt->selection_num, ppt->has_nc_rsd)) \
36  bin = index_tt - ptr->index_tt_d1; \
37  if (_index_tt_in_range_(ptr->index_tt_nc_g1, ppt->selection_num, ppt->has_nc_gr)) \
38  bin = index_tt - ptr->index_tt_nc_g1; \
39  if (_index_tt_in_range_(ptr->index_tt_nc_g2, ppt->selection_num, ppt->has_nc_gr)) \
40  bin = index_tt - ptr->index_tt_nc_g2; \
41  if (_index_tt_in_range_(ptr->index_tt_nc_g3, ppt->selection_num, ppt->has_nc_gr)) \
42  bin = index_tt - ptr->index_tt_nc_g3;
43 /* macro: bin number associated to particular redshift bin and selection function for integrated contributions*/
44 #define _get_bin_integrated_ncl_(index_tt) \
45  if (_index_tt_in_range_(ptr->index_tt_lensing, ppt->selection_num, ppt->has_cl_lensing_potential)) \
46  bin = index_tt - ptr->index_tt_lensing; \
47  if (_index_tt_in_range_(ptr->index_tt_nc_lens, ppt->selection_num, ppt->has_nc_lens)) \
48  bin = index_tt - ptr->index_tt_nc_lens; \
49  if (_index_tt_in_range_(ptr->index_tt_nc_g4, ppt->selection_num, ppt->has_nc_gr)) \
50  bin = index_tt - ptr->index_tt_nc_g4; \
51  if (_index_tt_in_range_(ptr->index_tt_nc_g5, ppt->selection_num, ppt->has_nc_gr)) \
52  bin = index_tt - ptr->index_tt_nc_g5;
76 struct transfer {
77 
83 
84  double lcmb_rescale;
87  double lcmb_tilt;
90  double lcmb_pivot;
96  short has_nz_file;
98  FileName nz_file_name;
99  int nz_size;
100  double * nz_z;
101  double * nz_nz;
102  double * nz_ddnz;
106  FileName nz_evo_file_name;
108  double * nz_evo_z;
109  double * nz_evo_nz;
110  double * nz_evo_dlog_nz;
111  double * nz_evo_dd_dlog_nz;
114 
118 
119  short has_cls;
122 
126 
127  int md_size;
148  int * tt_size;
151 
155 
156  int ** l_size_tt;
158  int * l_size;
162  int * l;
164  //int * l_size_bessel; /**< for each wavenumber, maximum value of l at which bessel functions must be evaluated */
165 
169 
173 
174  size_t q_size;
176  double * q;
178  double ** k;
183 
187 
188  double ** transfer;
191 
195 
198  ErrorMsg error_message;
201 };
202 
211 
215 
216  HyperInterpStruct HIS;
220  HyperInterpStruct * pBIS;
222  int l_size;
225 
229 
230  int tau_size;
236  double * sources;
241  double * tau0_minus_tau;
242  double * w_trapz;
243  double * chi;
247  double * cscKgen;
248  double * cotKgen;
251 
255 
256  double K;
257  int sgnK;
260 
263 };
264 
272 typedef enum {SCALAR_TEMPERATURE_0,
273  SCALAR_TEMPERATURE_1,
274  SCALAR_TEMPERATURE_2,
275  SCALAR_POLARISATION_E,
276  VECTOR_TEMPERATURE_1,
277  VECTOR_TEMPERATURE_2,
278  VECTOR_POLARISATION_E,
279  VECTOR_POLARISATION_B,
280  TENSOR_TEMPERATURE_2,
281  TENSOR_POLARISATION_E,
282  TENSOR_POLARISATION_B,
283  NC_RSD} radial_function_type;
284 
285 enum Hermite_Interpolation_Order {HERMITE3, HERMITE4, HERMITE6};
286 
287 /*************************************************************************************************************/
288 /* @cond INCLUDE_WITH_DOXYGEN */
289 /*
290  * Boilerplate for C++
291  */
292 #ifdef __cplusplus
293 extern "C" {
294 #endif
295 
297  struct transfer * ptr,
298  int index_md,
299  int index_ic,
300  int index_type,
301  int index_l,
302  double q,
303  double * ptransfer_local
304  );
305 
306  int transfer_init(
307  struct precision * ppr,
308  struct background * pba,
309  struct thermodynamics * pth,
310  struct perturbations * ppt,
311  struct fourier * pfo,
312  struct transfer * ptr
313  );
314 
315  int transfer_free(
316  struct transfer * ptr
317  );
318 
319  int transfer_indices(
320  struct precision * ppr,
321  struct perturbations * ppt,
322  struct transfer * ptr,
323  double q_period,
324  double K,
325  int sgnK
326  );
327 
328  int transfer_perturbation_copy_sources_and_nl_corrections(
329  struct perturbations * ppt,
330  struct fourier * pfo,
331  struct transfer * ptr,
332  double *** sources
333  );
334 
335  int transfer_perturbation_source_spline(
336  struct perturbations * ppt,
337  struct transfer * ptr,
338  double *** sources,
339  double *** sources_spline
340  );
341 
342  int transfer_perturbation_sources_free(
343  struct perturbations * ppt,
344  struct fourier * pfo,
345  struct transfer * ptr,
346  double *** sources
347  );
348 
349  int transfer_perturbation_sources_spline_free(
350  struct perturbations * ppt,
351  struct transfer * ptr,
352  double *** sources_spline
353  );
354 
356  struct precision * ppr,
357  struct perturbations * ppt,
358  struct transfer * ptr
359  );
360 
362  struct precision * ppr,
363  struct perturbations * ppt,
364  struct transfer * ptr,
365  double q_period,
366  double K,
367  int sgnK
368  );
369 
370  int transfer_get_q_list_v1(
371  struct precision * ppr,
372  struct perturbations * ppt,
373  struct transfer * ptr,
374  double q_period,
375  double K,
376  int sgnK
377  );
378 
380  struct perturbations * ppt,
381  struct transfer * ptr,
382  double K
383  );
384 
386  struct perturbations * ppt,
387  struct transfer * ptr,
388  int ** tp_of_tt
389  );
390 
391  int transfer_free_source_correspondence(
392  struct transfer * ptr,
393  int ** tp_of_tt
394  );
395 
396  int transfer_source_tau_size_max(
397  struct precision * ppr,
398  struct background * pba,
399  struct perturbations * ppt,
400  struct transfer * ptr,
401  double tau_rec,
402  double tau0,
403  int * tau_size_max
404  );
405 
407  struct precision * ppr,
408  struct background * pba,
409  struct perturbations * ppt,
410  struct transfer * ptr,
411  double tau_rec,
412  double tau0,
413  int index_md,
414  int index_tt,
415  int * tau_size
416  );
417 
419  struct precision * ppr,
420  struct background * pba,
421  struct perturbations * ppt,
422  struct transfer * ptr,
423  int ** tp_of_tt,
424  int index_q,
425  int tau_size_max,
426  double tau_rec,
427  double *** sources,
428  double *** sources_spline,
429  double * window,
430  struct transfer_workspace * ptw
431  );
432 
433  int transfer_radial_coordinates(
434  struct transfer * ptr,
435  struct transfer_workspace * ptw,
436  int index_md,
437  int index_q
438  );
439 
441  struct perturbations * ppt,
442  struct transfer * ptr,
443  int index_q,
444  int index_md,
445  int index_ic,
446  int index_type,
447  double * sources,
448  double * source_spline,
449  double * interpolated_sources
450  );
451 
452  int transfer_sources(
453  struct precision * ppr,
454  struct background * pba,
455  struct perturbations * ppt,
456  struct transfer * ptr,
457  double * interpolated_sources,
458  double tau_rec,
459  int index_q,
460  int index_md,
461  int index_tt,
462  double * sources,
463  double * window,
464  int tau_size_max,
465  double * tau0_minus_tau,
466  double * delta_tau,
467  int * tau_size_out
468  );
469 
471  struct precision * ppr,
472  struct perturbations * ppt,
473  struct transfer * ptr,
474  int bin,
475  double z,
476  double * selection);
477 
479  struct transfer * ptr,
480  double z,
481  double * dNdz,
482  double * dln_dNdz_dz);
483 
485  struct precision * ppr,
486  struct background * pba,
487  struct perturbations * ppt,
488  struct transfer * ptr,
489  int bin,
490  double * tau0_minus_tau,
491  int tau_size);
492 
494  struct precision * ppr,
495  struct background * pba,
496  struct perturbations * ppt,
497  struct transfer * ptr,
498  int bin,
499  double tau0,
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_minus_tau,
510  int tau_size,
511  int index_md,
512  double tau0,
513  double * interpolated_sources,
514  double * sources);
515 
517  struct precision * ppr,
518  struct background * pba,
519  struct perturbations * ppt,
520  struct transfer * ptr,
521  int bin,
522  double * tau_min,
523  double * tau_mean,
524  double * tau_max);
525 
527  struct precision * ppr,
528  struct background * pba,
529  struct perturbations * ppt,
530  struct transfer * ptr,
531  double * selection,
532  double * tau0_minus_tau,
533  double * delta_tau,
534  int tau_size,
535  double * pvecback,
536  double tau0,
537  int bin);
538 
540  struct transfer_workspace * ptw,
541  struct precision * ppr,
542  struct perturbations * ppt,
543  struct transfer * ptr,
544  int index_q,
545  int index_md,
546  int index_ic,
547  int index_tt,
548  int index_l,
549  double l,
550  double q_max_bessel,
551  radial_function_type radial_type
552  );
553 
554  int transfer_use_limber(
555  struct precision * ppr,
556  struct perturbations * ppt,
557  struct transfer * ptr,
558  double q_max_bessel,
559  int index_md,
560  int index_tt,
561  double q,
562  double l,
563  short * use_limber
564  );
565 
566  int transfer_integrate(
567  struct perturbations * ppt,
568  struct transfer * ptr,
569  struct transfer_workspace *ptw,
570  int index_q,
571  int index_md,
572  int index_tt,
573  double l,
574  int index_l,
575  double q,
576  radial_function_type radial_type,
577  double * trsf
578  );
579 
580  int transfer_limber(
581  struct transfer * ptr,
582  struct transfer_workspace * ptw,
583  int index_md,
584  int index_q,
585  double l,
586  double q,
587  radial_function_type radial_type,
588  double * trsf
589  );
590 
592  struct transfer * ptr,
593  double * tau0_minus_tau,
594  double * sources,
595  int tau_size,
596  double tau0_minus_tau_limber,
597  double * S
598  );
599 
600  int transfer_limber2(
601  int tau_size,
602  struct transfer * ptr,
603  int index_md,
604  int index_q,
605  double l,
606  double q,
607  double * tau0_minus_tau,
608  double * sources,
609  radial_function_type radial_type,
610  double * trsf
611  );
612 
613  int transfer_can_be_neglected(
614  struct precision * ppr,
615  struct perturbations * ppt,
616  struct transfer * ptr,
617  int index_md,
618  int index_ic,
619  int index_tt,
620  double ra_rec,
621  double q,
622  double l,
623  short * neglect
624  );
625 
626  int transfer_late_source_can_be_neglected(
627  struct precision * ppr,
628  struct perturbations * ppt,
629  struct transfer * ptr,
630  int index_md,
631  int index_tt,
632  double l,
633  short * neglect);
634 
635  int transfer_select_radial_function(
636  struct perturbations * ppt,
637  struct transfer * ptr,
638  int index_md,
639  int index_tt,
640  radial_function_type *radial_type
641  );
642 
643  int transfer_radial_function(
644  struct transfer_workspace * ptw,
645  struct perturbations * ppt,
646  struct transfer * ptr,
647  double k,
648  int index_q,
649  int index_l,
650  int x_size,
651  double * radial_function,
652  radial_function_type radial_type
653  );
654 
655  int transfer_init_HIS_from_bessel(
656  struct transfer * ptr,
657  HyperInterpStruct *pHIS
658  );
659 
660  int transfer_global_selection_read(
661  struct transfer * ptr
662  );
663 
664  int transfer_workspace_init(
665  struct transfer * ptr,
666  struct precision * ppr,
667  struct transfer_workspace **ptw,
668  int perturbations_tau_size,
669  int tau_size_max,
670  double K,
671  int sgnK,
672  double tau0_minus_tau_cut,
673  HyperInterpStruct * pBIS
674  );
675 
676  int transfer_workspace_free(
677  struct transfer * ptr,
678  struct transfer_workspace *ptw
679  );
680 
681  int transfer_update_HIS(
682  struct precision * ppr,
683  struct transfer * ptr,
684  struct transfer_workspace * ptw,
685  int index_q,
686  double tau0
687  );
688 
689  int transfer_get_lmax(int (*get_xmin_generic)(int sgnK,
690  int l,
691  double nu,
692  double xtol,
693  double phiminabs,
694  double *x_nonzero,
695  int *fevals),
696  int sgnK,
697  double nu,
698  int *lvec,
699  int lsize,
700  double phiminabs,
701  double xmax,
702  double xtol,
703  int *index_l_left,
704  int *index_l_right,
705  ErrorMsg error_message);
706 
708  struct precision * ppr,
709  struct background * pba,
710  struct perturbations * ppt,
711  struct transfer * ptr,
712  double tau_rec,
713  int tau_size_max,
714  double ** window
715  );
716 
717  int transfer_f_evo(
718  struct background* pba,
719  struct transfer * ptr,
720  double* pvecback,
721  int last_index,
722  double cotKgen,
723  double* f_evo
724  );
725 
726 #ifdef __cplusplus
727 }
728 #endif
729 
730 #endif
731 /* @endcond */
Definition: background.h:48
Definition: common.h:378
Definition: fourier.h:33
#define _SELECTION_NUM_MAX_
Definition: perturbations.h:70
Definition: perturbations.h:98
Definition: thermodynamics.h:59
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:3523
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)
Definition: transfer.c:2999
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:61
int transfer_get_source_correspondence(struct perturbations *ppt, struct transfer *ptr, int **tp_of_tt)
Definition: transfer.c:1326
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)
Definition: transfer.c:1683
int transfer_get_k_list(struct perturbations *ppt, struct transfer *ptr, double K)
Definition: transfer.c:1249
int transfer_indices(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, double q_period, double K, int sgnK)
Definition: transfer.c:492
int transfer_interpolate_sources(struct perturbations *ppt, struct transfer *ptr, int index_q, int index_md, int index_ic, int index_type, double *pert_source, double *pert_source_spline, double *interpolated_sources)
Definition: transfer.c:2001
int transfer_sources(struct precision *ppr, struct background *pba, struct perturbations *ppt, struct transfer *ptr, double *interpolated_sources, double tau_rec, int index_q, 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:2087
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:1510
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:3185
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:2607
int transfer_get_q_list(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, double q_period, double K, int sgnK)
Definition: transfer.c:1018
int transfer_dNdz_analytic(struct transfer *ptr, double z, double *dNdz, double *dln_dNdz_dz)
Definition: transfer.c:2563
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:3610
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:2799
int transfer_init(struct precision *ppr, struct background *pba, struct thermodynamics *pth, struct perturbations *ppt, struct fourier *pfo, struct transfer *ptr)
Definition: transfer.c:115
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:3366
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:2675
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:2732
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:4600
int transfer_selection_function(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr, int bin, double z, double *selection)
Definition: transfer.c:2418
int transfer_free(struct transfer *ptr)
Definition: transfer.c:435
int transfer_get_l_list(struct precision *ppr, struct perturbations *ppt, struct transfer *ptr)
Definition: transfer.c:819
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:2879
double angular_rescaling
Definition: transfer.h:166
int * l_size
Definition: transfer.h:158
double * cscKgen
Definition: transfer.h:247
double selection_magnification_bias[_SELECTION_NUM_MAX_]
Definition: transfer.h:94
double tau0_minus_tau_cut
Definition: transfer.h:261
double * nz_evo_nz
Definition: transfer.h:109
int index_tt_d1
Definition: transfer.h:140
double * nz_evo_z
Definition: transfer.h:108
double K
Definition: transfer.h:256
double * nz_z
Definition: transfer.h:100
short has_nz_file
Definition: transfer.h:96
int index_tt_nc_g3
Definition: transfer.h:144
double * chi
Definition: transfer.h:243
double * w_trapz
Definition: transfer.h:242
int index_tt_lcmb
Definition: transfer.h:134
short transfer_verbose
Definition: transfer.h:196
int nz_size
Definition: transfer.h:99
short has_nz_evo_analytic
Definition: transfer.h:105
double * nz_evo_dd_dlog_nz
Definition: transfer.h:111
int * l
Definition: transfer.h:162
int nz_evo_size
Definition: transfer.h:107
short has_nz_analytic
Definition: transfer.h:97
int index_tt_nc_g1
Definition: transfer.h:142
int index_tt_e
Definition: transfer.h:132
int tau_size
Definition: transfer.h:230
double lcmb_tilt
Definition: transfer.h:87
double * nz_evo_dlog_nz
Definition: transfer.h:110
double ** k
Definition: transfer.h:178
double selection_bias[_SELECTION_NUM_MAX_]
Definition: transfer.h:93
double * sources
Definition: transfer.h:236
FileName nz_evo_file_name
Definition: transfer.h:106
int l_size
Definition: transfer.h:222
int ** l_size_tt
Definition: transfer.h:156
double * nz_nz
Definition: transfer.h:101
int index_tt_nc_lens
Definition: transfer.h:141
int index_q_flat_approximation
Definition: transfer.h:180
int index_tt_nc_g4
Definition: transfer.h:145
HyperInterpStruct HIS
Definition: transfer.h:216
int HIS_allocated
Definition: transfer.h:218
ErrorMsg error_message
Definition: transfer.h:198
int index_tt_t0
Definition: transfer.h:129
int index_tt_density
Definition: transfer.h:135
double * q
Definition: transfer.h:176
HyperInterpStruct * pBIS
Definition: transfer.h:220
double * tau0_minus_tau
Definition: transfer.h:241
int index_tt_b
Definition: transfer.h:133
int md_size
Definition: transfer.h:127
int index_tt_lensing
Definition: transfer.h:136
short neglect_late_source
Definition: transfer.h:262
double * cotKgen
Definition: transfer.h:248
int * tt_size
Definition: transfer.h:148
int index_tt_nc_g2
Definition: transfer.h:143
int index_tt_t1
Definition: transfer.h:130
size_t q_size
Definition: transfer.h:174
int index_tt_d0
Definition: transfer.h:139
double * interpolated_sources
Definition: transfer.h:232
FileName nz_file_name
Definition: transfer.h:98
double lcmb_rescale
Definition: transfer.h:84
int index_tt_nc_g5
Definition: transfer.h:146
short has_cls
Definition: transfer.h:119
int tau_size_max
Definition: transfer.h:231
short has_nz_evo_file
Definition: transfer.h:104
double ** transfer
Definition: transfer.h:188
int index_tt_t2
Definition: transfer.h:131
radial_function_type
Definition: transfer.h:272
int sgnK
Definition: transfer.h:257
int index_tt_rsd
Definition: transfer.h:138
double lcmb_pivot
Definition: transfer.h:90
int l_size_max
Definition: transfer.h:160
double * nz_ddnz
Definition: transfer.h:102
Definition: transfer.h:76
Definition: transfer.h:210