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

Last change on this file since 3382 was 3077, checked in by cmv, 19 years ago

remplacement nbrandom.h (obsolete) -> srandgen.h cmv 14/09/2006

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