source: Sophya/trunk/SophyaLib/Samba/sphericaltransformserver.cc@ 2946

Last change on this file since 2946 was 2872, checked in by ansari, 20 years ago

Portage/compilation sur AIX-XlC (regatta): Ajout qualification de nom pour instanciation explicite de template - Reza 3 Jan 2006

File size: 46.3 KB
Line 
1#include "sopnamsp.h"
2#include "machdefs.h"
3#include <iostream>
4#include <math.h>
5#include <complex>
6#include "sphericaltransformserver.h"
7#include "tvector.h"
8#include "nbrandom.h"
9#include "nbmath.h"
10#include "timing.h"
11//#include "spherehealpix.h"
12
13
14/*!
15 \ingroup Samba
16 \class SOPHYA::SphericalTransformServer
17
18 \brief Analysis/synthesis in spherical harmonics server.
19
20 Class for performing analysis and synthesis of sky maps using spin-0 or spin-2 spherical harmonics.
21
22Maps must be SOPHYA SphericalMaps (SphereHEALPix or SphereThetaPhi or SphereECP).
23
24Temperature and polarization (Stokes parameters) can be developped on spherical harmonics :
25\f[
26\frac{\Delta T}{T}(\hat{n})=\sum_{lm}a_{lm}^TY_l^m(\hat{n})
27\f]
28\f[
29Q(\hat{n})=\frac{1}{\sqrt{2}}\sum_{lm}N_l\left(a_{lm}^EW_{lm}(\hat{n})+a_{lm}^BX_{lm}(\hat{n})\right)
30\f]
31\f[
32U(\hat{n})=-\frac{1}{\sqrt{2}}\sum_{lm}N_l\left(a_{lm}^EX_{lm}(\hat{n})-a_{lm}^BW_{lm}(\hat{n})\right)
33\f]
34\f[
35\left(Q \pm iU\right)(\hat{n})=\sum_{lm}a_{\pm 2lm}\, _{\pm 2}Y_l^m(\hat{n})
36\f]
37
38\f[
39Y_l^m(\hat{n})=\lambda_l^m(\theta)e^{im\phi}
40\f]
41\f[
42_{\pm}Y_l^m(\hat{n})=_{\pm}\lambda_l^m(\theta)e^{im\phi}
43\f]
44\f[
45W_{lm}(\hat{n})=\frac{1}{N_l}\,_{w}\lambda_l^m(\theta)e^{im\phi}
46\f]
47\f[
48X_{lm}(\hat{n})=\frac{-i}{N_l}\,_{x}\lambda_l^m(\theta)e^{im\phi}
49\f]
50
51(see LambdaLMBuilder, LambdaPMBuilder, LambdaWXBuilder classes)
52
53power spectra :
54
55\f[
56C_l^T=\frac{1}{2l+1}\sum_{m=0}^{+ \infty }\left|a_{lm}^T\right|^2=\langle\left|a_{lm}^T\right|^2\rangle
57\f]
58\f[
59C_l^E=\frac{1}{2l+1}\sum_{m=0}^{+\infty}\left|a_{lm}^E\right|^2=\langle\left|a_{lm}^E\right|^2\rangle
60\f]
61\f[
62C_l^B=\frac{1}{2l+1}\sum_{m=0}^{+\infty}\left|a_{lm}^B\right|^2=\langle\left|a_{lm}^B\right|^2\rangle
63\f]
64
65\arg
66\b Synthesis : Get temperature and polarization maps from \f$a_{lm}\f$ coefficients or from power spectra, (methods GenerateFrom...).
67
68\b Temperature:
69\f[
70\frac{\Delta T}{T}(\hat{n})=\sum_{lm}a_{lm}^TY_l^m(\hat{n}) = \sum_{-\infty}^{+\infty}b_m(\theta)e^{im\phi}
71\f]
72
73with
74\f[
75b_m(\theta)=\sum_{l=\left|m\right|}^{+\infty}a_{lm}^T\lambda_l^m(\theta)
76\f]
77
78\b Polarisation
79\f[
80Q \pm iU = \sum_{-\infty}^{+\infty}b_m^{\pm}(\theta)e^{im\phi}
81\f]
82
83where :
84\f[
85b_m^{\pm}(\theta) = \sum_{l=\left|m\right|}^{+\infty}a_{\pm 2lm}\,_{\pm}\lambda_l^m(\theta)
86\f]
87
88or :
89\f[
90Q = \sum_{-\infty}^{+\infty}b_m^{Q}(\theta)e^{im\phi}
91\f]
92\f[
93U = \sum_{-\infty}^{+\infty}b_m^{U}(\theta)e^{im\phi}
94\f]
95
96where:
97\f[
98b_m^{Q}(\theta) = \frac{1}{\sqrt{2}}\sum_{l=\left|m\right|}^{+\infty}\left(a_{lm}^E\,_{w}\lambda_l^m(\theta)-ia_{lm}^B\,_{x}\lambda_l^m(\theta)\right)
99\f]
100\f[
101b_m^{U}(\theta) = \frac{1}{\sqrt{2}}\sum_{l=\left|m\right|}^{+\infty}\left(ia_{lm}^E\,_{x}\lambda_l^m(\theta)+a_{lm}^B\,_{w}\lambda_l^m(\theta)\right)
102\f]
103
104Since the pixelization provides "slices" with constant \f$\theta\f$ and \f$\phi\f$ equally distributed on \f$2\pi\f$ \f$\frac{\Delta T}{T}\f$, \f$Q\f$,\f$U\f$ can be computed by FFT.
105
106
107\arg
108\b Analysis : Get \f$a_{lm}\f$ coefficients or power spectra from temperature and polarization maps (methods DecomposeTo...).
109
110\b Temperature:
111\f[
112a_{lm}^T=\int\frac{\Delta T}{T}(\hat{n})Y_l^{m*}(\hat{n})d\hat{n}
113\f]
114
115approximated as :
116\f[
117a_{lm}^T=\sum_{\theta_k}\omega_kC_m(\theta_k)\lambda_l^m(\theta_k)
118\f]
119where :
120\f[
121C_m (\theta _k)=\sum_{\phi _{k\prime}}\frac{\Delta T}{T}(\theta _k,\phi_{k\prime})e^{-im\phi _{k\prime}}
122\f]
123Since the pixelization provides "slices" with constant \f$\theta\f$ and \f$\phi\f$ equally distributed on \f$2\pi\f$ (\f$\omega_k\f$ is the solid angle of each pixel of the slice \f$\theta_k\f$) \f$C_m\f$ can be computed by FFT.
124
125\b polarisation:
126
127\f[
128a_{\pm 2lm}=\sum_{\theta_k}\omega_kC_m^{\pm}(\theta_k)\,_{\pm}\lambda_l^m(\theta_k)
129\f]
130where :
131\f[
132C_m^{\pm} (\theta _k)=\sum_{\phi _{k\prime}}\left(Q \pm iU\right)(\theta _k,\phi_{k\prime})e^{-im\phi _{k\prime}}
133\f]
134or :
135
136\f[
137a_{lm}^E=\frac{1}{\sqrt{2}}\sum_{\theta_k}\omega_k\left(C_m^{Q}(\theta_k)\,_{w}\lambda_l^m(\theta_k)-iC_m^{U}(\theta_k)\,_{x}\lambda_l^m(\theta_k)\right)
138\f]
139\f[
140a_{lm}^B=\frac{1}{\sqrt{2}}\sum_{\theta_k}\omega_k\left(iC_m^{Q}(\theta_k)\,_{x}\lambda_l^m(\theta_k)+C_m^{U}(\theta_k)\,_{w}\lambda_l^m(\theta_k)\right)
141\f]
142
143where :
144\f[
145C_m^{Q} (\theta _k)=\sum_{\phi _{k\prime}}Q(\theta _k,\phi_{k\prime})e^{-im\phi _{k\prime}}
146\f]
147\f[
148C_m^{U} (\theta _k)=\sum_{\phi _{k\prime}}U(\theta _k,\phi_{k\prime})e^{-im\phi _{k\prime}}
149\f]
150
151 */
152
153 /*! \fn void SOPHYA::SphericalTransformServer::GenerateFromAlm( SphericalMap<T>& map, int_4 pixelSizeIndex, const Alm<T>& alm) const
154
155 synthesis of a temperature map from Alm coefficients
156*/
157template<class T>
158void SphericalTransformServer<T>::GenerateFromAlm( SphericalMap<T>& map, int_4 pixelSizeIndex, const Alm<T>& alm) const
159{
160 /*=======================================================================
161 computes a map from its alm for the HEALPIX pixelisation
162 map(theta,phi) = sum_l_m a_lm Y_lm(theta,phi)
163 = sum_m {e^(i*m*phi) sum_l a_lm*lambda_lm(theta)}
164
165 where Y_lm(theta,phi) = lambda(theta) * e^(i*m*phi)
166
167 * the recurrence of Ylm is the standard one (cf Num Rec)
168 * the sum over m is done by FFT
169
170 =======================================================================*/
171 int_4 nlmax=alm.Lmax();
172 int_4 nmmax=nlmax;
173 int_4 nsmax=0;
174 // le Resize est suppose mettre a zero
175 map.Resize(pixelSizeIndex);
176 string sphere_type=map.TypeOfMap();
177 int premiereTranche = 0;
178 int derniereTranche = map.NbThetaSlices()-1;
179 if (sphere_type.substr(0,4) == "RING")
180 {
181 nsmax=map.SizeIndex();
182 }
183 else
184 {
185 // pour une sphere Gorski le nombre de pixels est 12*nsmax**2
186 // on calcule une quantite equivalente a nsmax pour la sphere-theta-phi
187 // en vue de l'application du critere Healpix : nlmax<=3*nsmax-1
188 // c'est approximatif ; a raffiner.
189 if (sphere_type.substr(0,6) == "TETAFI")
190 {
191 nsmax=(int_4)sqrt(map.NbPixels()/12.);
192 premiereTranche++;
193 derniereTranche--;
194 }
195 else
196 {
197 cout << " unknown type of sphere : " << sphere_type << endl;
198 throw IOExc(" unknown type of sphere: " + (string)sphere_type );
199 }
200 // cout << "GenerateFromAlm: the sphere is of type : " << sphere_type << endl;
201 // cout << "GenerateFromAlm: size index (nside) of the sphere= " << nsmax << endl;
202 // cout << "GenerateFromAlm: nlmax (from Alm) = " << nlmax << endl;
203 // if (nlmax>3*nsmax-1)
204 // {
205 // cout << "GenerateFromAlm: nlmax should be <= 3*nside-1" << endl;
206 // if (strncmp(sphere_type,"TETAFI",6) == 0)
207 // {
208 // cout << "GenerateFromAlm: nlmax should be <= 3*nside-1" << endl;
209 // cout << " (for this criterium, nsmax is computed as sqrt(nbPixels/12))" << endl;
210 // }
211 //}
212 }
213 Bm<complex<T> > b_m_theta(nmmax);
214
215 // map.Resize(nsmax);
216
217
218 // pour chaque tranche en theta
219 for (int_4 ith = premiereTranche; ith <= derniereTranche;ith++)
220 {
221 int_4 nph;
222 r_8 phi0;
223 r_8 theta;
224 TVector<int_4> pixNumber;
225 TVector<T> datan;
226
227 map.GetThetaSlice(ith,theta,phi0, pixNumber,datan);
228 nph = pixNumber.NElts();
229
230 // -----------------------------------------------------
231 // for each theta, and each m, computes
232 // b(m,theta) = sum_over_l>m (lambda_l_m(theta) * a_l_m)
233 // ------------------------------------------------------
234 LambdaLMBuilder lb(theta,nlmax,nmmax);
235 // somme sur m de 0 a l'infini
236 int m;
237 for (m = 0; m <= nmmax; m++)
238 {
239 b_m_theta(m) = (T)( lb.lamlm(m,m) ) * alm(m,m);
240 for (int l = m+1; l<= nlmax; l++)
241 {
242 b_m_theta(m) += (T)( lb.lamlm(l,m) ) * alm(l,m);
243 }
244 }
245 // obtains the negative m of b(m,theta) (= complex conjugate)
246
247 for (m=1;m<=nmmax;m++)
248 {
249 b_m_theta(-m) = conj(b_m_theta(m));
250 }
251 // ---------------------------------------------------------------
252 // sum_m b(m,theta)*exp(i*m*phi) -> f(phi,theta)
253 // ---------------------------------------------------------------*/
254
255
256 if (sphere_type.substr(0,4) == "RING")
257 {
258 TVector<T> Temp = RfourierSynthesisFromB(b_m_theta,nph,phi0);
259 for (int i=0;i< nph;i++) map(pixNumber(i))=Temp(i);
260 }
261 else
262 // pour des pixelisations quelconques (autres que HEALPix
263 // nph n'est pas toujours pair
264 // ca fait des problemes pour les transformees de Fourier
265 // car le server de TF ajuste la longueur du vecteur reel
266 // en sortie de TF, bref, la securite veut qu'on prenne une
267 // TF complexe
268 {
269 TVector<complex<T> > Temp = fourierSynthesisFromB(b_m_theta,nph,phi0);
270 for (int i=0;i< nph;i++) map(pixNumber(i))=Temp(i).real();
271 }
272 }
273}
274
275
276
277 /*! \fn TVector< complex<T> > SOPHYA::SphericalTransformServer::fourierSynthesisFromB(const Bm<complex<T> >& b_m, int_4 nph, r_8 phi0) const
278
279\return a vector with nph elements which are sums :\f$\sum_{m=-mmax}^{mmax}b_m(\theta)e^{im\varphi}\f$ for nph values of \f$\varphi\f$ regularly distributed in \f$[0,\pi]\f$ ( calculated by FFT)
280
281 The object b_m (\f$b_m\f$) of the class Bm is a special vector which index goes from -mmax to mmax.
282 */
283template<class T>
284TVector< complex<T> > SphericalTransformServer<T>::fourierSynthesisFromB(const Bm<complex<T> >& b_m, int_4 nph, r_8 phi0) const
285{
286 /*=======================================================================
287 dataout(j) = sum_m datain(m) * exp(i*m*phi(j))
288 with phi(j) = j*2pi/nph + kphi0*pi/nph and kphi0 =0 or 1
289
290 as the set of frequencies {m} is larger than nph,
291 we wrap frequencies within {0..nph-1}
292 ie m = k*nph + m' with m' in {0..nph-1}
293 then
294 noting bw(m') = exp(i*m'*phi0)
295 * sum_k (datain(k*nph+m') exp(i*k*pi*kphi0))
296 with bw(nph-m') = CONJ(bw(m')) (if datain(-m) = CONJ(datain(m)))
297 dataout(j) = sum_m' [ bw(m') exp (i*j*m'*2pi/nph) ]
298 = Fourier Transform of bw
299 is real
300
301 NB nph is not necessarily a power of 2
302
303=======================================================================*/
304 //**********************************************************************
305 // pour une valeur de phi (indexee par j) la temperature est la transformee
306 // de Fourier de bm (somme sur m de -nmax a +nmmax de bm*exp(i*m*phi)).
307 // on demande nph (nombre de pixels sur la tranche) valeurs de transformees, pour nph valeurs de phi, regulierement reparties sur 2*pi. On a:
308 // DT/T(j) = sum_m b(m) * exp(i*m*phi(j))
309 // sommation de -infini a +infini, en fait limitee a -nmamx, +nmmax
310 // On pose m=k*nph + m', avec m' compris entre 0 et nph-1. Alors :
311 // DT/T(j) = somme_k somme_m' b(k*nph + m')*exp(i*(k*nph + m')*phi(j))
312 // somme_k : de -infini a +infini
313 // somme_m' : de 0 a nph-1
314 // On echange les sommations :
315 // DT/T(j) = somme_m' (exp(i*m'*phi(j)) somme_k b(k*nph + m')*exp(i*(k*nph*phi(j))
316 // mais phi(j) est un multiple entier de 2*pi/nph, la seconde exponentielle
317 // vaut 1.
318 // Il reste a calculer les transformees de Fourier de somme_m' b(k*nph + m')
319 // si phi0 n'est pas nul, il y a juste un decalage a faire.
320 //**********************************************************************
321
322 TVector< complex<T> > bw(nph);
323 TVector< complex<T> > dataout(nph);
324 TVector< complex<T> > data(nph);
325
326
327 for (int kk=0; kk<bw.NElts(); kk++) bw(kk)=(T)0.;
328 int m;
329 for (m=-b_m.Mmax();m<=-1;m++)
330 {
331 int maux=m;
332 while (maux<0) maux+=nph;
333 int iw=maux%nph;
334 double aux=(m-iw)*phi0;
335 bw(iw) += b_m(m) * complex<T>( (T)cos(aux),(T)sin(aux) ) ;
336 }
337 for (m=0;m<=b_m.Mmax();m++)
338 {
339 // int iw=((m % nph) +nph) % nph; //between 0 and nph = m'
340 int iw=m%nph;
341 double aux=(m-iw)*phi0;
342 bw(iw)+=b_m(m) * complex<T>( (T)cos(aux),(T)sin(aux) );
343 }
344
345 // applies the shift in position <-> phase factor in Fourier space
346 for (int mprime=0; mprime < nph; mprime++)
347 {
348 complex<double> aux(cos(mprime*phi0),sin(mprime*phi0));
349 data(mprime)=bw(mprime)*
350 (complex<T>)(complex<double>(cos(mprime*phi0),sin(mprime*phi0)));
351 }
352
353 //sortie.ReSize(nph);
354 TVector< complex<T> > sortie(nph);
355
356 fftIntfPtr_-> FFTBackward(data, sortie);
357
358 return sortie;
359}
360
361//********************************************
362/*! \fn TVector<T> SOPHYA::SphericalTransformServer::RfourierSynthesisFromB(const Bm<complex<T> >& b_m, int_4 nph, r_8 phi0) const
363
364same as fourierSynthesisFromB, but return a real vector, taking into account the fact that b(-m) is conjugate of b(m) */
365template<class T>
366TVector<T> SphericalTransformServer<T>::RfourierSynthesisFromB(const Bm<complex<T> >& b_m, int_4 nph, r_8 phi0) const
367{
368 /*=======================================================================
369 dataout(j) = sum_m datain(m) * exp(i*m*phi(j))
370 with phi(j) = j*2pi/nph + kphi0*pi/nph and kphi0 =0 or 1
371
372 as the set of frequencies {m} is larger than nph,
373 we wrap frequencies within {0..nph-1}
374 ie m = k*nph + m' with m' in {0..nph-1}
375 then
376 noting bw(m') = exp(i*m'*phi0)
377 * sum_k (datain(k*nph+m') exp(i*k*pi*kphi0))
378 with bw(nph-m') = CONJ(bw(m')) (if datain(-m) = CONJ(datain(m)))
379 dataout(j) = sum_m' [ bw(m') exp (i*j*m'*2pi/nph) ]
380 = Fourier Transform of bw
381 is real
382
383 NB nph is not necessarily a power of 2
384
385=======================================================================*/
386 //**********************************************************************
387 // pour une valeur de phi (indexee par j) la temperature est la transformee
388 // de Fourier de bm (somme sur m de -nmax a +nmmax de bm*exp(i*m*phi)).
389 // on demande nph (nombre de pixels sur la tranche) valeurs de transformees, pour nph valeurs de phi, regulierement reparties sur 2*pi. On a:
390 // DT/T(j) = sum_m b(m) * exp(i*m*phi(j))
391 // sommation de -infini a +infini, en fait limitee a -nmamx, +nmmax
392 // On pose m=k*nph + m', avec m' compris entre 0 et nph-1. Alors :
393 // DT/T(j) = somme_k somme_m' b(k*nph + m')*exp(i*(k*nph + m')*phi(j))
394 // somme_k : de -infini a +infini
395 // somme_m' : de 0 a nph-1
396 // On echange les sommations :
397 // DT/T(j) = somme_m' (exp(i*m'*phi(j)) somme_k b(k*nph + m')*exp(i*(k*nph*phi(j))
398 // mais phi(j) est un multiple entier de 2*pi/nph, la seconde exponentielle
399 // vaut 1.
400 // Il reste a calculer les transformees de Fourier de somme_k b(k*nph + m')
401 // si phi0 n'est pas nul, il y a juste un decalage a faire.
402 //**********************************************************************
403 TVector< complex<T> > bw(nph);
404 TVector< complex<T> > dataout(nph);
405 TVector< complex<T> > data(nph/2+1);
406
407
408 for (int kk=0; kk<bw.NElts(); kk++) bw(kk)=(T)0.;
409 int m;
410 for (m=-b_m.Mmax();m<=-1;m++)
411 {
412 int maux=m;
413 while (maux<0) maux+=nph;
414 int iw=maux%nph;
415 double aux=(m-iw)*phi0;
416 bw(iw) += b_m(m) * complex<T>( (T)cos(aux),(T)sin(aux) ) ;
417 }
418 for (m=0;m<=b_m.Mmax();m++)
419 {
420 // int iw=((m % nph) +nph) % nph; //between 0 and nph = m'
421 int iw=m%nph;
422 double aux=(m-iw)*phi0;
423 bw(iw)+=b_m(m) * complex<T>( (T)cos(aux),(T)sin(aux) );
424 }
425
426 // applies the shift in position <-> phase factor in Fourier space
427 // cout << " TF : nph= " << nph << " vec. entree " << data.Size() << endl;
428 for (int mprime=0; mprime <= nph/2; mprime++)
429 {
430 complex<double> aux(cos(mprime*phi0),sin(mprime*phi0));
431 data(mprime)=bw(mprime)*
432 (complex<T>)(complex<double>(cos(mprime*phi0),sin(mprime*phi0)));
433 }
434
435 TVector<T> sortie;
436 fftIntfPtr_-> FFTBackward(data, sortie);
437
438 return sortie;
439}
440//*******************************************
441
442 /*! \fn Alm<T> SOPHYA::SphericalTransformServer::DecomposeToAlm(const SphericalMap<T>& map, int_4 nlmax, r_8 cos_theta_cut) const
443
444\return the Alm coefficients from analysis of a temperature map.
445
446 \param<nlmax> : maximum value of the l index
447
448 \param<cos_theta_cut> : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
449
450 */
451template<class T>
452void SphericalTransformServer<T>::DecomposeToAlm(const SphericalMap<T>& map, Alm<T>& alm, int_4 nlmax, r_8 cos_theta_cut) const
453{
454 DecomposeToAlm(const_cast< SphericalMap<T>& >(map), alm, nlmax, cos_theta_cut, 0);
455}
456//*******************************************
457
458 /*! \fn Alm<T> SOPHYA::SphericalTransformServer::DecomposeToAlm(const SphericalMap<T>& map, int_4 nlmax, r_8 cos_theta_cut, int iterationOrder) const
459
460\return the Alm coefficients from analysis of a temperature map. THE MAP CAN BE MODIFIED (if iterationOrder >0)
461
462 \param<nlmax> : maximum value of the l index
463
464 \param<cos_theta_cut> : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
465
466\param<iterationOrder> : 1,2,3,4.... order of an iterative analysis. (Default : 0 -> standard analysis). If iterationOrder is not null, the method works with SphereHEALPix but NOT WITH SphereThetaPhi maps ! */
467template<class T>
468void SphericalTransformServer<T>::DecomposeToAlm(SphericalMap<T>& map, Alm<T>& alm, int_4 nlmax, r_8 cos_theta_cut, int iterationOrder) const
469{
470 int_4 nmmax = nlmax;
471 // PrtTim("appel carteVersAlm");
472 carteVersAlm(map, nlmax, cos_theta_cut, alm);
473 // PrtTim("retour carteVersAlm");
474 if (iterationOrder > 0)
475 {
476 TVector<int_4> fact(iterationOrder+2);
477 fact(0) = 1;
478 int k;
479 for (k=1; k <= iterationOrder+1; k++)
480 {
481 fact(k) = fact(k-1)*k;
482 }
483 Alm<T> alm2(alm);
484 T Tzero = (T)0.;
485 complex<T> complexZero = complex<T>(Tzero, Tzero);
486 alm = complexZero;
487 int signe = 1;
488 int nbIteration = iterationOrder+1;
489 for (k=1; k <= nbIteration; k++)
490 {
491 T facMult = (T)(0.5*signe*fact(iterationOrder)*(2*nbIteration-k)/(fact(k)*fact(nbIteration-k)));
492 for (int m = 0; m <= nmmax; m++)
493 {
494 for (int l = m; l<= nlmax; l++)
495 {
496 alm(l,m) += facMult*alm2(l,m);
497 }
498 }
499 if (k == nbIteration) break;
500 signe = -signe;
501 for (int k=0; k< map.NbPixels(); k++) map(k) = (T)0.;
502 // synthetize a map from the estimated alm
503 // PrtTim("appel GenerateFromAlm");
504 GenerateFromAlm( map, map.SizeIndex(), alm2);
505 // PrtTim("retour GenerateFromAlm");
506 alm2 = complexZero;
507 // analyse the new map
508 // PrtTim("appel carteVersAlm");
509 carteVersAlm(map, nlmax, cos_theta_cut, alm2);
510 // PrtTim("retour carteVersAlm");
511 }
512 }
513}
514
515template<class T>
516 void SphericalTransformServer<T>::carteVersAlm(const SphericalMap<T>& map, int_4 nlmax, r_8 cos_theta_cut, Alm<T>& alm) const
517{
518
519 /*-----------------------------------------------------------------------
520 computes the integral in phi : phas_m(theta)
521 for each parallele from north to south pole
522 -----------------------------------------------------------------------*/
523 TVector<T> data;
524 TVector<int_4> pixNumber;
525 int_4 nmmax = nlmax;
526 TVector< complex<T> > phase(nmmax+1);
527
528 alm.ReSizeToLmax(nlmax);
529 for (int_4 ith = 0; ith < map.NbThetaSlices(); ith++)
530 {
531 r_8 phi0;
532 r_8 theta;
533 // PrtTim("debut 1ere tranche ");
534 map.GetThetaSlice(ith,theta,phi0,pixNumber ,data);
535 phase = complex<T>((T)0.,(T)0.);
536 double cth = cos(theta);
537
538 //part of the sky out of the symetric cut
539 bool keep_it = (fabs(cth) >= cos_theta_cut);
540
541 // PrtTim("fin 1ere tranche ");
542
543 if (keep_it)
544 {
545 // phase = CFromFourierAnalysis(nmmax,data,phi0);
546 // PrtTim("avant Fourier ");
547 CFromFourierAnalysis(nmmax,data,phase, phi0);
548 // PrtTim("apres Fourier ");
549
550 }
551
552// ---------------------------------------------------------------------
553// computes the a_lm by integrating over theta
554// lambda_lm(theta) * phas_m(theta)
555// for each m and l
556// -----------------------------------------------------------------------
557 // PrtTim("avant instanciation LM ");
558 LambdaLMBuilder lb(theta,nlmax,nmmax);
559 // PrtTim("apres instanciation LM ");
560 r_8 domega=map.PixSolAngle(map.PixIndexSph(theta,phi0));
561
562 // PrtTim("avant mise a jour Alm ");
563 complex<T> fi;
564 T facteur;
565 int index;
566 for (int m = 0; m <= nmmax; m++)
567 {
568 fi = phase(m);
569 for (int l = m; l<= nlmax; l++)
570 {
571 index = alm.indexOfElement(l,m);
572 // facteur = (T)(lb.lamlm(l,m) * domega);
573 facteur = (T)(lb.lamlm(index) * domega);
574 // alm(l,m) += facteur * fi ;
575 alm(index) += facteur * fi ;
576 }
577 }
578
579
580
581 //
582 //
583 // PrtTim("apres mise a jour Alm ");
584 }
585}
586 /*! \fn TVector< complex<T> > SOPHYA::SphericalTransformServer::CFromFourierAnalysis(int_4 nmmax, const TVector<complex<T> >datain, r_8 phi0) const
587
588\return a vector with mmax elements which are sums :
589\f$\sum_{k=0}^{nphi}datain(\theta,\varphi_k)e^{im\varphi_k}\f$ for (mmax+1) values of \f$m\f$ from 0 to mmax.
590 */
591template<class T>
592TVector< complex<T> > SphericalTransformServer<T>::CFromFourierAnalysis(int_4 nmmax, const TVector<complex<T> >datain, r_8 phi0) const
593{
594 /*=======================================================================
595 integrates (data * phi-dependence-of-Ylm) over phi
596 --> function of m can be computed by FFT
597
598 datain est modifie
599 =======================================================================*/
600 int_4 nph=datain.NElts();
601 if (nph <= 0)
602 {
603 throw PException("bizarre : vecteur datain de longueur nulle (CFromFourierAnalysis)");
604 }
605 TVector<complex<T> > transformedData(nph);
606 fftIntfPtr_-> FFTForward(datain, transformedData);
607
608 TVector< complex<T> > dataout(nmmax+1);
609
610 int im_max=min(nph,nmmax+1);
611 int i;
612 dataout = complex<T>((T)0.,(T)0.);
613 // for (i=0;i< dataout.NElts();i++) dataout(i)=complex<T>((T)0.,(T)0.);
614 for (i=0;i<im_max;i++) dataout(i)=transformedData(i);
615
616
617 for (int kk=nph; kk<dataout.NElts(); kk++) dataout(kk)=dataout(kk%nph);
618 for (i = 0;i <dataout.NElts();i++){
619 dataout(i)*= (complex<T>)(complex<double>(cos(-i*phi0),sin(-i*phi0)));
620 }
621 return dataout;
622}
623
624//&&&&&&&&& nouvelle version
625/* \fn TVector< complex<T> > SOPHYA::SphericalTransformServer::CFromFourierAnalysis(int_4 nmmax, const TVector<T> datain, r_8 phi0) const
626
627same as previous one, but with a "datain" which is real (not complex) */
628template<class T>
629void SphericalTransformServer<T>::CFromFourierAnalysis(int_4 nmmax, const TVector<T> datain, TVector< complex<T> >& dataout, r_8 phi0) const
630{
631 //=======================================================================
632 // integrates (data * phi-dependence-of-Ylm) over phi
633 // --> function of m can be computed by FFT
634 // ! with 0<= m <= npoints/2 (: Nyquist)
635 // ! because the data is real the negative m are the conjugate of the
636 // ! positive ones
637
638 // datain est modifie
639 //
640 // =======================================================================
641 int_4 nph=datain.NElts();
642 if (nph <= 0)
643 {
644 throw PException("bizarre : vecteur datain de longueur nulle (CFromFourierAnalysis)");
645 }
646 // if (nph%2 != 0 )
647 // {
648 // throw PException("SphericalTransformServer<T>::CFromFourierAnalysis : longueur de datain impair ?");
649 // }
650 TVector<complex<T> > transformedData;
651
652 // la taille du vecteur complexe retourne est nph/2+1 (si la taille
653 // du vecteur reel entre est nph)
654 // cout << " longueur de datain = " << nph << endl;
655 fftIntfPtr_-> FFTForward(datain, transformedData);
656 // cout << " taille de la transformee " << transformedData.Size() << endl;
657 // TVector< complex<T> > dataout(nmmax+1);
658 dataout.ReSize(nmmax+1);
659
660 // on transfere le resultat de la fft dans dataout.
661
662 int maxFreqAccessiblesParFFT = min(nph/2,nmmax);
663 int i;
664 for (i=0;i<=maxFreqAccessiblesParFFT;i++) dataout(i)=transformedData(i);
665
666
667 // si dataout n'est pas plein, on complete jusqu'a nph+1 valeurs (a moins
668 // que dataout ne soit plein avant d'atteindre nph)
669 if (maxFreqAccessiblesParFFT != nmmax )
670 {
671 int maxMfft = min(nph,nmmax);
672 for (i=maxFreqAccessiblesParFFT+1; i<=maxMfft; i++)
673 {
674 dataout(i) = conj(dataout(nph-i) );
675 }
676 // on conplete, si necessaire, par periodicite
677 if ( maxMfft != nmmax )
678 {
679 for (int kk=nph+1; kk <= nmmax; kk++)
680 {
681 dataout(kk)=dataout(kk%nph);
682 }
683 }
684 }
685 for (i = 0;i <dataout.NElts();i++)
686 {
687 dataout(i)*= (complex<T>)(complex<double>(cos(-i*phi0),sin(-i*phi0)));
688 }
689 // return dataout;
690}
691
692 /*! \fn void SOPHYA::SphericalTransformServer::GenerateFromAlm(SphericalMap<T>& mapq,
693 SphericalMap<T>& mapu,
694 int_4 pixelSizeIndex,
695 const Alm<T>& alme,
696 const Alm<T>& almb) const
697
698synthesis of a polarization map from Alm coefficients. The spheres mapq and mapu contain respectively the Stokes parameters. */
699template<class T>
700void SphericalTransformServer<T>::GenerateFromAlm(SphericalMap<T>& mapq,
701 SphericalMap<T>& mapu,
702 int_4 pixelSizeIndex,
703 const Alm<T>& alme,
704 const Alm<T>& almb) const
705{
706 /*=======================================================================
707 computes a map form its alm for the HEALPIX pixelisation
708 map(theta,phi) = sum_l_m a_lm Y_lm(theta,phi)
709 = sum_m {e^(i*m*phi) sum_l a_lm*lambda_lm(theta)}
710
711 where Y_lm(theta,phi) = lambda(theta) * e^(i*m*phi)
712
713 * the recurrence of Ylm is the standard one (cf Num Rec)
714 * the sum over m is done by FFT
715
716 =======================================================================*/
717 int_4 nlmax=alme.Lmax();
718 if (nlmax != almb.Lmax())
719 {
720 cout << " SphericalTransformServer: les deux tableaux alm n'ont pas la meme taille" << endl;
721 throw SzMismatchError("SphericalTransformServer: les deux tableaux alm n'ont pas la meme taille");
722 }
723 int_4 nmmax=nlmax;
724 int_4 nsmax=0;
725 mapq.Resize(pixelSizeIndex);
726 mapu.Resize(pixelSizeIndex);
727 string sphere_type=mapq.TypeOfMap();
728 if (sphere_type != mapu.TypeOfMap())
729 {
730 cout << " SphericalTransformServer: les deux spheres ne sont pas de meme type" << endl;
731 cout << " type 1 " << sphere_type << endl;
732 cout << " type 2 " << mapu.TypeOfMap() << endl;
733 throw SzMismatchError("SphericalTransformServer: les deux spheres ne sont pas de meme type");
734
735 }
736 bool healpix = true;
737 if (sphere_type.substr(0,4) == "RING")
738 {
739 nsmax=mapq.SizeIndex();
740 }
741 else
742 // pour une sphere Gorski le nombre de pixels est 12*nsmax**2
743 // on calcule une quantite equivalente a nsmax pour la sphere-theta-phi
744 // en vue de l'application du critere Healpix : nlmax<=3*nsmax-1
745 // c'est approximatif ; a raffiner.
746 healpix = false;
747 if (sphere_type.substr(0,6) == "TETAFI")
748 {
749 nsmax=(int_4)sqrt(mapq.NbPixels()/12.);
750 }
751 else
752 {
753 cout << " unknown type of sphere : " << sphere_type << endl;
754 throw IOExc(" unknown type of sphere ");
755 }
756 cout << "GenerateFromAlm: the spheres are of type : " << sphere_type << endl;
757 cout << "GenerateFromAlm: size indices (nside) of spheres= " << nsmax << endl;
758 cout << "GenerateFromAlm: nlmax (from Alm) = " << nlmax << endl;
759 if (nlmax>3*nsmax-1)
760 {
761 cout << "GenerateFromAlm: nlmax should be <= 3*nside-1" << endl;
762 if (sphere_type.substr(0,6) == "TETAFI")
763 {
764 cout << " (for this criterium, nsmax is computed as sqrt(nbPixels/12))" << endl;
765 }
766 }
767 if (alme.Lmax()!=almb.Lmax())
768 {
769 cout << "GenerateFromAlm: arrays Alme and Almb have not the same size ? " << endl;
770 throw SzMismatchError("SphericalTransformServer: arrays Alme and Almb have not the same size ? ");
771 }
772 mapFromWX(nlmax, nmmax, mapq, mapu, alme, almb, healpix);
773 // mapFromPM(nlmax, nmmax, mapq, mapu, alme, almb);
774}
775 /*! \fn void SOPHYA::SphericalTransformServer::DecomposeToAlm(const SphericalMap<T>& mapq,
776 const SphericalMap<T>& mapu,
777 Alm<T>& alme,
778 Alm<T>& almb,
779 int_4 nlmax,
780 r_8 cos_theta_cut) const
781
782analysis of a polarization map into Alm coefficients.
783
784 The spheres \c mapq and \c mapu contain respectively the Stokes parameters.
785
786 \c a2lme and \c a2lmb will receive respectively electric and magnetic Alm's
787 nlmax : maximum value of the l index
788
789 \c cos_theta_cut : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
790
791
792 */
793template<class T>
794void SphericalTransformServer<T>::DecomposeToAlm(const SphericalMap<T>& mapq,
795 const SphericalMap<T>& mapu,
796 Alm<T>& alme,
797 Alm<T>& almb,
798 int_4 nlmax,
799 r_8 cos_theta_cut) const
800{
801 DecomposeToAlm(const_cast< SphericalMap<T>& >(mapq), const_cast< SphericalMap<T>& >(mapu), alme, almb, nlmax, cos_theta_cut);
802}
803
804 /*! \fn void SOPHYA::SphericalTransformServer::DecomposeToAlm(const SphericalMap<T>& mapq,
805 const SphericalMap<T>& mapu,
806 Alm<T>& alme,
807 Alm<T>& almb,
808 int_4 nlmax,
809 r_8 cos_theta_cut,
810 int iterationOrder) const
811
812analysis of a polarization map into Alm coefficients.
813
814 The spheres \c mapq and \c mapu contain respectively the Stokes parameters.
815
816 \c a2lme and \c a2lmb will receive respectively electric and magnetic Alm's
817 nlmax : maximum value of the l index
818
819 \c cos_theta_cut : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
820
821\param<iterationOrder> : 1,2,3,4.... order of an iterative analysis. (Default : 0 -> standard analysis). If iterationOrder is not null, the method works with SphereHEALPix but NOT WITH SphereThetaPhi maps !
822
823THE INPUT MAPS CAN BE MODIFIED (only if iterationOrder >0)
824
825 */
826template<class T>
827void SphericalTransformServer<T>::DecomposeToAlm(SphericalMap<T>& mapq,
828 SphericalMap<T>& mapu,
829 Alm<T>& alme,
830 Alm<T>& almb,
831 int_4 nlmax,
832 r_8 cos_theta_cut,
833 int iterationOrder) const
834{
835 int_4 nmmax = nlmax;
836 carteVersAlm(mapq, mapu, alme, almb, nlmax, cos_theta_cut);
837 if (iterationOrder > 0)
838 {
839 TVector<int_4> fact(iterationOrder+2);
840 fact(0) = 1;
841 int k;
842 for (k=1; k <= iterationOrder+1; k++)
843 {
844 fact(k) = fact(k-1)*k;
845 }
846 Alm<T> alme2(alme);
847 Alm<T> almb2(almb);
848 T Tzero = (T)0.;
849 complex<T> complexZero = complex<T>(Tzero, Tzero);
850 alme = complexZero;
851 almb = complexZero;
852 int signe = 1;
853 int nbIteration = iterationOrder+1;
854 for (k=1; k <= nbIteration; k++)
855 {
856 T facMult = (T)(0.5*signe*fact(iterationOrder)*(2*nbIteration-k)/(fact(k)*fact(nbIteration-k)));
857 for (int m = 0; m <= nmmax; m++)
858 {
859 for (int l = m; l<= nlmax; l++)
860 {
861 alme(l,m) += facMult*alme2(l,m);
862 almb(l,m) += facMult*almb2(l,m);
863 }
864 }
865 if (k == nbIteration) break;
866 signe = -signe;
867 for (int k=0; k< mapq.NbPixels(); k++)
868 {
869 mapq(k) = (T)0.;
870 mapu(k) = (T)0.;
871 }
872 // synthetize a map from the estimated alm
873 GenerateFromAlm(mapq,mapu,mapq.SizeIndex(),alme2,almb2);
874 alme2 = complexZero;
875 almb2 = complexZero;
876 // analyse the new map
877 carteVersAlm(mapq, mapu, alme2, almb2, nlmax, cos_theta_cut);
878 }
879 }
880}
881
882template<class T>
883void SphericalTransformServer<T>::carteVersAlm(const SphericalMap<T>& mapq,
884 const SphericalMap<T>& mapu,
885 Alm<T>& alme,
886 Alm<T>& almb,
887 int_4 nlmax,
888 r_8 cos_theta_cut) const
889{
890 int_4 nmmax = nlmax;
891 // resize et remise a zero
892 alme.ReSizeToLmax(nlmax);
893 almb.ReSizeToLmax(nlmax);
894
895
896 TVector<T> dataq;
897 TVector<T> datau;
898 TVector<int_4> pixNumber;
899
900 string sphere_type=mapq.TypeOfMap();
901 if (sphere_type != mapu.TypeOfMap())
902 {
903 cout << " SphericalTransformServer: les deux spheres ne sont pas de meme type" << endl;
904 cout << " type 1 " << sphere_type << endl;
905 cout << " type 2 " << mapu.TypeOfMap() << endl;
906 throw SzMismatchError("SphericalTransformServer: les deux spheres ne sont pas de meme type");
907
908 }
909 if (mapq.NbPixels()!=mapu.NbPixels())
910 {
911 cout << " DecomposeToAlm: map Q and map U have not same size ?" << endl;
912 throw SzMismatchError("SphericalTransformServer::DecomposeToAlm: map Q and map U have not same size ");
913 }
914 for (int_4 ith = 0; ith < mapq.NbThetaSlices(); ith++)
915 {
916 r_8 phi0;
917 r_8 theta;
918 mapq.GetThetaSlice(ith,theta,phi0, pixNumber,dataq);
919 mapu.GetThetaSlice(ith,theta,phi0, pixNumber,datau);
920 if (dataq.NElts() != datau.NElts() )
921 {
922 throw SzMismatchError("the spheres have not the same pixelization");
923 }
924 r_8 domega=mapq.PixSolAngle(mapq.PixIndexSph(theta,phi0));
925 double cth = cos(theta);
926 //part of the sky out of the symetric cut
927 bool keep_it = (fabs(cth) >= cos_theta_cut);
928 if (keep_it)
929 {
930 // almFromPM(pixNumber.NElts(), nlmax, nmmax, phi0, domega, theta, dataq, datau, alme, almb);
931 almFromWX(nlmax, nmmax, phi0, domega, theta, dataq, datau, alme, almb);
932 }
933 }
934}
935
936
937 /*! \fn void SOPHYA::SphericalTransformServer::almFromWX(int_4 nlmax, int_4 nmmax,
938 r_8 phi0, r_8 domega,
939 r_8 theta,
940 const TVector<T>& dataq,
941 const TVector<T>& datau,
942 Alm<T>& alme,
943 Alm<T>& almb) const
944
945Compute polarized Alm's as :
946\f[
947a_{lm}^E=\frac{1}{\sqrt{2}}\sum_{slices}{\omega_{pix}\left(\,_{w}\lambda_l^m\tilde{Q}-i\,_{x}\lambda_l^m\tilde{U}\right)}
948\f]
949\f[
950a_{lm}^B=\frac{1}{\sqrt{2}}\sum_{slices}{\omega_{pix}\left(i\,_{x}\lambda_l^m\tilde{Q}+\,_{w}\lambda_l^m\tilde{U}\right)}
951\f]
952
953where \f$\tilde{Q}\f$ and \f$\tilde{U}\f$ are C-coefficients computed by FFT (method CFromFourierAnalysis, called by present method) from the Stokes parameters.
954
955\f$\omega_{pix}\f$ are solid angle of each pixel.
956
957dataq, datau : Stokes parameters.
958
959 */
960template<class T>
961void SphericalTransformServer<T>::almFromWX(int_4 nlmax, int_4 nmmax,
962 r_8 phi0, r_8 domega,
963 r_8 theta,
964 const TVector<T>& dataq,
965 const TVector<T>& datau,
966 Alm<T>& alme,
967 Alm<T>& almb) const
968{
969 TVector< complex<T> > phaseq(nmmax+1);
970 TVector< complex<T> > phaseu(nmmax+1);
971 // TVector<complex<T> > datain(nph);
972 for (int i=0;i< nmmax+1;i++)
973 {
974 phaseq(i)=0;
975 phaseu(i)=0;
976 }
977 // for(int kk=0; kk<nph; kk++) datain(kk)=complex<T>(dataq(kk),0.);
978
979 // phaseq = CFromFourierAnalysis(nmmax,dataq,phi0);
980 CFromFourierAnalysis(nmmax,dataq,phaseq, phi0);
981
982 // phaseu= CFromFourierAnalysis(nmmax,datau,phi0);
983 CFromFourierAnalysis(nmmax,datau,phaseu, phi0);
984
985 LambdaWXBuilder lwxb(theta,nlmax,nmmax);
986
987 r_8 sqr2inv=1/Rac2;
988 for (int m = 0; m <= nmmax; m++)
989 {
990 r_8 lambda_w=0.;
991 r_8 lambda_x=0.;
992 lwxb.lam_wx(m, m, lambda_w, lambda_x);
993 complex<T> zi_lam_x((T)0., (T)lambda_x);
994 alme(m,m) += ( (T)(lambda_w)*phaseq(m)-zi_lam_x*phaseu(m) )*(T)(domega*sqr2inv);
995 almb(m,m) += ( (T)(lambda_w)*phaseu(m)+zi_lam_x*phaseq(m) )*(T)(domega*sqr2inv);
996
997 for (int l = m+1; l<= nlmax; l++)
998 {
999 lwxb.lam_wx(l, m, lambda_w, lambda_x);
1000 zi_lam_x = complex<T>((T)0., (T)lambda_x);
1001 alme(l,m) += ( (T)(lambda_w)*phaseq(m)-zi_lam_x*phaseu(m) )*(T)(domega*sqr2inv);
1002 almb(l,m) += ( (T)(lambda_w)*phaseu(m)+zi_lam_x*phaseq(m) )*(T)(domega*sqr2inv);
1003 }
1004 }
1005}
1006
1007
1008 /*! \fn void SOPHYA::SphericalTransformServer::almFromPM(int_4 nph, int_4 nlmax,
1009 int_4 nmmax,
1010 r_8 phi0, r_8 domega,
1011 r_8 theta,
1012 const TVector<T>& dataq,
1013 const TVector<T>& datau,
1014 Alm<T>& alme,
1015 Alm<T>& almb) const
1016
1017Compute polarized Alm's as :
1018\f[
1019a_{lm}^E=-\frac{1}{2}\sum_{slices}{\omega_{pix}\left(\,_{+}\lambda_l^m\tilde{P^+}+\,_{-}\lambda_l^m\tilde{P^-}\right)}
1020\f]
1021\f[
1022a_{lm}^B=\frac{i}{2}\sum_{slices}{\omega_{pix}\left(\,_{+}\lambda_l^m\tilde{P^+}-\,_{-}\lambda_l^m\tilde{P^-}\right)}
1023\f]
1024
1025where \f$\tilde{P^{\pm}}=\tilde{Q}\pm\tilde{U}\f$ computed by FFT (method CFromFourierAnalysis, called by present method) from the Stokes parameters,\f$Q\f$ and \f$U\f$ .
1026
1027\f$\omega_{pix}\f$ are solid angle of each pixel.
1028
1029dataq, datau : Stokes parameters.
1030
1031 */
1032template<class T>
1033void SphericalTransformServer<T>::almFromPM(int_4 nph, int_4 nlmax,
1034 int_4 nmmax,
1035 r_8 phi0, r_8 domega,
1036 r_8 theta,
1037 const TVector<T>& dataq,
1038 const TVector<T>& datau,
1039 Alm<T>& alme,
1040 Alm<T>& almb) const
1041{
1042 TVector< complex<T> > phasep(nmmax+1);
1043 TVector< complex<T> > phasem(nmmax+1);
1044 TVector<complex<T> > datain(nph);
1045 for (int i=0;i< nmmax+1;i++)
1046 {
1047 phasep(i)=0;
1048 phasem(i)=0;
1049 }
1050 int kk;
1051 for(kk=0; kk<nph; kk++) datain(kk)=complex<T>(dataq(kk),datau(kk));
1052
1053 phasep = CFromFourierAnalysis(nmmax,datain,phi0);
1054
1055 for(kk=0; kk<nph; kk++) datain(kk)=complex<T>(dataq(kk),-datau(kk));
1056 phasem = CFromFourierAnalysis(nmmax,datain,phi0);
1057 LambdaPMBuilder lpmb(theta,nlmax,nmmax);
1058
1059 for (int m = 0; m <= nmmax; m++)
1060 {
1061 r_8 lambda_p=0.;
1062 r_8 lambda_m=0.;
1063 complex<T> im((T)0.,(T)1.);
1064 lpmb.lam_pm(m, m, lambda_p, lambda_m);
1065
1066 alme(m,m) += -( (T)(lambda_p)*phasep(m) + (T)(lambda_m)*phasem(m) )*(T)(domega*0.5);
1067 almb(m,m) += im*( (T)(lambda_p)*phasep(m) - (T)(lambda_m)*phasem(m) )*(T)(domega*0.5);
1068 for (int l = m+1; l<= nlmax; l++)
1069 {
1070 lpmb.lam_pm(l, m, lambda_p, lambda_m);
1071 alme(l,m) += -( (T)(lambda_p)*phasep(m) + (T)(lambda_m)*phasem(m) )*(T)(domega*0.5);
1072 almb(l,m) += im* ( (T)(lambda_p)*phasep(m) - (T)(lambda_m)*phasem(m) )*(T)(domega*0.5);
1073 }
1074 }
1075}
1076
1077
1078/*! \fn void SOPHYA::SphericalTransformServer::mapFromWX(int_4 nlmax, int_4 nmmax,
1079 SphericalMap<T>& mapq,
1080 SphericalMap<T>& mapu,
1081 const Alm<T>& alme,
1082 const Alm<T>& almb, bool healpix) const
1083
1084synthesis of Stokes parameters following formulae :
1085
1086\f[
1087Q=\sum_{m=-mmax}^{mmax}b_m^qe^{im\varphi}
1088\f]
1089\f[
1090U=\sum_{m=-mmax}^{mmax}b_m^ue^{im\varphi}
1091\f]
1092
1093computed by FFT (method fourierSynthesisFromB called by the present one)
1094
1095with :
1096
1097\f[
1098b_m^q=-\frac{1}{\sqrt{2}}\sum_{l=|m|}^{lmax}{\left(\,_{w}\lambda_l^ma_{lm}^E-i\,_{x}\lambda_l^ma_{lm}^B\right) }
1099\f]
1100\f[
1101b_m^u=\frac{1}{\sqrt{2}}\sum_{l=|m|}^{lmax}{\left(i\,_{x}\lambda_l^ma_{lm}^E+\,_{w}\lambda_l^ma_{lm}^B\right) }
1102\f]
1103 */
1104template<class T>
1105void SphericalTransformServer<T>::mapFromWX(int_4 nlmax, int_4 nmmax,
1106 SphericalMap<T>& mapq,
1107 SphericalMap<T>& mapu,
1108 const Alm<T>& alme,
1109 const Alm<T>& almb, bool healpix) const
1110{
1111 int i;
1112
1113 Bm<complex<T> > b_m_theta_q(nmmax);
1114 Bm<complex<T> > b_m_theta_u(nmmax);
1115
1116 for (int_4 ith = 0; ith < mapq.NbThetaSlices();ith++)
1117 {
1118 int_4 nph;
1119 r_8 phi0;
1120 r_8 theta;
1121 TVector<int_4> pixNumber;
1122 TVector<T> datan;
1123
1124 mapq.GetThetaSlice(ith,theta,phi0, pixNumber,datan);
1125 nph = pixNumber.NElts();
1126 // -----------------------------------------------------
1127 // for each theta, and each m, computes
1128 // b(m,theta) = sum_over_l>m (lambda_l_m(theta) * a_l_m)
1129 // ------------------------------------------------------
1130 LambdaWXBuilder lwxb(theta,nlmax,nmmax);
1131 // LambdaPMBuilder lpmb(theta,nlmax,nmmax);
1132 r_8 sqr2inv=1/Rac2;
1133 int m;
1134 for (m = 0; m <= nmmax; m++)
1135 {
1136 r_8 lambda_w=0.;
1137 r_8 lambda_x=0.;
1138 lwxb.lam_wx(m, m, lambda_w, lambda_x);
1139 complex<T> zi_lam_x((T)0., (T)lambda_x);
1140
1141 b_m_theta_q(m) = ( (T)(lambda_w) * alme(m,m) - zi_lam_x * almb(m,m))*(T)sqr2inv ;
1142 b_m_theta_u(m) = ( (T)(lambda_w) * almb(m,m) + zi_lam_x * alme(m,m))*(T)sqr2inv;
1143
1144
1145 for (int l = m+1; l<= nlmax; l++)
1146 {
1147
1148 lwxb.lam_wx(l, m, lambda_w, lambda_x);
1149 zi_lam_x= complex<T>((T)0., (T)lambda_x);
1150
1151 b_m_theta_q(m) += ((T)(lambda_w)*alme(l,m)-zi_lam_x *almb(l,m))*(T)sqr2inv;
1152 b_m_theta_u(m) += ((T)(lambda_w)*almb(l,m)+zi_lam_x *alme(l,m))*(T)sqr2inv;
1153
1154 }
1155 }
1156 // obtains the negative m of b(m,theta) (= complex conjugate)
1157 for (m=1;m<=nmmax;m++)
1158 {
1159 b_m_theta_q(-m) = conj(b_m_theta_q(m));
1160 b_m_theta_u(-m) = conj(b_m_theta_u(m));
1161 }
1162 if (healpix)
1163 {
1164 TVector<T> Tempq = RfourierSynthesisFromB(b_m_theta_q,nph,phi0);
1165 TVector<T> Tempu = RfourierSynthesisFromB(b_m_theta_u,nph,phi0);
1166 for (i=0;i< nph;i++)
1167 {
1168 mapq(pixNumber(i))=Tempq(i);
1169 mapu(pixNumber(i))=Tempu(i);
1170 }
1171 }
1172 else
1173 // pour des pixelisations quelconques (autres que HEALPix
1174 // nph n'est pas toujours pair
1175 // ca fait des problemes pour les transformees de Fourier
1176 // car le server de TF ajuste la longueur du vecteur reel
1177 // en sortie de TF, bref, la securite veut qu'on prenne une
1178 // TF complexe
1179 {
1180 TVector<complex<T> > Tempq = fourierSynthesisFromB(b_m_theta_q,nph,phi0);
1181 TVector<complex<T> > Tempu = fourierSynthesisFromB(b_m_theta_u,nph,phi0);
1182 for (i=0;i< nph;i++)
1183 {
1184 mapq(pixNumber(i))=Tempq(i).real();
1185 mapu(pixNumber(i))=Tempu(i).real();
1186 }
1187 }
1188 }
1189}
1190/*! \fn void SOPHYA::SphericalTransformServer::mapFromPM(int_4 nlmax, int_4 nmmax,
1191 SphericalMap<T>& mapq,
1192 SphericalMap<T>& mapu,
1193 const Alm<T>& alme,
1194 const Alm<T>& almb) const
1195
1196synthesis of polarizations following formulae :
1197
1198\f[
1199P^+ = \sum_{m=-mmax}^{mmax} {b_m^+e^{im\varphi} }
1200\f]
1201\f[
1202P^- = \sum_{m=-mmax}^{mmax} {b_m^-e^{im\varphi} }
1203\f]
1204
1205computed by FFT (method fourierSynthesisFromB called by the present one)
1206
1207with :
1208
1209\f[
1210b_m^+=-\sum_{l=|m|}^{lmax}{\,_{+}\lambda_l^m \left( a_{lm}^E+ia_{lm}^B \right) }
1211\f]
1212\f[
1213b_m^-=-\sum_{l=|m|}^{lmax}{\,_{+}\lambda_l^m \left( a_{lm}^E-ia_{lm}^B \right) }
1214\f]
1215 */
1216template<class T>
1217void SphericalTransformServer<T>::mapFromPM(int_4 nlmax, int_4 nmmax,
1218 SphericalMap<T>& mapq,
1219 SphericalMap<T>& mapu,
1220 const Alm<T>& alme,
1221 const Alm<T>& almb) const
1222{
1223 Bm<complex<T> > b_m_theta_p(nmmax);
1224 Bm<complex<T> > b_m_theta_m(nmmax);
1225 for (int_4 ith = 0; ith < mapq.NbThetaSlices();ith++)
1226 {
1227 int_4 nph;
1228 r_8 phi0;
1229 r_8 theta;
1230 TVector<int_4> pixNumber;
1231 TVector<T> datan;
1232
1233 mapq.GetThetaSlice(ith,theta,phi0, pixNumber,datan);
1234 nph = pixNumber.NElts();
1235
1236 // -----------------------------------------------------
1237 // for each theta, and each m, computes
1238 // b(m,theta) = sum_over_l>m (lambda_l_m(theta) * a_l_m)
1239 //------------------------------------------------------
1240
1241 LambdaPMBuilder lpmb(theta,nlmax,nmmax);
1242 int m;
1243 for (m = 0; m <= nmmax; m++)
1244 {
1245 r_8 lambda_p=0.;
1246 r_8 lambda_m=0.;
1247 lpmb.lam_pm(m, m, lambda_p, lambda_m);
1248 complex<T> im((T)0.,(T)1.);
1249
1250 b_m_theta_p(m) = (T)(lambda_p )* (-alme(m,m) - im * almb(m,m));
1251 b_m_theta_m(m) = (T)(lambda_m) * (-alme(m,m) + im * almb(m,m));
1252
1253
1254 for (int l = m+1; l<= nlmax; l++)
1255 {
1256 lpmb.lam_pm(l, m, lambda_p, lambda_m);
1257 b_m_theta_p(m) += (T)(lambda_p)*(-alme(l,m)-im *almb(l,m));
1258 b_m_theta_m(m) += (T)(lambda_m)*(-alme(l,m)+im *almb(l,m));
1259 }
1260 }
1261
1262 // obtains the negative m of b(m,theta) (= complex conjugate)
1263 for (m=1;m<=nmmax;m++)
1264 {
1265 b_m_theta_p(-m) = conj(b_m_theta_m(m));
1266 b_m_theta_m(-m) = conj(b_m_theta_p(m));
1267 }
1268
1269 TVector<complex<T> > Tempp = fourierSynthesisFromB(b_m_theta_p,nph,phi0);
1270 TVector<complex<T> > Tempm = fourierSynthesisFromB(b_m_theta_m,nph,phi0);
1271
1272 for (int i=0;i< nph;i++)
1273 {
1274 mapq(pixNumber(i))=0.5*(Tempp(i)+Tempm(i)).real();
1275 mapu(pixNumber(i))=0.5*(Tempp(i)-Tempm(i)).imag();
1276 }
1277 }
1278}
1279
1280
1281 /*! \fn void SOPHYA::SphericalTransformServer::GenerateFromCl(SphericalMap<T>& sphq,
1282 SphericalMap<T>& sphu,
1283 int_4 pixelSizeIndex,
1284 const TVector<T>& Cle,
1285 const TVector<T>& Clb,
1286 const r_8 fwhm) const
1287
1288synthesis of a polarization map from power spectra electric-Cl and magnetic-Cl (Alm's are generated randomly, following a gaussian distribution).
1289 \param fwhm FWHM in arcmin for random generation of Alm's (eg. 5)
1290*/
1291template<class T>
1292void SphericalTransformServer<T>::GenerateFromCl(SphericalMap<T>& sphq,
1293 SphericalMap<T>& sphu,
1294 int_4 pixelSizeIndex,
1295 const TVector<T>& Cle,
1296 const TVector<T>& Clb,
1297 const r_8 fwhm) const
1298{
1299 if (Cle.NElts() != Clb.NElts())
1300 {
1301 cout << " SphericalTransformServer: les deux tableaux Cl n'ont pas la meme taille" << endl;
1302 throw SzMismatchError("SphericalTransformServer::GenerateFromCl : two Cl arrays have not same size");
1303 }
1304
1305 // Alm<T> a2lme,a2lmb;
1306 // almFromCl(a2lme, Cle, fwhm);
1307 // almFromCl(a2lmb, Clb, fwhm);
1308 // Alm<T> a2lme = almFromCl(Cle, fwhm);
1309 // Alm<T> a2lmb = almFromCl(Clb, fwhm);
1310 Alm<T> a2lme(Cle, fwhm);
1311 Alm<T> a2lmb(Clb, fwhm);
1312
1313 GenerateFromAlm(sphq,sphu,pixelSizeIndex,a2lme,a2lmb);
1314}
1315 /*! \fn void SOPHYA::SphericalTransformServer::GenerateFromCl(SphericalMap<T>& sph,
1316 int_4 pixelSizeIndex,
1317 const TVector<T>& Cl,
1318 const r_8 fwhm) const
1319
1320synthesis of a temperature map from power spectrum Cl (Alm's are generated randomly, following a gaussian distribution). */
1321template<class T>
1322void SphericalTransformServer<T>::GenerateFromCl(SphericalMap<T>& sph,
1323 int_4 pixelSizeIndex,
1324 const TVector<T>& Cl,
1325 const r_8 fwhm) const
1326{
1327
1328 Alm<T> alm(Cl, fwhm);
1329 GenerateFromAlm(sph,pixelSizeIndex, alm );
1330}
1331
1332
1333
1334/*! \fn TVector<T> SOPHYA::SphericalTransformServer::DecomposeToCl(SphericalMap<T>& sph, int_4 nlmax, r_8 cos_theta_cut, int iterationOrder) const
1335
1336\return power spectrum from analysis of a temperature map. THE MAP CAN BE MODIFIED (if iterationOrder >0)
1337
1338 \param<nlmax> : maximum value of the l index
1339
1340 \param<cos_theta_cut> : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
1341
1342\param<iterationOrder> : 1,2,3,4.... order of an iterative analysis. If iterationOrder is not null, the method works with SphereHEALPix but NOT WITH SphereThetaPhi maps !
1343
1344 */
1345template <class T>
1346TVector<T> SphericalTransformServer<T>::DecomposeToCl(SphericalMap<T>& sph, int_4 nlmax, r_8 cos_theta_cut, int iterationOrder) const
1347{
1348 Alm<T> alm;
1349 DecomposeToAlm( sph, alm, nlmax, cos_theta_cut, iterationOrder);
1350 // power spectrum
1351 return alm.powerSpectrum();
1352}
1353
1354
1355/*! \fn TVector<T> SOPHYA::SphericalTransformServer::DecomposeToCl(const SphericalMap<T>& sph, int_4 nlmax, r_8 cos_theta_cut) const
1356
1357\return power spectrum from analysis of a temperature map.
1358
1359 \param<nlmax> : maximum value of the l index
1360
1361 \param<cos_theta_cut> : cosinus of the symmetric cut EULER angle theta : cos_theta_cut=0 means no cut ; cos_theta_cut=1 all the sphere is cut.
1362
1363
1364 */
1365
1366
1367template <class T>
1368TVector<T> SphericalTransformServer<T>::DecomposeToCl(const SphericalMap<T>& sph, int_4 nlmax, r_8 cos_theta_cut) const
1369{
1370 Alm<T> alm;
1371 DecomposeToAlm( sph, alm, nlmax, cos_theta_cut);
1372 // power spectrum
1373 return alm.powerSpectrum();
1374}
1375
1376#ifdef __CXX_PRAGMA_TEMPLATES__
1377#pragma define_template SphericalTransformServer<r_8>
1378#pragma define_template SphericalTransformServer<r_4>
1379#endif
1380#if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES)
1381template class SOPHYA::SphericalTransformServer<r_8>;
1382template class SOPHYA::SphericalTransformServer<r_4>;
1383#endif
Note: See TracBrowser for help on using the repository browser.