CLASS MANUAL
lensing.h
Go to the documentation of this file.
1 
3 #ifndef __LENSING__
4 #define __LENSING__
5 
6 #include "harmonic.h"
7 
17 struct lensing {
18 
25 
29 
33 
34  int has_tt;
35  int has_ee;
36  int has_te;
37  int has_bb;
38  int has_pp;
39  int has_tp;
40  int has_dd;
41  int has_td;
42  int has_ll;
43  int has_tl;
56  int lt_size;
59 
63 
68  /* interpolable version: */
69 
70  int l_size;
72  int * l_max_lt;
75  double * l;
76  double * cl_lens;
80  double * ddcl_lens;
83 
87 
90  ErrorMsg error_message;
93 };
94 
95 /*************************************************************************************************************/
96 /* @cond INCLUDE_WITH_DOXYGEN */
97 /*
98  * Boilerplate for C++
99  */
100 #ifdef __cplusplus
101 extern "C" {
102 #endif
103 
104  int lensing_cl_at_l(
105  struct lensing * ple,
106  int l,
107  double * cl_lensed
108  );
109 
110  int lensing_init(
111  struct precision * ppr,
112  struct perturbations * ppt,
113  struct harmonic * phr,
114  struct fourier * pfo,
115  struct lensing * ple
116  );
117 
118  int lensing_free(
119  struct lensing * ple
120  );
121 
122  int lensing_indices(
123  struct precision * ppr,
124  struct harmonic * phr,
125  struct lensing * ple
126  );
127 
129  double *ksi,
130  double **d00,
131  double *w8,
132  int nmu,
133  struct lensing * ple
134  );
135 
137  double *ksiX,
138  double **d20,
139  double *w8,
140  int nmu,
141  struct lensing * ple
142  );
143 
145  double *ksip,
146  double *ksim,
147  double **d22,
148  double **d2m2,
149  double *w8,
150  int nmu,
151  struct lensing * ple
152  );
154  struct lensing *ple,
155  double *cl_tt
156  );
157 
159  struct lensing *ple,
160  double *cl_te
161  );
162 
164  struct lensing *ple,
165  double *cl_ee,
166  double *cl_bb
167  );
168 
169 
170  int lensing_X000(
171  double * mu,
172  int num_mu,
173  int lmax,
174  double * sigma2,
175  double ** X000
176  );
177 
178  int lensing_Xp000(
179  double * mu,
180  int num_mu,
181  int lmax,
182  double * sigma2,
183  double ** Xp000
184  );
185 
186  int lensing_X220(
187  double * mu,
188  int num_mu,
189  int lmax,
190  double * sigma2,
191  double ** X220
192  );
193 
194  int lensing_X022(
195  double * mu,
196  int num_mu,
197  int lmax,
198  double * sigma2,
199  double ** X022
200  );
201 
202  int lensing_Xp022(
203  double * mu,
204  int num_mu,
205  int lmax,
206  double * sigma2,
207  double ** Xp022
208  );
209 
210  int lensing_X121(
211  double * mu,
212  int num_mu,
213  int lmax,
214  double * sigma2,
215  double ** X121
216  );
217 
218  int lensing_X132(
219  double * mu,
220  int num_mu,
221  int lmax,
222  double * sigma2,
223  double ** X132
224  );
225 
226  int lensing_X242(
227  double * mu,
228  int num_mu,
229  int lmax,
230  double * sigma2,
231  double ** X242
232  );
233 
234  int lensing_d00(
235  double * mu,
236  int num_mu,
237  int lmax,
238  double ** d00
239  );
240 
241  int lensing_d11(
242  double * mu,
243  int num_mu,
244  int lmax,
245  double ** d11
246  );
247 
248  int lensing_d1m1(
249  double * mu,
250  int num_mu,
251  int lmax,
252  double ** d1m1
253  );
254 
255  int lensing_d2m2(
256  double * mu,
257  int num_mu,
258  int lmax,
259  double ** d2m2
260  );
261 
262  int lensing_d22(
263  double * mu,
264  int num_mu,
265  int lmax,
266  double ** d22
267  );
268 
269  int lensing_d20(
270  double * mu,
271  int num_mu,
272  int lmax,
273  double ** d20
274  );
275 
276  int lensing_d31(
277  double * mu,
278  int num_mu,
279  int lmax,
280  double ** d3m1
281  );
282 
283  int lensing_d3m1(
284  double * mu,
285  int num_mu,
286  int lmax,
287  double ** d3m1
288  );
289 
290  int lensing_d3m3(
291  double * mu,
292  int num_mu,
293  int lmax,
294  double ** d3m3
295  );
296 
297  int lensing_d40(
298  double * mu,
299  int num_mu,
300  int lmax,
301  double ** d40
302  );
303 
304  int lensing_d4m2(
305  double * mu,
306  int num_mu,
307  int lmax,
308  double ** d4m2
309  );
310 
311  int lensing_d4m4(
312  double * mu,
313  int num_mu,
314  int lmax,
315  double ** d4m4
316  );
317 
318 #ifdef __cplusplus
319 }
320 #endif
321 
322 #endif
323 /* @endcond */
Definition: common.h:378
Definition: fourier.h:33
Definition: harmonic.h:17
int lensing_indices(struct precision *ppr, struct harmonic *phr, struct lensing *ple)
Definition: lensing.c:833
int lensing_d1m1(double *mu, int num_mu, int lmax, double **d1m1)
Definition: lensing.c:1360
int lensing_d40(double *mu, int num_mu, int lmax, double **d40)
Definition: lensing.c:1760
int lensing_d20(double *mu, int num_mu, int lmax, double **d20)
Definition: lensing.c:1531
int lensing_lensed_cl_te(double *ksiX, double **d20, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1109
int lensing_addback_cl_te(struct lensing *ple, double *cl_te)
Definition: lensing.c:1147
int lensing_lensed_cl_tt(double *ksi, double **d00, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1046
int lensing_d2m2(double *mu, int num_mu, int lmax, double **d2m2)
Definition: lensing.c:1417
int lensing_d4m2(double *mu, int num_mu, int lmax, double **d4m2)
Definition: lensing.c:1817
int lensing_addback_cl_ee_bb(struct lensing *ple, double *cl_ee, double *cl_bb)
Definition: lensing.c:1217
int lensing_d11(double *mu, int num_mu, int lmax, double **d11)
Definition: lensing.c:1303
int lensing_d00(double *mu, int num_mu, int lmax, double **d00)
Definition: lensing.c:1246
int lensing_addback_cl_tt(struct lensing *ple, double *cl_tt)
Definition: lensing.c:1084
int lensing_d31(double *mu, int num_mu, int lmax, double **d31)
Definition: lensing.c:1586
int lensing_free(struct lensing *ple)
Definition: lensing.c:807
int lensing_d3m1(double *mu, int num_mu, int lmax, double **d3m1)
Definition: lensing.c:1644
int lensing_d3m3(double *mu, int num_mu, int lmax, double **d3m3)
Definition: lensing.c:1702
int lensing_d22(double *mu, int num_mu, int lmax, double **d22)
Definition: lensing.c:1474
int lensing_cl_at_l(struct lensing *ple, int l, double *cl_lensed)
Definition: lensing.c:39
int lensing_lensed_cl_ee_bb(double *ksip, double *ksim, double **d22, double **d2m2, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1174
int lensing_init(struct precision *ppr, struct perturbations *ppt, struct harmonic *phr, struct fourier *pfo, struct lensing *ple)
Definition: lensing.c:84
int lensing_d4m4(double *mu, int num_mu, int lmax, double **d4m4)
Definition: lensing.c:1876
int has_td
Definition: lensing.h:41
int has_te
Definition: lensing.h:36
int has_tt
Definition: lensing.h:34
int l_size
Definition: lensing.h:70
int index_lt_dd
Definition: lensing.h:51
int index_lt_td
Definition: lensing.h:52
int has_bb
Definition: lensing.h:37
int index_lt_te
Definition: lensing.h:47
ErrorMsg error_message
Definition: lensing.h:90
int index_lt_tt
Definition: lensing.h:45
double * l
Definition: lensing.h:75
int has_ee
Definition: lensing.h:35
short lensing_verbose
Definition: lensing.h:88
int index_lt_ll
Definition: lensing.h:53
int has_dd
Definition: lensing.h:40
int index_lt_tp
Definition: lensing.h:50
int has_tl
Definition: lensing.h:43
int * l_max_lt
Definition: lensing.h:72
int index_lt_ee
Definition: lensing.h:46
int index_lt_tl
Definition: lensing.h:54
double * ddcl_lens
Definition: lensing.h:80
int lt_size
Definition: lensing.h:56
int index_lt_bb
Definition: lensing.h:48
int has_tp
Definition: lensing.h:39
int has_pp
Definition: lensing.h:38
double * cl_lens
Definition: lensing.h:76
int index_lt_pp
Definition: lensing.h:49
int l_unlensed_max
Definition: lensing.h:64
short has_lensed_cls
Definition: lensing.h:26
int l_lensed_max
Definition: lensing.h:66
int has_ll
Definition: lensing.h:42
Definition: lensing.h:17
Definition: perturbations.h:98