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

Last change on this file since 3725 was 3613, checked in by ansari, 16 years ago

Adaptation a l'introduction de la suite des classes RandomGeneratorInterface et Cie , Reza 30/04/2009

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