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;
92  short is_allocated;
95 };
96 
97 /*************************************************************************************************************/
98 /* @cond INCLUDE_WITH_DOXYGEN */
99 /*
100  * Boilerplate for C++
101  */
102 #ifdef __cplusplus
103 extern "C" {
104 #endif
105 
106  int lensing_cl_at_l(
107  struct lensing * ple,
108  int l,
109  double * cl_lensed
110  );
111 
112  int lensing_init(
113  struct precision * ppr,
114  struct perturbations * ppt,
115  struct harmonic * phr,
116  struct fourier * pfo,
117  struct lensing * ple
118  );
119 
120  int lensing_free(
121  struct lensing * ple
122  );
123 
124  int lensing_indices(
125  struct precision * ppr,
126  struct harmonic * phr,
127  struct lensing * ple
128  );
129 
131  double *ksi,
132  double **d00,
133  double *w8,
134  int nmu,
135  struct lensing * ple
136  );
137 
139  double *ksiX,
140  double **d20,
141  double *w8,
142  int nmu,
143  struct lensing * ple
144  );
145 
147  double *ksip,
148  double *ksim,
149  double **d22,
150  double **d2m2,
151  double *w8,
152  int nmu,
153  struct lensing * ple
154  );
156  struct lensing *ple,
157  double *cl_tt
158  );
159 
161  struct lensing *ple,
162  double *cl_te
163  );
164 
166  struct lensing *ple,
167  double *cl_ee,
168  double *cl_bb
169  );
170 
171 
172  int lensing_X000(
173  double * mu,
174  int num_mu,
175  int lmax,
176  double * sigma2,
177  double ** X000
178  );
179 
180  int lensing_Xp000(
181  double * mu,
182  int num_mu,
183  int lmax,
184  double * sigma2,
185  double ** Xp000
186  );
187 
188  int lensing_X220(
189  double * mu,
190  int num_mu,
191  int lmax,
192  double * sigma2,
193  double ** X220
194  );
195 
196  int lensing_X022(
197  double * mu,
198  int num_mu,
199  int lmax,
200  double * sigma2,
201  double ** X022
202  );
203 
204  int lensing_Xp022(
205  double * mu,
206  int num_mu,
207  int lmax,
208  double * sigma2,
209  double ** Xp022
210  );
211 
212  int lensing_X121(
213  double * mu,
214  int num_mu,
215  int lmax,
216  double * sigma2,
217  double ** X121
218  );
219 
220  int lensing_X132(
221  double * mu,
222  int num_mu,
223  int lmax,
224  double * sigma2,
225  double ** X132
226  );
227 
228  int lensing_X242(
229  double * mu,
230  int num_mu,
231  int lmax,
232  double * sigma2,
233  double ** X242
234  );
235 
236  int lensing_d00(
237  double * mu,
238  int num_mu,
239  int lmax,
240  double ** d00
241  );
242 
243  int lensing_d11(
244  double * mu,
245  int num_mu,
246  int lmax,
247  double ** d11
248  );
249 
250  int lensing_d1m1(
251  double * mu,
252  int num_mu,
253  int lmax,
254  double ** d1m1
255  );
256 
257  int lensing_d2m2(
258  double * mu,
259  int num_mu,
260  int lmax,
261  double ** d2m2
262  );
263 
264  int lensing_d22(
265  double * mu,
266  int num_mu,
267  int lmax,
268  double ** d22
269  );
270 
271  int lensing_d20(
272  double * mu,
273  int num_mu,
274  int lmax,
275  double ** d20
276  );
277 
278  int lensing_d31(
279  double * mu,
280  int num_mu,
281  int lmax,
282  double ** d3m1
283  );
284 
285  int lensing_d3m1(
286  double * mu,
287  int num_mu,
288  int lmax,
289  double ** d3m1
290  );
291 
292  int lensing_d3m3(
293  double * mu,
294  int num_mu,
295  int lmax,
296  double ** d3m3
297  );
298 
299  int lensing_d40(
300  double * mu,
301  int num_mu,
302  int lmax,
303  double ** d40
304  );
305 
306  int lensing_d4m2(
307  double * mu,
308  int num_mu,
309  int lmax,
310  double ** d4m2
311  );
312 
313  int lensing_d4m4(
314  double * mu,
315  int num_mu,
316  int lmax,
317  double ** d4m4
318  );
319 
320 #ifdef __cplusplus
321 }
322 #endif
323 
324 #endif
325 /* @endcond */
Definition: common.h:406
Definition: fourier.h:33
Definition: harmonic.h:17
int lensing_indices(struct precision *ppr, struct harmonic *phr, struct lensing *ple)
Definition: lensing.c:836
int lensing_d1m1(double *mu, int num_mu, int lmax, double **d1m1)
Definition: lensing.c:1375
int lensing_d40(double *mu, int num_mu, int lmax, double **d40)
Definition: lensing.c:1795
int lensing_d20(double *mu, int num_mu, int lmax, double **d20)
Definition: lensing.c:1554
int lensing_lensed_cl_te(double *ksiX, double **d20, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1115
int lensing_addback_cl_te(struct lensing *ple, double *cl_te)
Definition: lensing.c:1155
int lensing_lensed_cl_tt(double *ksi, double **d00, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1049
int lensing_d2m2(double *mu, int num_mu, int lmax, double **d2m2)
Definition: lensing.c:1434
int lensing_d4m2(double *mu, int num_mu, int lmax, double **d4m2)
Definition: lensing.c:1855
int lensing_addback_cl_ee_bb(struct lensing *ple, double *cl_ee, double *cl_bb)
Definition: lensing.c:1227
int lensing_d11(double *mu, int num_mu, int lmax, double **d11)
Definition: lensing.c:1315
int lensing_d00(double *mu, int num_mu, int lmax, double **d00)
Definition: lensing.c:1256
int lensing_addback_cl_tt(struct lensing *ple, double *cl_tt)
Definition: lensing.c:1090
int lensing_d31(double *mu, int num_mu, int lmax, double **d31)
Definition: lensing.c:1612
int lensing_free(struct lensing *ple)
Definition: lensing.c:808
int lensing_d3m1(double *mu, int num_mu, int lmax, double **d3m1)
Definition: lensing.c:1673
int lensing_d3m3(double *mu, int num_mu, int lmax, double **d3m3)
Definition: lensing.c:1734
int lensing_d22(double *mu, int num_mu, int lmax, double **d22)
Definition: lensing.c:1494
int lensing_cl_at_l(struct lensing *ple, int l, double *cl_lensed)
Definition: lensing.c:40
int lensing_lensed_cl_ee_bb(double *ksip, double *ksim, double **d22, double **d2m2, double *w8, int nmu, struct lensing *ple)
Definition: lensing.c:1182
int lensing_init(struct precision *ppr, struct perturbations *ppt, struct harmonic *phr, struct fourier *pfo, struct lensing *ple)
Definition: lensing.c:85
int lensing_d4m4(double *mu, int num_mu, int lmax, double **d4m4)
Definition: lensing.c:1917
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
short is_allocated
Definition: lensing.h:92
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