Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
alm.hh
Go to the documentation of this file.
1 /*! \file alm.hh
2  * Class for storing spherical harmonic coefficients.
3  *
4  * alm wat interface to healpix alm.h (aligned with healpix 3.31)
5  */
6 
7 #ifndef ALM_HH
8 #define ALM_HH
9 
10 #ifndef __CINT__
11 #include "alm.h"
12 #include "xcomplex.h" // starting from healpix 3.30 the complex class is the std::complex !!!
13 #include "alm_powspec_tools.h"
14 #endif
15 #include <complex>
16 #include "TMath.h"
17 
18 using namespace std;
19 
20 namespace wat {
21 
22 
23 /*! Base class for calculating the storage layout of spherical harmonic
24  coefficients. */
25 class Alm_Base
26  {
27  protected:
28  int lmax, mmax, tval;
29 
30  public:
31  /*! Constructs an Alm_Base object with given \a lmax and \a mmax. */
32  Alm_Base (int lmax_=0, int mmax_=0)
33  : lmax(lmax_), mmax(mmax_), tval(2*lmax+1) {}
34 
35  /*! Returns the total number of coefficients for maximum quantum numbers
36  \a l and \a m. */
37  int Num_Alms (int l, int m)
38  {
39  planck_assert(m<=l,"mmax must not be larger than lmax");
40  return ((m+1)*(m+2))/2 + (m+1)*(l-m);
41  }
42 
43  /*! Changes the object's maximum quantum numbers to \a lmax and \a mmax. */
44  void Set (int lmax_, int mmax_)
45  {
46  lmax=lmax_;
47  mmax=mmax_;
48  tval=2*lmax+1;
49  }
50 
51  /*! Returns the maximum \a l */
52  int Lmax() const { return lmax; }
53  /*! Returns the maximum \a m */
54  int Mmax() const { return mmax; }
55 
56  /*! Returns an array index for a given m, from which the index of a_lm
57  can be obtained by adding l. */
58  int index_l0 (int m) const
59  { return ((m*(tval-m))>>1); }
60 
61  /*! Returns the array index of the specified coefficient. */
62  int index (int l, int m) const
63  { return index_l0(m) + l; }
64 
65  /*! Returns \a true, if both objects have the same \a lmax and \a mmax,
66  else \a false. */
67  bool conformable (const Alm_Base &other) const
68  { return ((lmax==other.lmax) && (mmax==other.mmax)); }
69 
70  /*! Swaps the contents of two Alm_Base objects. */
71  void swap (Alm_Base &other);
72  };
73 
74 
75 /*! Class for storing spherical harmonic coefficients. */
76 template<typename T> class Alm_Template: public Alm_Base
77  {
78  private:
79 #ifndef __CINT__
80  arr<T> alm;
81 #endif
82 
83  public:
84  /*! Constructs an Alm_Template object with given \a lmax and \a mmax. */
85  Alm_Template (int lmax_=0, int mmax_=0)
86  : Alm_Base(lmax_,mmax_), alm (Num_Alms(lmax,mmax)) {}
87 
88  /*! Deletes the old coefficients and allocates storage according to
89  \a lmax and \a mmax. */
90  void Set (int lmax_, int mmax_)
91  {
92  Alm_Base::Set(lmax_, mmax_);
93  alm.alloc(Num_Alms(lmax,mmax));
94  }
95 
96  /*! Deallocates the old coefficients and uses the content of \a data
97  for storage. \a data is deallocated during the call. */
98 #ifndef __CINT__
99  void Set (arr<T> &data, int lmax_, int mmax_)
100  {
101  planck_assert (Num_Alms(lmax_,mmax_)==data.size(),"wrong array size");
102  Alm_Base::Set(lmax_, mmax_);
103  alm.transfer(data);
104  }
105 #endif
106 
107  /*! Sets all coefficients to zero. */
108  void SetToZero ()
109  { alm.fill (0); }
110 
111  /*! Multiplies all coefficients by \a factor. */
112  template<typename T2> void Scale (const T2 &factor)
113  { for (tsize m=0; m<alm.size(); ++m) alm[m]*=factor; }
114 #ifndef __CINT__
115  /*! \a a(l,m) *= \a factor[l] for all \a l,m. */
116  template<typename T2> void ScaleL (const arr<T2> &factor)
117  {
118  planck_assert(factor.size()>tsize(lmax),
119  "alm.ScaleL: factor array too short");
120  for (int m=0; m<=mmax; ++m)
121  for (int l=m; l<=lmax; ++l)
122  operator()(l,m)*=factor[l];
123  }
124  /*! \a a(l,m) *= \a factor[m] for all \a l,m. */
125  template<typename T2> void ScaleM (const arr<T2> &factor)
126  {
127  planck_assert(factor.size()>tsize(mmax),
128  "alm.ScaleM: factor array too short");
129  for (int m=0; m<=mmax; ++m)
130  for (int l=m; l<=lmax; ++l)
131  operator()(l,m)*=factor[m];
132  }
133 #endif
134  /*! Adds \a num to a_00. */
135  template<typename T2> void Add (const T2 &num)
136  { alm[0]+=num; }
137 
138  /*! Returns a reference to the specified coefficient. */
139  T &operator() (int l, int m)
140  { return alm[index(l,m)]; }
141  /*! Returns a constant reference to the specified coefficient. */
142  const T &operator() (int l, int m) const
143  { return alm[index(l,m)]; }
144 
145  /*! Returns a pointer for a given m, from which the address of a_lm
146  can be obtained by adding l. */
147  T *mstart (int m)
148  { return &alm[index_l0(m)]; }
149  /*! Returns a pointer for a given m, from which the address of a_lm
150  can be obtained by adding l. */
151  const T *mstart (int m) const
152  { return &alm[index_l0(m)]; }
153 
154  /*! Returns a constant reference to the a_lm data. */
155 #ifndef __CINT__
156  const arr<T> &Alms () const { return alm; }
157 #endif
158 
159  /*! Swaps the contents of two Alm_Template objects. */
160  void swap (Alm_Template &other)
161  {
162  Alm_Base::swap(other);
163  alm.swap(other.alm);
164  }
165 
166  /*! Adds all coefficients from \a other to the own coefficients. */
167  void Add (const Alm_Template &other)
168  {
169  planck_assert (conformable(other), "A_lm are not conformable");
170  for (tsize m=0; m<alm.size(); ++m)
171  alm[m] += other.alm[m];
172  }
173  };
174 
175  class Alm: public Alm_Template<complex<double> >
176  {
177  public:
178  Alm (int lmax_=0, int mmax_=0)
179  : Alm_Template<complex<double> >(lmax_,mmax_) {}
180 
181  //: Copy constructors
182  Alm(const Alm& alm) {*this = alm;}
183 
184 #ifndef __CINT__
185  Alm(const ::Alm<xcomplex<double> >& alm) {*this = alm;}
186 #endif
187 
188  // applies gaussian smoothing to alm
189  void smoothWithGauss(double fwhm) {
190 #ifndef __CINT__
191  double deg2rad = PI/180.;
192  fwhm*=deg2rad;
193 
194  ::Alm<xcomplex<double> > alm(Lmax(),Mmax());
195  for(int l=0;l<=Lmax();l++) {
196  int limit = TMath::Min(l,Mmax());
197  for (int m=0; m<=limit; m++)
198  alm(l,m)=complex<double>((*this)(l,m).real(),(*this)(l,m).imag());
199  }
200  //::Alm<xcomplex<double> > alm = *this;
201  ::smoothWithGauss(alm, fwhm);
202  *this = alm;
203 #endif
204  }
205 
206 
207  // Rotates alm through the Euler angles psi, theta and phi.
208  // The Euler angle convention is right handed, rotations are active.
209  // - psi is the first rotation about the z-axis (vertical)
210  // - then theta about the ORIGINAL (unrotated) y-axis
211  // - then phi about the ORIGINAL (unrotated) z-axis (vertical)
212  // relates Alm */
213  void rotate(double psi, double theta, double phi) {
214 #ifndef __CINT__
215  ::Alm<xcomplex<double> > alm(Lmax(),Mmax());
216  for(int l=0;l<=Lmax();l++) {
217  int limit = TMath::Min(l,Mmax());
218  for (int m=0; m<=limit; m++)
219  alm(l,m)=complex<double>((*this)(l,m).real(),(*this)(l,m).imag());
220  }
221  rotate_alm(alm, psi, theta, phi);
222  *this = alm;
223 #endif
224  }
225 #ifndef __CINT__
226  void operator=(const ::Alm<xcomplex<double> >& alm) {
227  Set(alm.Lmax(),alm.Mmax());
228  for(int l=0;l<=Lmax();l++) {
229  int limit = TMath::Min(l,Mmax());
230  for (int m=0; m<=limit; m++)
231  (*this)(l,m) = complex<double>( alm(l,m).real(), alm(l,m).imag());
232  }
233  }
234 /*
235  ::Alm<xcomplex<double> >& operator=(const Alm& a) {
236  ::Alm<xcomplex<double> > alm(Lmax(),Mmax());
237  for(int l=0;l<=Lmax();l++) {
238  int limit = TMath::Min(l,Mmax());
239  for (int m=0; m<=limit; m++)
240  alm(l,m)=complex<double>((*this)(l,m).real(),(*this)(l,m).imag());
241  }
242  return alm;
243  }
244 */
245 #endif
246 
247  };
248 
249 } // end namespace
250 
251 #endif
void smoothWithGauss(double fwhm)
Definition: alm.hh:189
void ScaleM(const arr< T2 > &factor)
a(l,m) *= factor[m] for all l,m.
Definition: alm.hh:125
int index_l0(int m) const
Returns an array index for a given m, from which the index of a_lm can be obtained by adding l...
Definition: alm.hh:58
int tval
Definition: alm.hh:28
const arr< T > & Alms() const
Returns a constant reference to the a_lm data.
Definition: alm.hh:156
void Add(const Alm_Template &other)
Adds all coefficients from other to the own coefficients.
Definition: alm.hh:167
bool conformable(const Alm_Base &other) const
Returns true, if both objects have the same lmax and mmax, else false.
Definition: alm.hh:67
int Lmax() const
Returns the maximum l.
Definition: alm.hh:52
void operator=(const ::Alm< xcomplex< double > > &alm)
Definition: alm.hh:226
void Set(int lmax_, int mmax_)
Deletes the old coefficients and allocates storage according to lmax and mmax.
Definition: alm.hh:90
float theta
STL namespace.
int m
Definition: cwb_net.C:10
void Set(arr< T > &data, int lmax_, int mmax_)
Deallocates the old coefficients and uses the content of data for storage.
Definition: alm.hh:99
Definition: alm.hh:20
#define PI
Definition: watfun.hh:14
void Set(int lmax_, int mmax_)
Changes the object's maximum quantum numbers to lmax and mmax.
Definition: alm.hh:44
void swap(Alm_Template &other)
Swaps the contents of two Alm_Template objects.
Definition: alm.hh:160
void Add(const T2 &num)
Adds num to a_00.
Definition: alm.hh:135
float phi
double factor
void SetToZero()
Sets all coefficients to zero.
Definition: alm.hh:108
float psi
int Mmax() const
Returns the maximum m.
Definition: alm.hh:54
Alm(const ::Alm< xcomplex< double > > &alm)
Definition: alm.hh:185
const T * mstart(int m) const
Returns a pointer for a given m, from which the address of a_lm can be obtained by adding l...
Definition: alm.hh:151
double deg2rad
Alm_Template(int lmax_=0, int mmax_=0)
Constructs an Alm_Template object with given lmax and mmax.
Definition: alm.hh:85
int mmax
Definition: alm.hh:28
arr< T > alm
Definition: alm.hh:80
Alm(int lmax_=0, int mmax_=0)
Definition: alm.hh:178
T * mstart(int m)
Returns a pointer for a given m, from which the address of a_lm can be obtained by adding l...
Definition: alm.hh:147
int index(int l, int m) const
Returns the array index of the specified coefficient.
Definition: alm.hh:62
int num
void ScaleL(const arr< T2 > &factor)
a(l,m) *= factor[l] for all l,m.
Definition: alm.hh:116
Alm_Base(int lmax_=0, int mmax_=0)
Constructs an Alm_Base object with given lmax and mmax.
Definition: alm.hh:32
Base class for calculating the storage layout of spherical harmonic coefficients. ...
Definition: alm.hh:25
Definition: alm.hh:175
double T
Definition: testWDM_4.C:11
void rotate(double psi, double theta, double phi)
Definition: alm.hh:213
int Num_Alms(int l, int m)
Returns the total number of coefficients for maximum quantum numbers l and m.
Definition: alm.hh:37
wavearray< int > index
Class for storing spherical harmonic coefficients.
Definition: alm.hh:76
Alm(const Alm &alm)
Definition: alm.hh:182
int l
Definition: cbc_plots.C:434
double mmax
void Scale(const T2 &factor)
Multiplies all coefficients by factor.
Definition: alm.hh:112
int lmax
Definition: alm.hh:28