1 | #include "sopnamsp.h"
|
---|
2 | #include "machdefs.h"
|
---|
3 | #include <stdio.h>
|
---|
4 | #include <stdlib.h>
|
---|
5 | #include <iostream>
|
---|
6 | #include <math.h>
|
---|
7 |
|
---|
8 | #include <typeinfo>
|
---|
9 | #include <vector>
|
---|
10 | #include <string>
|
---|
11 |
|
---|
12 | #include "piacmd.h"
|
---|
13 | #include "nobjmgr.h"
|
---|
14 | #include "pistdimgapp.h"
|
---|
15 | #include "servnobjm.h"
|
---|
16 | #include "pidrawer.h"
|
---|
17 | #include "nomgadapter.h"
|
---|
18 |
|
---|
19 | #include "tvector.h"
|
---|
20 | #include "ntuple.h"
|
---|
21 | #include "slinparbuff.h"
|
---|
22 | #include "localmap.h"
|
---|
23 | #include "spherehealpix.h"
|
---|
24 | #include "spherethetaphi.h"
|
---|
25 | #include "sphericaltransformserver.h"
|
---|
26 |
|
---|
27 | extern "C" {
|
---|
28 | void skymapmodule_init();
|
---|
29 | void skymapmodule_end();
|
---|
30 | }
|
---|
31 |
|
---|
32 | class skymapmoduleExecutor : public CmdExecutor {
|
---|
33 | public:
|
---|
34 | enum {UnKnown=(uint_2) 0,
|
---|
35 | TypeR8=(uint_2) 16, TypeR4=(uint_2) 32,
|
---|
36 | Spherical=(uint_2) 1, HealPix=(uint_2) 2, ThetaPhi=(uint_2) 4,
|
---|
37 | Spherical8=(uint_2) Spherical|TypeR8,
|
---|
38 | Spherical4=(uint_2) Spherical|TypeR4,
|
---|
39 | HealPix8 =(uint_2) Spherical|HealPix|TypeR8,
|
---|
40 | HealPix4 =(uint_2) Spherical|HealPix|TypeR4,
|
---|
41 | ThetaPhi8 =(uint_2) Spherical|ThetaPhi|TypeR8,
|
---|
42 | ThetaPhi4 =(uint_2) Spherical|ThetaPhi|TypeR4};
|
---|
43 |
|
---|
44 | skymapmoduleExecutor();
|
---|
45 | virtual ~skymapmoduleExecutor();
|
---|
46 | virtual int Execute(string& keyw, vector<string>& args, string& toks);
|
---|
47 |
|
---|
48 | void Resol2SizeIndex(vector<string>& tokens);
|
---|
49 | void TypeMap(vector<string>& tokens);
|
---|
50 | void Map2Double(vector<string>& tokens);
|
---|
51 | void Map2Float(vector<string>& tokens);
|
---|
52 | void Map2Map(vector<string>& tokens);
|
---|
53 | void MapMult(vector<string>& tokens);
|
---|
54 | void MapCover(vector<string>& tokens);
|
---|
55 | void Map2Cl(vector<string>& tokens);
|
---|
56 | void Cl2Map(vector<string>& tokens);
|
---|
57 | void Map2Alm(vector<string>& tokens);
|
---|
58 | void Alm2Map(vector<string>& tokens);
|
---|
59 | void CrMapMask(vector<string>& tokens);
|
---|
60 | void CrMaskFrMap(vector<string>& tokens);
|
---|
61 | void MaskMap(vector<string>& tokens);
|
---|
62 | void MapProj(vector<string>& tokens);
|
---|
63 | void Map2Local(vector<string>& tokens);
|
---|
64 | void MapOper(vector<string>& tokens);
|
---|
65 | void MapStat(vector<string>& tokens);
|
---|
66 | void Alm2Cl(vector<string>& tokens);
|
---|
67 | void Cl2llCl(vector<string>& tokens);
|
---|
68 | void ClMean(vector<string>& tokens);
|
---|
69 | void ClMult(vector<string>& tokens);
|
---|
70 | void ClOper(vector<string>& tokens);
|
---|
71 | void ClRebin(vector<string>& tokens);
|
---|
72 | protected:
|
---|
73 | void SetTypeMap(vector<string>& tokens);
|
---|
74 | bool IsNSideGood(int_4 nside);
|
---|
75 | void GetNSideGood(int_4 &nside);
|
---|
76 | uint_2 typemap(AnyDataObj* obj);
|
---|
77 |
|
---|
78 | uint_2 DefTypMap;
|
---|
79 | };
|
---|
80 |
|
---|
81 | /* --Methode-- */
|
---|
82 | skymapmoduleExecutor::skymapmoduleExecutor()
|
---|
83 | {
|
---|
84 | PIACmd * mpiac;
|
---|
85 | NamedObjMgr omg;
|
---|
86 | mpiac = omg.GetImgApp()->CmdInterpreter();
|
---|
87 |
|
---|
88 | string hgrp = "SkyMapCmd";
|
---|
89 | string gdesc = "4 Pi spherical, local maps manipulation commands \n";
|
---|
90 | gdesc += "and spherical harmonics analysis Ylm. \n This is mainly ";
|
---|
91 | gdesc += "Based on SOPHYA SkyMap and Samba modules. ";
|
---|
92 | mpiac->AddHelpGroup(hgrp, gdesc);
|
---|
93 | string kw,usage;
|
---|
94 |
|
---|
95 | kw = "settypemap";
|
---|
96 | usage = "Set le type de map(double) par default";
|
---|
97 | usage += "\n Usage: settypemap type";
|
---|
98 | usage += "\n type= H for Healpix (default)";
|
---|
99 | usage += "\n T for ThetaPhi";
|
---|
100 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
101 |
|
---|
102 | kw = "resol2szidx";
|
---|
103 | usage = "Compute SizeIndex value (=nside for HEALPix) for a";
|
---|
104 | usage += "\n given resolution, (resol in arcminutes)";
|
---|
105 | usage += "\n Usage: resol2szidx resol";
|
---|
106 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
107 |
|
---|
108 | kw = "typemap";
|
---|
109 | usage = "Imprime le type de map";
|
---|
110 | usage += "\n Usage: typemap map";
|
---|
111 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
112 |
|
---|
113 | kw = "map2double";
|
---|
114 | usage = "Convert a float map to a double map";
|
---|
115 | usage += "\n Usage: map2double map(float)";
|
---|
116 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
117 |
|
---|
118 | kw = "map2float";
|
---|
119 | usage = "Convert a double map to a float map";
|
---|
120 | usage += "\n Usage: map2float map(double)";
|
---|
121 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
122 |
|
---|
123 | kw = "map2map";
|
---|
124 | usage = "Convert a map(double) into an other map type";
|
---|
125 | usage += "\n Usage: map2map map(double) type [nside]";
|
---|
126 | usage += "\n type= H for Healpix";
|
---|
127 | usage += "\n T for ThetaPhi";
|
---|
128 | usage += "\n if nside <=1 use nside from map";
|
---|
129 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
130 |
|
---|
131 | kw = "mapmult";
|
---|
132 | usage = "Multiply a map(double) by a constant";
|
---|
133 | usage += "\n Usage: mapmult map(double) cste";
|
---|
134 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
135 |
|
---|
136 | kw = "mapcover";
|
---|
137 | usage = "Print the covering of a map(double)";
|
---|
138 | usage += "\n Usage: mapcover map(double) v1,v2 [varname]";
|
---|
139 | usage += "\n v1,v2: good pixels have v1<=val<=v2 (def=0.99,1.01)";
|
---|
140 | usage += "\n (use default (0.99,1.01): \"!\")";
|
---|
141 | usage += "\n $varname: percent of sky covering";
|
---|
142 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
143 |
|
---|
144 | kw = "map2cl";
|
---|
145 | usage = "SkyMap map(double) to Cl";
|
---|
146 | usage += "\n Usage: map2cl map(double) clvec [lmax] [niter]";
|
---|
147 | usage += "\n lmax: if <=0 or \"!\" use default (3*nside-1)";
|
---|
148 | usage += "\n niter: number of iterations (<=0 or def=3)";
|
---|
149 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
150 |
|
---|
151 | kw = "cl2map";
|
---|
152 | usage = "Generate SkyMap map(double) from Cl";
|
---|
153 | usage += "\n Usage: cl2map clvec map(double) nside [fwhm]";
|
---|
154 | usage += "\n (see command settypemap)";
|
---|
155 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
156 |
|
---|
157 | kw = "map2alm";
|
---|
158 | usage = "SkyMap map(double) to Alm";
|
---|
159 | usage += "\n Usage: map2alm map(double) almmat [lmax] [niter]";
|
---|
160 | usage += "\n lmax: if <=0 or \"!\" use default (3*nside-1)";
|
---|
161 | usage += "\n niter: number of iterations (<=0 or def=3)";
|
---|
162 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
163 |
|
---|
164 | kw = "alm2map";
|
---|
165 | usage = "Generate SkyMap map(double) from Alm";
|
---|
166 | usage += "\n Usage: alm2map almmat map(double) nside";
|
---|
167 | usage += "\n (see command settypemap)";
|
---|
168 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
169 |
|
---|
170 | kw = "crmapmask";
|
---|
171 | usage = "Create a map mask(double) (nside) between latitudes lat1,lat2 longitudes lon1,lon2";
|
---|
172 | usage += "\n Usage: crmapmask mapmsk(double) nside lat1,lat2 lon1,lon2 [valmsk,valnomsk]";
|
---|
173 | usage += "\n mask pixels between [lat1,lat2] ([-90,90] in degrees)";
|
---|
174 | usage += "\n and [lon1,lon2] ([0,360[ in degrees)";
|
---|
175 | usage += "\n (for no mask \"!\")";
|
---|
176 | usage += "\n valmsk=value to be written into masked pixels (def=0)";
|
---|
177 | usage += "\n valnomsk=value to be written into unmasked pixels (def=1)";
|
---|
178 | usage += "\n (see command settypemap)";
|
---|
179 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
180 |
|
---|
181 | kw = "crmaskfrmap";
|
---|
182 | usage = "Create a map mask(double) (nside) relative to an other map(double) pixels values";
|
---|
183 | usage += "\n Usage: crmaskfrmap mapmsk(double) nside map(double) v1,v2 [valmsk,valnomsk]";
|
---|
184 | usage += "\n mask if v1<mapmsk(i)<v2";
|
---|
185 | usage += "\n valmsk=value to be written into masked pixels (def=0)";
|
---|
186 | usage += "\n valnomsk=value to be written into unmasked pixels (def=1)";
|
---|
187 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
188 |
|
---|
189 | kw = "maskmap";
|
---|
190 | usage = "Mask a map(double) with a mask map(double)";
|
---|
191 | usage += "\n Usage: maskmap map msk";
|
---|
192 | usage += "\n operation is map(i) *= msk(theta,phi)";
|
---|
193 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
194 |
|
---|
195 | kw = "maproj";
|
---|
196 | usage = "Project a map(double) into an other one";
|
---|
197 | usage += "\n Usage: maproj map(double) mapproj(double) [nside]";
|
---|
198 | usage += "\n Create mapproj(nside) and project map into mapproj";
|
---|
199 | usage += "\n (use \"maproj map mapproj\" to make a copy)";
|
---|
200 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
201 |
|
---|
202 | kw = "map2local";
|
---|
203 | usage = "Project a map(double) into a local map(double)";
|
---|
204 | usage += "\n Usage: map2local map(double) localmap(double) nx,ny angleX,angleY phi0,theta0 [x0,y0] [angle]";
|
---|
205 | usage += "\n nx,ny: number of pixels in x(col),y(row) direction";
|
---|
206 | usage += "\n X-axis==phi, Y-axis==theta";
|
---|
207 | usage += "\n angleX,angleY: total angle aperture in x,y direction (degrees)";
|
---|
208 | usage += "\n phi0,theta0: define origin in phi,theta (degrees)";
|
---|
209 | usage += "\n x0,y0: map phi0,theta0 to pixel x0,y0 (\"!\" or def=middle of the localmap)";
|
---|
210 | usage += "\n angle: angle (degrees) of rotation between x-axis of map-coordinate system";
|
---|
211 | usage += "\n and the tangent to parallel on the sphere (def: 0.)";
|
---|
212 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
213 |
|
---|
214 | kw = "mapop";
|
---|
215 | usage = "operation on a map(double)";
|
---|
216 | usage += "\n Usage: mapop map(double) op1 map1(double) [op2 map2(double)] [op2 map3(double)] [...]";
|
---|
217 | usage += "\n Do \"map op1= map1\", \"map op2=map2\", ...";
|
---|
218 | usage += "\n where op = +,-,*,/";
|
---|
219 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
220 |
|
---|
221 | kw = "mapstat";
|
---|
222 | usage = "Statistics from a map(double)";
|
---|
223 | usage += "\n Usage: mapstat map(double) [msk(double)] [meanvarname] [sigmavarname]";
|
---|
224 | usage += "\n msk = mask sphere used as weight";
|
---|
225 | usage += "\n if msk(i)<=0 do not use that pixel";
|
---|
226 | usage += "\n if msk(i)>0 use that pixel with weight msk(i)";
|
---|
227 | usage += "\n msk = \"!\" means no mask sphere";
|
---|
228 | usage += "\n meanvarname = variable name to store mean";
|
---|
229 | usage += "\n sigmavarname = variable name to store sigma";
|
---|
230 | usage += "\n (\"!\" means do not store)";
|
---|
231 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
232 |
|
---|
233 | kw = "alm2cl";
|
---|
234 | usage = "Project alm to cl";
|
---|
235 | usage += "\n Usage: alm2cl alm cl";
|
---|
236 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
237 |
|
---|
238 | kw = "cl2llcl";
|
---|
239 | usage = "Project cl to l*(l+1)*cl";
|
---|
240 | usage += "\n Usage: cl2llcl cl llcl";
|
---|
241 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
242 |
|
---|
243 | kw = "clmean";
|
---|
244 | usage = "Compute mean for a cl vector";
|
---|
245 | usage += "\n Usage: clmean cl imin,imax [meanvarname]";
|
---|
246 | usage += "\n imin,imax = compute between these indices";
|
---|
247 | usage += "\n meanvarname = variable name to store mean";
|
---|
248 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
249 |
|
---|
250 | kw = "clmult";
|
---|
251 | usage = "Multiply a cl by a constant";
|
---|
252 | usage += "\n Usage: clmult cl val";
|
---|
253 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
254 |
|
---|
255 | kw = "clop";
|
---|
256 | usage = "operation on a cl";
|
---|
257 | usage += "\n Usage: clop cl(double) op1 cl1(double) [op2 cl2(double)] [op2 cl3(double)] [...]";
|
---|
258 | usage += "\n Do \"cl op1= cl1\", \"cl op2=cl2\", ...";
|
---|
259 | usage += "\n where op = +,-,*,/";
|
---|
260 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
261 |
|
---|
262 | kw = "clrebin";
|
---|
263 | usage = "rebin a cl into an ntuple";
|
---|
264 | usage += "\n Usage: clrebin cl ntuple nbin,bin0";
|
---|
265 | usage += "\n nbin: rebin by nbin";
|
---|
266 | usage += "\n bin0: begin rebinning at bin bin0 (def=0)";
|
---|
267 | mpiac->RegisterCommand(kw, usage, this, hgrp);
|
---|
268 |
|
---|
269 | DefTypMap = HealPix;
|
---|
270 | }
|
---|
271 |
|
---|
272 | /* --Methode-- */
|
---|
273 | skymapmoduleExecutor::~skymapmoduleExecutor()
|
---|
274 | {
|
---|
275 | }
|
---|
276 |
|
---|
277 | /* --Methode-- */
|
---|
278 | int skymapmoduleExecutor::Execute(string& kw, vector<string>& tokens, string& toks)
|
---|
279 | {
|
---|
280 |
|
---|
281 | if(kw == "settypemap") {
|
---|
282 | SetTypeMap(tokens);
|
---|
283 | } else if(kw == "resol2szidx") {
|
---|
284 | if(tokens.size()<1) {
|
---|
285 | cout<<"Usage: resol2szidx resol"<<endl;
|
---|
286 | return(0);
|
---|
287 | }
|
---|
288 | Resol2SizeIndex(tokens);
|
---|
289 | } else if(kw == "typemap") {
|
---|
290 | if(tokens.size()<1) {
|
---|
291 | cout<<"Usage: typemap map"<<endl;
|
---|
292 | return(0);
|
---|
293 | }
|
---|
294 | TypeMap(tokens);
|
---|
295 | } else if(kw == "map2double") {
|
---|
296 | if(tokens.size()<1) {
|
---|
297 | cout<<"Usage: map2double map"<<endl;
|
---|
298 | return(0);
|
---|
299 | }
|
---|
300 | Map2Double(tokens);
|
---|
301 | } else if(kw == "map2float") {
|
---|
302 | if(tokens.size()<1) {
|
---|
303 | cout<<"Usage: map2float map"<<endl;
|
---|
304 | return(0);
|
---|
305 | }
|
---|
306 | Map2Float(tokens);
|
---|
307 | } else if(kw == "map2map") {
|
---|
308 | if(tokens.size()<2) {
|
---|
309 | cout<<"Usage: map2map map type [nside]"<<endl;
|
---|
310 | return(0);
|
---|
311 | }
|
---|
312 | Map2Map(tokens);
|
---|
313 | } else if(kw == "mapmult") {
|
---|
314 | if(tokens.size()<2) {
|
---|
315 | cout<<"Usage: mapmult map cste"<<endl;
|
---|
316 | return(0);
|
---|
317 | }
|
---|
318 | MapMult(tokens);
|
---|
319 | } else if(kw == "mapcover") {
|
---|
320 | if(tokens.size()<1) {
|
---|
321 | cout<<"Usage: mapcover map v1,v2 [varname]"<<endl;
|
---|
322 | return(0);
|
---|
323 | }
|
---|
324 | MapCover(tokens);
|
---|
325 | } else if(kw == "map2cl") {
|
---|
326 | if(tokens.size()<2) {
|
---|
327 | cout<<"Usage: map2cl map clvec [lmax] [niter]"<<endl;
|
---|
328 | return(0);
|
---|
329 | }
|
---|
330 | Map2Cl(tokens);
|
---|
331 | } else if(kw == "cl2map") {
|
---|
332 | if(tokens.size()<2) {
|
---|
333 | cout<<"Usage: cl2map clvec map nside [fwhm]"<<endl;
|
---|
334 | return(0);
|
---|
335 | }
|
---|
336 | Cl2Map(tokens);
|
---|
337 | } else if(kw == "map2alm") {
|
---|
338 | if(tokens.size()<2) {
|
---|
339 | cout<<"Usage: map2alm map(double) almmat [lmax] [niter]"<<endl;
|
---|
340 | return(0);
|
---|
341 | }
|
---|
342 | Map2Alm(tokens);
|
---|
343 | } else if(kw == "alm2map") {
|
---|
344 | if(tokens.size()<3) {
|
---|
345 | cout<<"Usage: alm2map almmat map nside"<<endl;
|
---|
346 | return(0);
|
---|
347 | }
|
---|
348 | Alm2Map(tokens);
|
---|
349 | } else if(kw == "crmapmask") {
|
---|
350 | if(tokens.size()<2) {
|
---|
351 | cout<<"Usage: crmapmask mapmsk nside lat1,lat2 lon1,lon2 [valmsk,valnomsk]"<<endl;
|
---|
352 | return(0);
|
---|
353 | }
|
---|
354 | CrMapMask(tokens);
|
---|
355 | } else if(kw == "crmaskfrmap") {
|
---|
356 | if(tokens.size()<3) {
|
---|
357 | cout<<"Usage: crmaskfrmap mapmsk nside map(double) v1,v2 [valmsk,valnomsk]"<<endl;
|
---|
358 | return(0);
|
---|
359 | }
|
---|
360 | CrMaskFrMap(tokens);
|
---|
361 | } else if(kw == "maskmap") {
|
---|
362 | if(tokens.size()<2) {
|
---|
363 | cout<<"Usage: maskmap map msk"<<endl;
|
---|
364 | return(0);
|
---|
365 | }
|
---|
366 | MaskMap(tokens);
|
---|
367 | } else if(kw == "maproj") {
|
---|
368 | if(tokens.size()<2) {
|
---|
369 | cout<<"Usage: maproj map mapproj [nside]"<<endl;
|
---|
370 | return(0);
|
---|
371 | }
|
---|
372 | MapProj(tokens);
|
---|
373 | } else if(kw == "map2local") {
|
---|
374 | if(tokens.size()<5) {
|
---|
375 | cout<<"Usage: map2local map localmap nx,ny angleX,angleY phi0,theta0 [x0,y0] [angle]"<<endl;
|
---|
376 | return(0);
|
---|
377 | }
|
---|
378 | Map2Local(tokens);
|
---|
379 | } else if(kw == "mapop") {
|
---|
380 | if(tokens.size()<3) {
|
---|
381 | cout<<"Usage: mapop map op map1 [op map2] [op map3] [...]"<<endl;
|
---|
382 | return(0);
|
---|
383 | }
|
---|
384 | MapOper(tokens);
|
---|
385 | } else if(kw == "mapstat") {
|
---|
386 | if(tokens.size()<1) {
|
---|
387 | cout<<"Usage: Usage: mapstat map [msk] [meanvarname] [sigmavarname]"<<endl;
|
---|
388 | return(0);
|
---|
389 | }
|
---|
390 | MapStat(tokens);
|
---|
391 | } else if(kw == "alm2cl") {
|
---|
392 | if(tokens.size()<2) {
|
---|
393 | cout<<"Usage: alm2cl alm cl"<<endl;
|
---|
394 | return(0);
|
---|
395 | }
|
---|
396 | Alm2Cl(tokens);
|
---|
397 | } else if(kw == "cl2llcl") {
|
---|
398 | if(tokens.size()<2) {
|
---|
399 | cout<<"Usage: cl2llcl cl llcl"<<endl;
|
---|
400 | return(0);
|
---|
401 | }
|
---|
402 | Cl2llCl(tokens);
|
---|
403 | } else if(kw == "clmean") {
|
---|
404 | if(tokens.size()<1) {
|
---|
405 | cout<<"Usage: clmean cl imin,imax [meanvarname]"<<endl;
|
---|
406 | return(0);
|
---|
407 | }
|
---|
408 | ClMean(tokens);
|
---|
409 | } else if(kw == "clmult") {
|
---|
410 | if(tokens.size()<1) {
|
---|
411 | cout<<"Usage: clmult cl val"<<endl;
|
---|
412 | return(0);
|
---|
413 | }
|
---|
414 | ClMult(tokens);
|
---|
415 | } else if(kw == "clop") {
|
---|
416 | if(tokens.size()<3) {
|
---|
417 | cout<<"Usage: clop cl op1 cl1 [op2 cl2] [op2 cl3] [...]"<<endl;
|
---|
418 | return(0);
|
---|
419 | }
|
---|
420 | ClOper(tokens);
|
---|
421 | } else if(kw == "clrebin") {
|
---|
422 | if(tokens.size()<2) {
|
---|
423 | cout<<"Usage: clrebin cl ntuple nbin,bin0"<<endl;
|
---|
424 | return(0);
|
---|
425 | }
|
---|
426 | ClRebin(tokens);
|
---|
427 | }
|
---|
428 |
|
---|
429 | return(0);
|
---|
430 | }
|
---|
431 |
|
---|
432 | /* --Methode-- */
|
---|
433 | void skymapmoduleExecutor::Resol2SizeIndex(vector<string>& tokens)
|
---|
434 | {
|
---|
435 | double resamin = atof(tokens[0].c_str());
|
---|
436 | double resrad = Angle(resamin,Angle::ArcMin).ToRadian();
|
---|
437 | int_4 mh = SphereHEALPix<r_4>::ResolToSizeIndex(resrad);
|
---|
438 | int_4 mt = SphereThetaPhi<r_4>::ResolToSizeIndex(resrad);
|
---|
439 | cout << "Resol2SizeIndex: Resol= " << resamin << " arcmin ->"
|
---|
440 | << resrad << " radian NPix ~= " << 4.*M_PI/(resrad*resrad)
|
---|
441 | << " \n ===> nside=m_HEALPix= " << mh << " m_ThetaPhi= " << mt << endl;
|
---|
442 | return;
|
---|
443 | }
|
---|
444 |
|
---|
445 | /* --Methode-- */
|
---|
446 | void skymapmoduleExecutor::TypeMap(vector<string>& tokens)
|
---|
447 | {
|
---|
448 | NamedObjMgr omg;
|
---|
449 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
450 | uint_2 t = typemap(obj);
|
---|
451 | if(t==UnKnown) {cout<<"TypeMap: UnKnown"<<endl; return;}
|
---|
452 |
|
---|
453 | int nside = 0;
|
---|
454 | if(t&TypeR8) {
|
---|
455 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
456 | nside = map->SizeIndex();
|
---|
457 | } else if(t&TypeR4) {
|
---|
458 | SphericalMap<r_4>* map = dynamic_cast<SphericalMap<r_4>*>(obj);
|
---|
459 | nside = map->SizeIndex();
|
---|
460 | }
|
---|
461 |
|
---|
462 | string dum = "";
|
---|
463 | if(t==HealPix8) dum = "SphereHEALPix<r_8>";
|
---|
464 | else if(t==HealPix4) dum = "SphereHEALPix<r_4>";
|
---|
465 | else if(t==ThetaPhi8) dum = "SphereThetaPhi<r_8>";
|
---|
466 | else if(t==ThetaPhi4) dum = "SphereThetaPhi<r_4>";
|
---|
467 | else if(t==Spherical8) dum = "SphericalMap<r_8>";
|
---|
468 | else if(t==Spherical4) dum = "SphericalMap<r_4>";
|
---|
469 |
|
---|
470 | cout<<"TypeMap: "<<dum<<" nside="<<nside<<endl;
|
---|
471 | }
|
---|
472 |
|
---|
473 | /* --Methode-- */
|
---|
474 | void skymapmoduleExecutor::Map2Double(vector<string>& tokens)
|
---|
475 | {
|
---|
476 | NamedObjMgr omg;
|
---|
477 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
478 | if(obj==NULL) {
|
---|
479 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
480 | return;
|
---|
481 | }
|
---|
482 | SphericalMap<r_4>* map = dynamic_cast<SphericalMap<r_4>*>(obj);
|
---|
483 | if(map==NULL) {
|
---|
484 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_4>"<<endl;
|
---|
485 | return;
|
---|
486 | }
|
---|
487 | uint_2 t = typemap(obj);
|
---|
488 | int_4 nside = map->SizeIndex();
|
---|
489 |
|
---|
490 | SphericalMap<r_8>* map8 = NULL;
|
---|
491 | if(t&HealPix) map8 = new SphereHEALPix<r_8>(nside);
|
---|
492 | else if(t&ThetaPhi) map8 = new SphereThetaPhi<r_8>(nside);
|
---|
493 | else return;
|
---|
494 |
|
---|
495 | cout<<"Map2Double: converting to double "<<tokens[0]<<" nside="<<nside<<endl;
|
---|
496 | for(int_4 i=0;i<map8->NbPixels();i++) (*map8)(i) = (r_8) (*map)(i);
|
---|
497 |
|
---|
498 | omg.AddObj(map8,tokens[0]);
|
---|
499 | }
|
---|
500 |
|
---|
501 | /* --Methode-- */
|
---|
502 | void skymapmoduleExecutor::Map2Float(vector<string>& tokens)
|
---|
503 | {
|
---|
504 | NamedObjMgr omg;
|
---|
505 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
506 | if(obj==NULL) {
|
---|
507 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
508 | return;
|
---|
509 | }
|
---|
510 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
511 | if(map==NULL) {
|
---|
512 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
513 | return;
|
---|
514 | }
|
---|
515 | uint_2 t = typemap(obj);
|
---|
516 | int_4 nside = map->SizeIndex();
|
---|
517 |
|
---|
518 | SphericalMap<r_4>* map4 = NULL;
|
---|
519 | if(t&HealPix) map4 = new SphereHEALPix<r_4>(nside);
|
---|
520 | else if(t&ThetaPhi) map4 = new SphereThetaPhi<r_4>(nside);
|
---|
521 | else return;
|
---|
522 |
|
---|
523 | cout<<"Map2Float: converting to float "<<tokens[0]<<" nside="<<nside<<endl;
|
---|
524 | for(int_4 i=0;i<map4->NbPixels();i++) (*map4)(i) = (r_4) (*map)(i);
|
---|
525 |
|
---|
526 | omg.AddObj(map4,tokens[0]);
|
---|
527 | }
|
---|
528 |
|
---|
529 | /* --Methode-- */
|
---|
530 | void skymapmoduleExecutor::Map2Map(vector<string>& tokens)
|
---|
531 | {
|
---|
532 | NamedObjMgr omg;
|
---|
533 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
534 | if(obj==NULL) {
|
---|
535 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
536 | return;
|
---|
537 | }
|
---|
538 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
539 | if(map==NULL) {
|
---|
540 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
541 | return;
|
---|
542 | }
|
---|
543 | uint_2 t = typemap(obj);
|
---|
544 | int_4 nside = map->SizeIndex();
|
---|
545 |
|
---|
546 | uint_2 tomap = UnKnown;
|
---|
547 | if(toupper(tokens[1][0])=='H') tomap = HealPix;
|
---|
548 | else if(toupper(tokens[1][0])=='T') tomap = ThetaPhi;
|
---|
549 | else {cout<<"unkown map type "<<(char)toupper(tokens[1][0])<<" (H or T)!"<<endl; return;}
|
---|
550 | if(t&tomap) {cout<<"map of same type, no conversion done"<<endl; return;}
|
---|
551 |
|
---|
552 | int_4 tonside = nside;
|
---|
553 | if(tokens.size()>2) tonside = atoi(tokens[2].c_str());
|
---|
554 | else {
|
---|
555 | if(tomap&HealPix) {tonside=int(nside/1.5); GetNSideGood(tonside);}
|
---|
556 | else tonside=int(map->NbThetaSlices()/2.5);
|
---|
557 | }
|
---|
558 | if(tomap&HealPix && !IsNSideGood(tonside))
|
---|
559 | {cout<<"Bad nside for Healpix "<<tonside<<endl; return;}
|
---|
560 |
|
---|
561 | SphericalMap<r_8>* map8 = NULL;
|
---|
562 | if(tomap&HealPix) map8 = new SphereHEALPix<r_8>(tonside);
|
---|
563 | else if(tomap&ThetaPhi) map8 = new SphereThetaPhi<r_8>(tonside);
|
---|
564 | else return;
|
---|
565 |
|
---|
566 | cout<<"Map2Map: convert "<<tokens[0]<<" nside="<<nside
|
---|
567 | <<" to type "<<tokens[1]<<" tonside="<<tonside
|
---|
568 | <<" ntheta="<<map->NbThetaSlices()<<" -> "<<map8->NbThetaSlices()<<endl;
|
---|
569 |
|
---|
570 | for(int_4 i=0;i<map8->NbPixels();i++) {
|
---|
571 | double theta,phi; map8->PixThetaPhi(i,theta,phi);
|
---|
572 | int_4 ip = map->PixIndexSph(theta,phi);
|
---|
573 | if(ip<0 || ip>=map->NbPixels()) continue;
|
---|
574 | (*map8)(i)=(*map)(ip);
|
---|
575 | }
|
---|
576 |
|
---|
577 | omg.AddObj(map8,tokens[0]);
|
---|
578 | }
|
---|
579 |
|
---|
580 | /* --Methode-- */
|
---|
581 | void skymapmoduleExecutor::MapMult(vector<string>& tokens)
|
---|
582 | {
|
---|
583 | NamedObjMgr omg;
|
---|
584 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
585 | if(obj==NULL) {
|
---|
586 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
587 | return;
|
---|
588 | }
|
---|
589 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
590 | if(map==NULL) {
|
---|
591 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
592 | return;
|
---|
593 | }
|
---|
594 | double v = atof(tokens[1].c_str());
|
---|
595 | cout<<"MapMult: Multiply "<<tokens[0]<<" by "<<v<<endl;
|
---|
596 | for(int_4 i=0;i<map->NbPixels();i++) (*map)(i) *= v;
|
---|
597 | }
|
---|
598 |
|
---|
599 | /* --Methode-- */
|
---|
600 | void skymapmoduleExecutor::MapCover(vector<string>& tokens)
|
---|
601 | {
|
---|
602 | NamedObjMgr omg;
|
---|
603 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
604 | if(obj==NULL) {
|
---|
605 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
606 | return;
|
---|
607 | }
|
---|
608 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
609 | if(map==NULL) {
|
---|
610 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
611 | return;
|
---|
612 | }
|
---|
613 |
|
---|
614 | double v1=0.99, v2=1.01;
|
---|
615 | if(tokens.size()>1) if(tokens[1][0]!='!')
|
---|
616 | sscanf(tokens[1].c_str(),"%lf,%lf",&v1,&v2);
|
---|
617 |
|
---|
618 | cout<<"MapCover: "<<tokens[0]<<" good px=["<<v1<<","<<v2<<"]"<<endl;
|
---|
619 |
|
---|
620 | int_4 ngood=0;
|
---|
621 | for(int_4 i=0;i<map->NbPixels();i++)
|
---|
622 | if(v1<=(*map)(i) && (*map)(i)<=v2) ngood++;
|
---|
623 | double per = (double)ngood/(double)map->NbPixels();
|
---|
624 | cout<<"NGood="<<ngood<<" NTot="<<map->NbPixels()
|
---|
625 | <<" -> "<<100.*per<<" %";
|
---|
626 | if(tokens.size()>2) {
|
---|
627 | if(omg.HasVar(tokens[2])) omg.DeleteVar(tokens[2]);
|
---|
628 | char str[512]; sprintf(str,"%g",per);
|
---|
629 | omg.SetVar(tokens[2],(string)str);
|
---|
630 | cout<<" stored into $"<<tokens[2];
|
---|
631 | }
|
---|
632 | cout<<endl;
|
---|
633 | }
|
---|
634 |
|
---|
635 | /* --Methode-- */
|
---|
636 | void skymapmoduleExecutor::Map2Cl(vector<string>& tokens)
|
---|
637 | {
|
---|
638 | NamedObjMgr omg;
|
---|
639 |
|
---|
640 | int_4 lmax=0, niter=3;
|
---|
641 | if(tokens.size()>2) if(tokens[2][0]!='!')
|
---|
642 | lmax = atoi(tokens[2].c_str());
|
---|
643 | if(tokens.size()>3) niter = atoi(tokens[3].c_str());
|
---|
644 | if(niter<=0) niter = 3;
|
---|
645 |
|
---|
646 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
647 | if(obj==NULL) {
|
---|
648 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
649 | return;
|
---|
650 | }
|
---|
651 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
652 | if(map==NULL) {
|
---|
653 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
654 | return;
|
---|
655 | }
|
---|
656 | uint_2 t = typemap(obj);
|
---|
657 |
|
---|
658 | SphericalTransformServer<r_8> ylmserver;
|
---|
659 |
|
---|
660 | int_4 nside = map->SizeIndex();
|
---|
661 | if(lmax<=0) lmax = 3*nside-1;
|
---|
662 |
|
---|
663 | SphericalMap<r_8>* tmpmap = NULL;
|
---|
664 | if(t&HealPix) tmpmap = new SphereHEALPix<r_8>(*((SphereHEALPix<r_8>*)map),false);
|
---|
665 | else if(t&ThetaPhi) tmpmap = new SphereThetaPhi<r_8>(*((SphereThetaPhi<r_8>*)map),false);
|
---|
666 | else return;
|
---|
667 |
|
---|
668 | cout<<"Map2Cl: "<<tokens[0]<<"(nside="<<nside
|
---|
669 | <<") lmax="<<lmax<<" niter="<<niter<<" -> "<<tokens[1];
|
---|
670 |
|
---|
671 | TVector<r_8> cltmp = ylmserver.DecomposeToCl(*tmpmap,lmax,0.,niter);
|
---|
672 | cout<<"("<<cltmp.Size()<<")"<<endl;
|
---|
673 | delete tmpmap;
|
---|
674 |
|
---|
675 | omg.AddObj(cltmp,tokens[1]);
|
---|
676 | }
|
---|
677 |
|
---|
678 | /* --Methode-- */
|
---|
679 | void skymapmoduleExecutor::Cl2Map(vector<string>& tokens)
|
---|
680 | {
|
---|
681 | NamedObjMgr omg;
|
---|
682 |
|
---|
683 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
684 | if(obj==NULL) {
|
---|
685 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
686 | return;
|
---|
687 | }
|
---|
688 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
689 | if(cl==NULL) {
|
---|
690 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
691 | return;
|
---|
692 | }
|
---|
693 |
|
---|
694 | int nside = 128;
|
---|
695 | if(tokens.size()>2) nside = atoi(tokens[2].c_str());
|
---|
696 | if(DefTypMap&HealPix && !IsNSideGood(nside))
|
---|
697 | {cout<<"Cl2Map: Bad nside for HealPix "<<nside<<endl; return;}
|
---|
698 |
|
---|
699 | double fwhm = 0.;
|
---|
700 | if(tokens.size()>3) fwhm = atof(tokens[3].c_str());
|
---|
701 |
|
---|
702 | SphericalMap<r_8>* map = NULL;
|
---|
703 | if(DefTypMap&HealPix) map = new SphereHEALPix<r_8>;
|
---|
704 | else if(DefTypMap&ThetaPhi) map = new SphereThetaPhi<r_8>;
|
---|
705 | else return;
|
---|
706 |
|
---|
707 | SphericalTransformServer<r_8> ylmserver;
|
---|
708 | cout<<"Cl2Map: Generate map "<<tokens[1]<<" with nside="<<nside
|
---|
709 | <<" from cl "<<tokens[0]<<"("<<cl->Size()<<")"
|
---|
710 | <<" with fwhm="<<fwhm<<endl;
|
---|
711 | ylmserver.GenerateFromCl(*map,nside,*cl,fwhm);
|
---|
712 |
|
---|
713 | omg.AddObj(map,tokens[1]);
|
---|
714 | }
|
---|
715 |
|
---|
716 | /* --Methode-- */
|
---|
717 | void skymapmoduleExecutor::Map2Alm(vector<string>& tokens)
|
---|
718 | {
|
---|
719 | NamedObjMgr omg;
|
---|
720 |
|
---|
721 | int_4 lmax=0, niter=3;
|
---|
722 | if(tokens.size()>2) if(tokens[2][0]!='!')
|
---|
723 | lmax = atoi(tokens[2].c_str());
|
---|
724 | if(tokens.size()>3) niter = atoi(tokens[3].c_str());
|
---|
725 | if(niter<=0) niter = 3;
|
---|
726 |
|
---|
727 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
728 | if(obj==NULL) {
|
---|
729 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
730 | return;
|
---|
731 | }
|
---|
732 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
733 | if(map==NULL) {
|
---|
734 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
735 | return;
|
---|
736 | }
|
---|
737 | uint_2 t = typemap(obj);
|
---|
738 |
|
---|
739 | SphericalTransformServer<r_8> ylmserver;
|
---|
740 |
|
---|
741 | int_4 nside = map->SizeIndex();
|
---|
742 | if(lmax<=0) lmax = 3*nside-1;
|
---|
743 |
|
---|
744 | SphericalMap<r_8>* tmpmap = NULL;
|
---|
745 | if(t&HealPix) tmpmap = new SphereHEALPix<r_8>(*((SphereHEALPix<r_8>*)map),false);
|
---|
746 | else if(t&ThetaPhi) tmpmap = new SphereThetaPhi<r_8>(*((SphereThetaPhi<r_8>*)map),false);
|
---|
747 | else return;
|
---|
748 |
|
---|
749 | cout<<"Map2Alm: "<<tokens[0]<<"(nside="<<nside
|
---|
750 | <<") lmax="<<lmax<<" niter="<<niter<<" -> "<<tokens[1];
|
---|
751 |
|
---|
752 | Alm<r_8> almtmp;
|
---|
753 | ylmserver.DecomposeToAlm(*tmpmap,almtmp,lmax,0.,niter);
|
---|
754 | cout<<"("<<almtmp.rowNumber()<<")"<<endl;
|
---|
755 | delete tmpmap;
|
---|
756 |
|
---|
757 | int n = almtmp.rowNumber();
|
---|
758 | TMatrix< complex<r_8> >* alm = new TMatrix< complex<r_8> >(n,n);
|
---|
759 | *alm = (complex<r_8>)0.;
|
---|
760 | for(int i=0;i<n;i++) for(int j=0;j<=i;j++) (*alm)(i,j)=almtmp(i,j);
|
---|
761 | omg.AddObj(alm,tokens[1]);
|
---|
762 | }
|
---|
763 |
|
---|
764 | /* --Methode-- */
|
---|
765 | void skymapmoduleExecutor::Alm2Map(vector<string>& tokens)
|
---|
766 | {
|
---|
767 | NamedObjMgr omg;
|
---|
768 |
|
---|
769 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
770 | if(obj==NULL) {
|
---|
771 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
772 | return;
|
---|
773 | }
|
---|
774 | TMatrix< complex<r_8> >* almmat = dynamic_cast<TMatrix< complex<r_8> >*>(obj);
|
---|
775 | if(almmat==NULL) {
|
---|
776 | cout<<"Error: "<<tokens[0]<<" is not a TMatrix< complex<r_8> >"<<endl;
|
---|
777 | return;
|
---|
778 | }
|
---|
779 |
|
---|
780 | int lmax = almmat->NRows()-1;
|
---|
781 | Alm<r_8> alm(lmax);
|
---|
782 | alm = (complex<r_8>)0.;
|
---|
783 | for(int i=0;i<lmax+1;i++) for(int j=0;j<=i;j++) alm(i,j)=(*almmat)(i,j);
|
---|
784 |
|
---|
785 | int nside = 128;
|
---|
786 | if(tokens.size()>2) nside = atoi(tokens[2].c_str());
|
---|
787 | if(DefTypMap&HealPix && !IsNSideGood(nside))
|
---|
788 | {cout<<"Alm2Map: Bad nside for HealPix "<<nside<<endl; return;}
|
---|
789 |
|
---|
790 | SphericalMap<r_8>* map = NULL;
|
---|
791 | if(DefTypMap&HealPix) map = new SphereHEALPix<r_8>;
|
---|
792 | else if(DefTypMap&ThetaPhi) map = new SphereThetaPhi<r_8>;
|
---|
793 | else return;
|
---|
794 |
|
---|
795 | SphericalTransformServer<r_8> ylmserver;
|
---|
796 | cout<<"Alm2Map: Generate map "<<tokens[1]<<" with nside="<<nside
|
---|
797 | <<" from alm "<<tokens[0]<<"("<<alm.rowNumber()<<")"<<endl;
|
---|
798 | ylmserver.GenerateFromAlm(*map,nside,alm);
|
---|
799 |
|
---|
800 | omg.AddObj(map,tokens[1]);
|
---|
801 | }
|
---|
802 |
|
---|
803 | /* --Methode-- */
|
---|
804 | void skymapmoduleExecutor::CrMapMask(vector<string>& tokens)
|
---|
805 | {
|
---|
806 | NamedObjMgr omg;
|
---|
807 |
|
---|
808 | int_4 nside = atoi(tokens[1].c_str());
|
---|
809 | if(DefTypMap&HealPix && !IsNSideGood(nside))
|
---|
810 | {cout<<"CrMapMask: Bad nside "<<nside<<endl; return;}
|
---|
811 |
|
---|
812 | SphericalMap<r_8>* map = NULL;
|
---|
813 | if(DefTypMap&HealPix) map = new SphereHEALPix<r_8>(nside);
|
---|
814 | else if(DefTypMap&ThetaPhi) map = new SphereThetaPhi<r_8>(nside);
|
---|
815 | else return;
|
---|
816 |
|
---|
817 | bool lat=false, lon=false;
|
---|
818 | double lat1=-99., lat2=99., lon1=-999., lon2=999.;
|
---|
819 | if(tokens.size()>2) if(tokens[2][0]!='!')
|
---|
820 | {lat=true; sscanf(tokens[2].c_str(),"%lf,%lf",&lat1,&lat2);}
|
---|
821 | if(tokens.size()>3) if(tokens[3][0]!='!')
|
---|
822 | {lon=true; sscanf(tokens[3].c_str(),"%lf,%lf",&lon1,&lon2);}
|
---|
823 |
|
---|
824 | double valmsk=0., unvalmsk=1.;
|
---|
825 | if(tokens.size()>4) sscanf(tokens[4].c_str(),"%lf,%lf",&valmsk,&unvalmsk);
|
---|
826 |
|
---|
827 | cout<<"CrMapMask "<<tokens[0]<<" nside="<<nside;
|
---|
828 | if(lat) cout<<" for latitude in ["<<lat1<<","<<lat2<<"]";
|
---|
829 | if(lon) cout<<" for longitude in ["<<lon1<<","<<lon2<<"]";
|
---|
830 | cout<<" valmsk="<<valmsk<<" unvalmsk="<<unvalmsk<<endl;
|
---|
831 |
|
---|
832 | //Attention: conversion en co-latitude ==> inversion: [lat2,lat1]
|
---|
833 | lat1 = (90.-lat1)*M_PI/180.; lat2 = (90.-lat2)*M_PI/180.;
|
---|
834 | lon1 *= M_PI/180.; lon2 *= M_PI/180.;
|
---|
835 | for(int_4 i=0;i<map->NbPixels();i++) {
|
---|
836 | (*map)(i)=unvalmsk;
|
---|
837 | if(!lat && !lon) continue;
|
---|
838 | double theta,phi; map->PixThetaPhi(i,theta,phi);
|
---|
839 | if(lat && (theta<lat2 || lat1<theta)) continue;
|
---|
840 | if(lon && ( phi<lon1 || lon2<phi )) continue;
|
---|
841 | (*map)(i)=valmsk;
|
---|
842 | }
|
---|
843 |
|
---|
844 | omg.AddObj(map,tokens[0]);
|
---|
845 | }
|
---|
846 |
|
---|
847 | /* --Methode-- */
|
---|
848 | void skymapmoduleExecutor::CrMaskFrMap(vector<string>& tokens)
|
---|
849 | {
|
---|
850 | NamedObjMgr omg;
|
---|
851 |
|
---|
852 | AnyDataObj* obj=omg.GetObj(tokens[2]);
|
---|
853 | if(obj==NULL) {
|
---|
854 | cout<<"Error: Pas d'objet de nom "<<tokens[2]<<endl;
|
---|
855 | return;
|
---|
856 | }
|
---|
857 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
858 | if(map==NULL) {
|
---|
859 | cout<<"Error: "<<tokens[2]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
860 | return;
|
---|
861 | }
|
---|
862 | uint_2 t = typemap(obj);
|
---|
863 |
|
---|
864 | int_4 nside = atoi(tokens[1].c_str());
|
---|
865 | if(nside<=1) return;
|
---|
866 | if(t&HealPix && !IsNSideGood(nside))
|
---|
867 | {cout<<"CrMapMask: Bad nside for HealPix "<<nside<<endl; return;}
|
---|
868 |
|
---|
869 | SphericalMap<r_8>* msk = NULL;
|
---|
870 | if(t&HealPix) msk = new SphereHEALPix<r_8>(nside);
|
---|
871 | else if(t&ThetaPhi) msk = new SphereThetaPhi<r_8>(nside);
|
---|
872 | else return;
|
---|
873 |
|
---|
874 | double v1=-1.e100, v2=1.e100;
|
---|
875 | if(tokens.size()>3) if(tokens[3][0]!='!')
|
---|
876 | sscanf(tokens[3].c_str(),"%lf,%lf",&v1,&v2);
|
---|
877 |
|
---|
878 | double valmsk=0., unvalmsk=1.;
|
---|
879 | if(tokens.size()>4) sscanf(tokens[4].c_str(),"%lf,%lf",&valmsk,&unvalmsk);
|
---|
880 |
|
---|
881 | cout<<"CrMaskFrMap "<<tokens[0]<<" nside="<<nside<<" with "<<tokens[2]<<endl
|
---|
882 | <<" for value in ["<<v1<<","<<v2<<"]"
|
---|
883 | <<" valmsk="<<valmsk<<" unvalmsk="<<unvalmsk<<endl;
|
---|
884 |
|
---|
885 | for(int_4 i=0;i<msk->NbPixels();i++) {
|
---|
886 | double theta,phi; msk->PixThetaPhi(i,theta,phi);
|
---|
887 | int_4 ip = map->PixIndexSph(theta,phi);
|
---|
888 | if(ip<0 || ip>=map->NbPixels()) continue;
|
---|
889 | double v = (*map)(ip);
|
---|
890 | if(v>v1 && v<v2) (*msk)(i)=valmsk;
|
---|
891 | else (*msk)(i)=unvalmsk;
|
---|
892 | }
|
---|
893 |
|
---|
894 | omg.AddObj(msk,tokens[0]);
|
---|
895 | }
|
---|
896 |
|
---|
897 | /* --Methode-- */
|
---|
898 | void skymapmoduleExecutor::MaskMap(vector<string>& tokens)
|
---|
899 | {
|
---|
900 | NamedObjMgr omg;
|
---|
901 |
|
---|
902 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
903 | if(obj==NULL) {
|
---|
904 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
905 | return;
|
---|
906 | }
|
---|
907 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
908 | if(map==NULL) {
|
---|
909 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
910 | return;
|
---|
911 | }
|
---|
912 |
|
---|
913 | AnyDataObj* objm=omg.GetObj(tokens[1]);
|
---|
914 | if(obj==NULL) {
|
---|
915 | cout<<"Error: Pas d'objet de nom "<<tokens[1]<<endl;
|
---|
916 | return;
|
---|
917 | }
|
---|
918 | SphericalMap<r_8>* msk = dynamic_cast<SphericalMap<r_8>*>(objm);
|
---|
919 | if(msk==NULL) {
|
---|
920 | cout<<"Error: "<<tokens[1]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
921 | return;
|
---|
922 | }
|
---|
923 |
|
---|
924 | cout<<"MskMap: "<<tokens[0]<<" with mask "<<tokens[1]<<endl;
|
---|
925 |
|
---|
926 | for(int_4 i=0;i<map->NbPixels();i++) {
|
---|
927 | double theta,phi; map->PixThetaPhi(i,theta,phi);
|
---|
928 | int_4 ip = msk->PixIndexSph(theta,phi);
|
---|
929 | if(ip<0 || ip>=msk->NbPixels()) continue;
|
---|
930 | (*map)(i) *= (*msk)(ip);
|
---|
931 | }
|
---|
932 |
|
---|
933 | }
|
---|
934 |
|
---|
935 | /* --Methode-- */
|
---|
936 | void skymapmoduleExecutor::MapProj(vector<string>& tokens)
|
---|
937 | {
|
---|
938 | NamedObjMgr omg;
|
---|
939 |
|
---|
940 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
941 | if(obj==NULL) {
|
---|
942 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
943 | return;
|
---|
944 | }
|
---|
945 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
946 | if(map==NULL) {
|
---|
947 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
948 | return;
|
---|
949 | }
|
---|
950 | int_4 nsidemap = map->SizeIndex();
|
---|
951 | uint_2 t = typemap(obj);
|
---|
952 |
|
---|
953 | int_4 nside = nsidemap;
|
---|
954 | if(tokens.size()>2) nside = atoi(tokens[2].c_str());
|
---|
955 | if(t&HealPix && !IsNSideGood(nside))
|
---|
956 | {cout<<"MapProj: Bad nside for Healpix "<<nside<<endl; return;}
|
---|
957 |
|
---|
958 | SphericalMap<r_8>* mapproj = NULL;
|
---|
959 | if(t&HealPix) mapproj = new SphereHEALPix<r_8>(nside);
|
---|
960 | else if(t&ThetaPhi) mapproj = new SphereThetaPhi<r_8>(nside);
|
---|
961 | else return;
|
---|
962 |
|
---|
963 | cout<<"MapProj: project "<<tokens[0]<<" n="<<nsidemap
|
---|
964 | <<" to "<<tokens[1]<<" n="<<nside<<endl;
|
---|
965 |
|
---|
966 | if(nsidemap==nside) { // tailles egales
|
---|
967 | for(int_4 i=0;i<mapproj->NbPixels();i++)
|
---|
968 | (*mapproj)(i) = (*map)(i);
|
---|
969 | } else if(nsidemap<nside) { // mapproj>map
|
---|
970 | for(int_4 i=0;i<mapproj->NbPixels();i++) {
|
---|
971 | double theta,phi; mapproj->PixThetaPhi(i,theta,phi);
|
---|
972 | int_4 ip = map->PixIndexSph(theta,phi);
|
---|
973 | if(ip<0 || ip>=map->NbPixels()) continue;
|
---|
974 | (*mapproj)(i) = (*map)(ip);
|
---|
975 | }
|
---|
976 | } else { // mapproj<map
|
---|
977 | SphericalMap<int_4>* nproj= NULL;
|
---|
978 | if(t&HealPix) nproj= new SphereHEALPix<int_4>(nside);
|
---|
979 | else if(t&ThetaPhi) nproj= new SphereThetaPhi<int_4>(nside);
|
---|
980 | int_4 i;
|
---|
981 | for(i=0;i<nproj->NbPixels();i++) {(*nproj)(i)=0;(*mapproj)(i)=0.;}
|
---|
982 | for(i=0;i<map->NbPixels();i++) {
|
---|
983 | double theta,phi; map->PixThetaPhi(i,theta,phi);
|
---|
984 | int_4 ip = mapproj->PixIndexSph(theta,phi);
|
---|
985 | if(ip<0 || ip>=mapproj->NbPixels()) continue;
|
---|
986 | (*mapproj)(ip) += (*map)(i);
|
---|
987 | (*nproj)(ip)++;
|
---|
988 | }
|
---|
989 | for(i=0;i<mapproj->NbPixels();i++)
|
---|
990 | (*mapproj)(i) /= (r_8) (*nproj)(i);
|
---|
991 | delete nproj;
|
---|
992 | }
|
---|
993 |
|
---|
994 | omg.AddObj(mapproj,tokens[1]);
|
---|
995 | }
|
---|
996 |
|
---|
997 | /* --Methode-- */
|
---|
998 | void skymapmoduleExecutor::Map2Local(vector<string>& tokens)
|
---|
999 | {
|
---|
1000 | NamedObjMgr omg;
|
---|
1001 |
|
---|
1002 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1003 | if(obj==NULL) {
|
---|
1004 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1005 | return;
|
---|
1006 | }
|
---|
1007 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
1008 | if(map==NULL) {
|
---|
1009 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
1010 | return;
|
---|
1011 | }
|
---|
1012 |
|
---|
1013 | int_4 nx=10,ny=10;
|
---|
1014 | if(tokens.size()>2) sscanf(tokens[2].c_str(),"%d,%d",&nx,&ny);
|
---|
1015 | if(nx<=0)
|
---|
1016 | {cout<<"Error: nx="<<nx<<" value >0"<<endl; return;}
|
---|
1017 | if(ny<=0)
|
---|
1018 | {cout<<"Error: ny="<<ny<<" value >0"<<endl; return;}
|
---|
1019 |
|
---|
1020 | double angleX=1., angleY=1.;
|
---|
1021 | if(tokens.size()>3) sscanf(tokens[3].c_str(),"%lf,%lf",&angleX,&angleY);
|
---|
1022 | if(angleX<=0. || angleX>180.)
|
---|
1023 | {cout<<"Error: bad angleX="<<angleX<<" value #]0,180]"<<endl; return;}
|
---|
1024 | if(angleY<=0. || angleY>180.)
|
---|
1025 | {cout<<"Error: bad angleY="<<angleY<<" value #]0,180]"<<endl; return;}
|
---|
1026 |
|
---|
1027 | double phi0=0., theta0=0.;
|
---|
1028 | if(tokens.size()>4) sscanf(tokens[4].c_str(),"%lf,%lf",&phi0,&theta0);
|
---|
1029 | if(phi0<0. || phi0>=360.)
|
---|
1030 | {cout<<"Error: bad phi0="<<phi0<<" value #[0,360["<<endl; return;}
|
---|
1031 | if(theta0<-90. || theta0>90.)
|
---|
1032 | {cout<<"Error: bad theta0="<<theta0<<" value #[-90,90]"<<endl; return;}
|
---|
1033 |
|
---|
1034 | int_4 x0=nx/2, y0=ny/2;
|
---|
1035 | if(tokens.size()>5) if(tokens[5][0]!='!')
|
---|
1036 | sscanf(tokens[5].c_str(),"%d,%d",&x0,&y0);
|
---|
1037 |
|
---|
1038 | double angle=0.;
|
---|
1039 | if(tokens.size()>6) sscanf(tokens[6].c_str(),"%lf",&angle);
|
---|
1040 |
|
---|
1041 | cout<<"Map2Local: proj "<<tokens[0]<<" to local "<<tokens[1]<<endl
|
---|
1042 | <<" nx="<<nx<<" ny="<<ny<<" angleX="<<angleX<<" angleY="<<angleY
|
---|
1043 | <<" phi0="<<phi0<<" theta0="<<theta0<<" x0="<<x0<<" y0="<<y0
|
---|
1044 | <<" angle="<<angle<<endl;
|
---|
1045 | if(angle<-90. || angle>90.)
|
---|
1046 | {cout<<"Error: bad angle="<<angle<<" value #[-90,90]"<<endl; return;}
|
---|
1047 |
|
---|
1048 | theta0 = (90.-theta0);
|
---|
1049 | LocalMap<r_8>* lmap = new LocalMap<r_8>(nx,ny,angleX,angleY,theta0,phi0,x0,y0,angle);
|
---|
1050 | *lmap = 0.; //lmap->print(cout);
|
---|
1051 |
|
---|
1052 | int_4 nbad=0;
|
---|
1053 | double tmax=-9999., pmax=-9999., tmin=9999., pmin=9999.;
|
---|
1054 | for(int_4 i=0;i<lmap->NbPixels();i++) {
|
---|
1055 | double theta,phi;
|
---|
1056 | lmap->PixThetaPhi(i,theta,phi);
|
---|
1057 | int_4 ip = map->PixIndexSph(theta,phi);
|
---|
1058 | if(ip<0 || ip>=map->NbPixels()) {nbad++; continue;}
|
---|
1059 | if(theta<tmin) tmin=theta; if(theta>tmax) tmax=theta;
|
---|
1060 | if(phi<pmin) pmin=phi; if(phi>pmax) pmax=phi;
|
---|
1061 | (*lmap)(i) = (*map)(ip);
|
---|
1062 | }
|
---|
1063 | if(nbad) cout<<"Warning: "<<nbad<<" bad pixels numbers"<<endl;
|
---|
1064 | cout<<" phi=["<<pmin/M_PI*180.<<","<<pmax/M_PI*180.<<"]="
|
---|
1065 | <<(pmax-pmin)/M_PI*180.
|
---|
1066 | <<" theta=["<<tmin/M_PI*180.<<","<<tmax/M_PI*180.<<"]="
|
---|
1067 | <<(tmax-tmin)/M_PI*180.<<endl;
|
---|
1068 |
|
---|
1069 | omg.AddObj(lmap,tokens[1]);
|
---|
1070 | }
|
---|
1071 |
|
---|
1072 | /* --Methode-- */
|
---|
1073 | void skymapmoduleExecutor::MapOper(vector<string>& tokens)
|
---|
1074 | {
|
---|
1075 | NamedObjMgr omg;
|
---|
1076 |
|
---|
1077 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1078 | if(obj==NULL) {
|
---|
1079 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1080 | return;
|
---|
1081 | }
|
---|
1082 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
1083 | if(map==NULL) {
|
---|
1084 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
1085 | return;
|
---|
1086 | }
|
---|
1087 |
|
---|
1088 | cout<<"MapOper: "<<tokens[0];
|
---|
1089 | for(unsigned short iop=1;iop<tokens.size()-1;iop+=2) {
|
---|
1090 | string op = tokens[iop];
|
---|
1091 | if(op!="+" && op!="-" && op!="*" && op!="/") continue;
|
---|
1092 | AnyDataObj* obj1=omg.GetObj(tokens[iop+1]);
|
---|
1093 | if(obj1==NULL)
|
---|
1094 | {cout<<"Error: Pas d'objet de nom "<<tokens[iop+1]<<endl;
|
---|
1095 | continue;}
|
---|
1096 | SphericalMap<r_8>* map1 = dynamic_cast<SphericalMap<r_8>*>(obj1);
|
---|
1097 | if(map1==NULL)
|
---|
1098 | {cout<<"Error: "<<tokens[iop+1]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
1099 | continue;}
|
---|
1100 | cout<<" "<<op<<"= "<<tokens[iop+1];
|
---|
1101 | for(int_4 i=0;i<map->NbPixels();i++) {
|
---|
1102 | double theta,phi; map->PixThetaPhi(i,theta,phi);
|
---|
1103 | int_4 ip = map1->PixIndexSph(theta,phi);
|
---|
1104 | if(ip<0 || ip>=map1->NbPixels()) continue;
|
---|
1105 | if(op=="+") (*map)(i) += (*map1)(ip);
|
---|
1106 | else if(op=="-") (*map)(i) -= (*map1)(ip);
|
---|
1107 | else if(op=="*") (*map)(i) *= (*map1)(ip);
|
---|
1108 | else if(op=="/")
|
---|
1109 | {if((*map1)(ip)!=0.) (*map)(i)/=(*map1)(ip); else (*map)(i)=0.;}
|
---|
1110 | }
|
---|
1111 | }
|
---|
1112 | cout<<endl;
|
---|
1113 | }
|
---|
1114 |
|
---|
1115 | /* --Methode-- */
|
---|
1116 | void skymapmoduleExecutor::MapStat(vector<string>& tokens)
|
---|
1117 | {
|
---|
1118 | NamedObjMgr omg;
|
---|
1119 |
|
---|
1120 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1121 | if(obj==NULL) {
|
---|
1122 | cout<<"Error: Pas d'objet de nom "<<tokens[1]<<endl;
|
---|
1123 | return;
|
---|
1124 | }
|
---|
1125 | SphericalMap<r_8>* map = dynamic_cast<SphericalMap<r_8>*>(obj);
|
---|
1126 | if(map==NULL) {
|
---|
1127 | cout<<"Error: "<<tokens[0]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
1128 | return;
|
---|
1129 | }
|
---|
1130 |
|
---|
1131 | SphericalMap<r_8>* msk = NULL;
|
---|
1132 | if(tokens.size()>1) {
|
---|
1133 | if(tokens[1][0]!='!') {
|
---|
1134 | AnyDataObj* objm=omg.GetObj(tokens[1]);
|
---|
1135 | if(objm==NULL) {
|
---|
1136 | cout<<"Error: Pas d'objet de nom "<<tokens[1]<<endl;
|
---|
1137 | return;
|
---|
1138 | }
|
---|
1139 | msk = dynamic_cast<SphericalMap<r_8>*>(objm);
|
---|
1140 | if(msk==NULL) {
|
---|
1141 | cout<<"Error: "<<tokens[1]<<" is not a SphericalMap<r_8>"<<endl;
|
---|
1142 | return;
|
---|
1143 | }
|
---|
1144 | }
|
---|
1145 | }
|
---|
1146 |
|
---|
1147 | cout<<"MapStat for map "<<tokens[0];
|
---|
1148 | if(msk) cout<<" using mask "<<tokens[1];
|
---|
1149 | cout<<endl;
|
---|
1150 |
|
---|
1151 | double sum=0., sum2=0., sumw=0.;
|
---|
1152 | int npix = map->NbPixels(), npixgood=0;
|
---|
1153 | for(int_4 i=0;i<npix;i++) {
|
---|
1154 | double w = 1.;
|
---|
1155 | if(msk) {
|
---|
1156 | double theta,phi; map->PixThetaPhi(i,theta,phi);
|
---|
1157 | int_4 ip = msk->PixIndexSph(theta,phi);
|
---|
1158 | if(ip<0 || ip>=msk->NbPixels()) w=0.;
|
---|
1159 | else w=(*msk)(ip);
|
---|
1160 | }
|
---|
1161 | if(w<=0.) continue;
|
---|
1162 | npixgood++; sumw += w;
|
---|
1163 | sum += (*map)(i)*w;
|
---|
1164 | sum2 += (*map)(i)*(*map)(i)*w*w;
|
---|
1165 | }
|
---|
1166 | if(sumw<=0. || npixgood==0) {
|
---|
1167 | sum = sum2 = sumw=0.;
|
---|
1168 | } else {
|
---|
1169 | sum /= sumw;
|
---|
1170 | sum2 = sum2/sumw - sum*sum;
|
---|
1171 | if(sum2<=0.) sum2=0.; else sum2=sqrt(sum2);
|
---|
1172 | }
|
---|
1173 | cout<<"Number of good px="<<npixgood<<" / "<<npix
|
---|
1174 | <<" ("<<100.*npixgood/npix<<" %)"<<endl;
|
---|
1175 | cout<<" mean="<<sum<<" sigma="<<sum2<<" sigma^2="<<sum2*sum2<<endl;
|
---|
1176 |
|
---|
1177 | if(tokens.size()>2) if(tokens[2][0]!='!') {
|
---|
1178 | if(omg.HasVar(tokens[2])) omg.DeleteVar(tokens[2]);
|
---|
1179 | char str[512]; sprintf(str,"%.8f",sum);
|
---|
1180 | omg.SetVar(tokens[2],(string)str);
|
---|
1181 | cout<<" mean stored into $"<<tokens[2];
|
---|
1182 | }
|
---|
1183 | if(tokens.size()>3) if(tokens[3][0]!='!') {
|
---|
1184 | if(omg.HasVar(tokens[3])) omg.DeleteVar(tokens[3]);
|
---|
1185 | char str[512]; sprintf(str,"%.8f",sum2);
|
---|
1186 | omg.SetVar(tokens[3],(string)str);
|
---|
1187 | cout<<" sigma stored into $"<<tokens[3];
|
---|
1188 | }
|
---|
1189 | if(tokens.size()>2) cout<<endl;
|
---|
1190 |
|
---|
1191 | }
|
---|
1192 |
|
---|
1193 | /* --Methode-- */
|
---|
1194 | void skymapmoduleExecutor::Alm2Cl(vector<string>& tokens)
|
---|
1195 | {
|
---|
1196 | NamedObjMgr omg;
|
---|
1197 |
|
---|
1198 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1199 | if(obj==NULL) {
|
---|
1200 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1201 | return;
|
---|
1202 | }
|
---|
1203 | TMatrix< complex<r_8> >* almmat = dynamic_cast<TMatrix< complex<r_8> >*>(obj);
|
---|
1204 | if(almmat==NULL) {
|
---|
1205 | cout<<"Error: "<<tokens[0]<<" is not a TMatrix< complex<r_8> >"<<endl;
|
---|
1206 | return;
|
---|
1207 | }
|
---|
1208 |
|
---|
1209 | int lmax = almmat->NRows()-1;
|
---|
1210 | Alm<r_8> alm(lmax);
|
---|
1211 | alm = (complex<r_8>)0.;
|
---|
1212 | for(int i=0;i<lmax+1;i++) for(int j=0;j<=i;j++) alm(i,j)=(*almmat)(i,j);
|
---|
1213 |
|
---|
1214 | TVector<r_8> cl(almmat->NRows());
|
---|
1215 | cl = 0.;
|
---|
1216 |
|
---|
1217 | cout<<"Alm2Cl: project alm "<<tokens[0]<<"("<<alm.rowNumber()
|
---|
1218 | <<") to cl "<<tokens[1]<<"("<<cl.Size()<<")"<<endl;
|
---|
1219 |
|
---|
1220 | for(int_4 l=0;l<alm.rowNumber();l++) {
|
---|
1221 | cl(l) = alm(l,0).real()*alm(l,0).real()+alm(l,0).imag()*alm(l,0).imag();
|
---|
1222 | if(l>0) for(int_4 m=1;m<=l;m++)
|
---|
1223 | cl(l) += 2.*(alm(l,m).real()*alm(l,m).real()+alm(l,m).imag()*alm(l,m).imag());
|
---|
1224 | cl(l) /= 2.*l+1;
|
---|
1225 | }
|
---|
1226 |
|
---|
1227 | omg.AddObj(cl,tokens[1]);
|
---|
1228 | }
|
---|
1229 |
|
---|
1230 | /* --Methode-- */
|
---|
1231 | void skymapmoduleExecutor::Cl2llCl(vector<string>& tokens)
|
---|
1232 | {
|
---|
1233 | NamedObjMgr omg;
|
---|
1234 |
|
---|
1235 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1236 | if(obj==NULL) {
|
---|
1237 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1238 | return;
|
---|
1239 | }
|
---|
1240 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
1241 | if(cl==NULL) {
|
---|
1242 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
1243 | return;
|
---|
1244 | }
|
---|
1245 |
|
---|
1246 | TVector<r_8> llcl(*cl,false);
|
---|
1247 |
|
---|
1248 | cout<<"Cl2llCl: project "<<tokens[0]<<"("<<cl->Size()
|
---|
1249 | <<") to l*(l+1)*cl "<<tokens[1]<<"("<<llcl.Size()<<")"<<endl;
|
---|
1250 |
|
---|
1251 | for(int_4 l=0;l<llcl.Size();l++) llcl(l) *= (double)l*(l+1.);
|
---|
1252 |
|
---|
1253 | omg.AddObj(llcl,tokens[1]);
|
---|
1254 | }
|
---|
1255 |
|
---|
1256 | /* --Methode-- */
|
---|
1257 | void skymapmoduleExecutor::ClMean(vector<string>& tokens)
|
---|
1258 | {
|
---|
1259 | NamedObjMgr omg;
|
---|
1260 |
|
---|
1261 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1262 | if(obj==NULL) {
|
---|
1263 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1264 | return;
|
---|
1265 | }
|
---|
1266 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
1267 | if(cl==NULL) {
|
---|
1268 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
1269 | return;
|
---|
1270 | }
|
---|
1271 |
|
---|
1272 | int lmin=0, lmax=-1;
|
---|
1273 | if(tokens.size()>1) if(tokens[1][0]!='!')
|
---|
1274 | sscanf(tokens[1].c_str(),"%d,%d",&lmin,&lmax);
|
---|
1275 | if(lmin<0) lmin=0; if(lmax>=cl->Size()) lmax=cl->Size()-1;
|
---|
1276 | if(lmin>lmax) {lmin=0; lmax=cl->Size()-1;}
|
---|
1277 |
|
---|
1278 | cout<<"ClMean: compute mean "<<tokens[0]<<"("<<cl->Size()
|
---|
1279 | <<") between ["<<lmin<<","<<lmax<<"]";
|
---|
1280 |
|
---|
1281 | double mean = 0.;
|
---|
1282 | for(int_4 l=lmin;l<=lmax;l++) mean += (*cl)(l);
|
---|
1283 | mean /= lmax-lmin+1;
|
---|
1284 |
|
---|
1285 | cout<<" mean="<<mean;
|
---|
1286 |
|
---|
1287 | if(tokens.size()>2) {
|
---|
1288 | if(omg.HasVar(tokens[2])) omg.DeleteVar(tokens[2]);
|
---|
1289 | char str[512]; sprintf(str,"%.8f",mean);
|
---|
1290 | omg.SetVar(tokens[2],(string)str);
|
---|
1291 | cout<<" stored into $"<<tokens[2];
|
---|
1292 | }
|
---|
1293 | cout<<endl;
|
---|
1294 |
|
---|
1295 | }
|
---|
1296 |
|
---|
1297 | /* --Methode-- */
|
---|
1298 | void skymapmoduleExecutor::ClMult(vector<string>& tokens)
|
---|
1299 | {
|
---|
1300 | NamedObjMgr omg;
|
---|
1301 |
|
---|
1302 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1303 | if(obj==NULL) {
|
---|
1304 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1305 | return;
|
---|
1306 | }
|
---|
1307 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
1308 | if(cl==NULL) {
|
---|
1309 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
1310 | return;
|
---|
1311 | }
|
---|
1312 |
|
---|
1313 | double val = 1.;
|
---|
1314 | if(tokens.size()>1) val=atof(tokens[1].c_str());
|
---|
1315 |
|
---|
1316 | cout<<"ClMult: multiply "<<tokens[0]<<"("<<cl->Size()<<") by "<<val<<endl;
|
---|
1317 |
|
---|
1318 | *cl *= val;
|
---|
1319 |
|
---|
1320 | }
|
---|
1321 |
|
---|
1322 | /* --Methode-- */
|
---|
1323 | void skymapmoduleExecutor::ClOper(vector<string>& tokens)
|
---|
1324 | {
|
---|
1325 | NamedObjMgr omg;
|
---|
1326 |
|
---|
1327 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1328 | if(obj==NULL) {
|
---|
1329 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1330 | return;
|
---|
1331 | }
|
---|
1332 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
1333 | if(cl==NULL) {
|
---|
1334 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
1335 | return;
|
---|
1336 | }
|
---|
1337 |
|
---|
1338 | cout<<"ClOper: "<<tokens[0]<<"("<<cl->Size()<<")";
|
---|
1339 | for(unsigned short iop=1;iop<tokens.size()-1;iop+=2) {
|
---|
1340 | string op = tokens[iop];
|
---|
1341 | if(op!="+" && op!="-" && op!="*" && op!="/") continue;
|
---|
1342 | AnyDataObj* obj1=omg.GetObj(tokens[iop+1]);
|
---|
1343 | if(obj1==NULL)
|
---|
1344 | {cout<<"Error: Pas d'objet de nom "<<tokens[iop+1]<<endl;
|
---|
1345 | continue;}
|
---|
1346 | TVector<r_8>* cl1 = dynamic_cast<TVector<r_8>*>(obj1);
|
---|
1347 | if(cl1==NULL)
|
---|
1348 | {cout<<"Error: "<<tokens[iop+1]<<" is not a TVector<r_8>"<<endl;
|
---|
1349 | continue;}
|
---|
1350 | cout<<" "<<op<<"= "<<tokens[iop+1]<<"("<<cl1->Size()<<")";
|
---|
1351 | int_4 n=cl->Size(); if(n>cl1->Size()) n=cl1->Size();
|
---|
1352 | for(int_4 i=0;i<n;i++) {
|
---|
1353 | if(op=="+") (*cl)(i) += (*cl1)(i);
|
---|
1354 | else if(op=="-") (*cl)(i) -= (*cl1)(i);
|
---|
1355 | else if(op=="*") (*cl)(i) *= (*cl1)(i);
|
---|
1356 | else if(op=="/")
|
---|
1357 | {if((*cl1)(i)!=0.) (*cl)(i)/=(*cl1)(i); else (*cl)(i)=0.;}
|
---|
1358 | }
|
---|
1359 | }
|
---|
1360 | cout<<endl;
|
---|
1361 | }
|
---|
1362 | /* --Methode-- */
|
---|
1363 | void skymapmoduleExecutor::ClRebin(vector<string>& tokens)
|
---|
1364 | {
|
---|
1365 | NamedObjMgr omg;
|
---|
1366 |
|
---|
1367 | AnyDataObj* obj=omg.GetObj(tokens[0]);
|
---|
1368 | if(obj==NULL) {
|
---|
1369 | cout<<"Error: Pas d'objet de nom "<<tokens[0]<<endl;
|
---|
1370 | return;
|
---|
1371 | }
|
---|
1372 | TVector<r_8>* cl = dynamic_cast<TVector<r_8>*>(obj);
|
---|
1373 | if(cl==NULL) {
|
---|
1374 | cout<<"Error: "<<tokens[0]<<" is not a TVector<r_8>"<<endl;
|
---|
1375 | return;
|
---|
1376 | }
|
---|
1377 | int_4 nn = cl->Size();
|
---|
1378 | if(nn<=0) return;
|
---|
1379 |
|
---|
1380 | int_4 nrebin=1, ibin0=0;
|
---|
1381 | if(tokens.size()>2) sscanf(tokens[2].c_str(),"%d,%d",&nrebin,&ibin0);
|
---|
1382 | if(ibin0<0 || ibin0>=nn) ibin0=0;
|
---|
1383 |
|
---|
1384 | cout<<"ClRebin: rebin "<<tokens[0]<<"("<<nn<<") by nbin="<<nrebin
|
---|
1385 | <<" begin at "<<ibin0<<" -> "<<tokens[1]<<endl;
|
---|
1386 |
|
---|
1387 | const int ndim = 8;
|
---|
1388 | char *nament[ndim] =
|
---|
1389 | {"l","n","clmean","cllin","clpar","sclmean","scllin","sclpar"};
|
---|
1390 | NTuple clbin(ndim,nament);
|
---|
1391 | r_4 xnt[ndim];
|
---|
1392 |
|
---|
1393 | SLinParBuff slp(nrebin+2);
|
---|
1394 |
|
---|
1395 | for(int_4 i=ibin0;i<nn;i+=nrebin) {
|
---|
1396 | r_8 ll=0.;
|
---|
1397 | slp.Reset();
|
---|
1398 | int_4 j;
|
---|
1399 | for(j=0;j<ndim;j++) xnt[j]=0.;
|
---|
1400 | for(j=i;j<nn;j++) {
|
---|
1401 | ll += (r_8) j;
|
---|
1402 | slp.Push((r_8)j,(*cl)(j));
|
---|
1403 | if((int_4)slp.NPoints()>=nrebin) break;
|
---|
1404 | }
|
---|
1405 | if(slp.NPoints()<=0) continue;
|
---|
1406 | ll /= (r_8) slp.NPoints();
|
---|
1407 | xnt[0] = ll;
|
---|
1408 | xnt[1] = slp.NPoints();
|
---|
1409 | r_8 s,a0,a1,a2;
|
---|
1410 | s = slp.Compute(a0);
|
---|
1411 | if(s>0.) {xnt[2]=a0; xnt[5]=s;}
|
---|
1412 | s = slp.Compute(a0,a1);
|
---|
1413 | if(s>0.) {xnt[3]=a0+a1*ll; xnt[6]=s;}
|
---|
1414 | s = slp.Compute(a0,a1,a2);
|
---|
1415 | if(s>0.) {xnt[4]=a0+a1*ll+a2*ll*ll; xnt[7]=s;}
|
---|
1416 | clbin.Fill(xnt);
|
---|
1417 | }
|
---|
1418 |
|
---|
1419 | omg.AddObj(clbin,tokens[1]);
|
---|
1420 | }
|
---|
1421 |
|
---|
1422 | ////////////////////////////////////////////////////////////
|
---|
1423 |
|
---|
1424 | /* --Methode-- */
|
---|
1425 | void skymapmoduleExecutor::SetTypeMap(vector<string>& tokens)
|
---|
1426 | {
|
---|
1427 | if(tokens.size()<1) {
|
---|
1428 | cout<<"SetTypeMap: Usage: settypemap type"<<endl;
|
---|
1429 | } else {
|
---|
1430 | if(toupper(tokens[0][0])=='H') DefTypMap = HealPix;
|
---|
1431 | else if(toupper(tokens[0][0])=='T') DefTypMap = ThetaPhi;
|
---|
1432 | else cout<<"unkown map type, should be (H or T)!"<<endl;
|
---|
1433 | }
|
---|
1434 | string dum;
|
---|
1435 | if(DefTypMap&HealPix) dum="HealPix";
|
---|
1436 | else if(DefTypMap&ThetaPhi) dum="ThetaPhi";
|
---|
1437 | cout<<"SetTypeMap: type map is "<<dum<<" ("<<DefTypMap<<")"<<endl;
|
---|
1438 | }
|
---|
1439 |
|
---|
1440 | /* --Methode-- */
|
---|
1441 | bool skymapmoduleExecutor::IsNSideGood(int_4 nside)
|
---|
1442 | {
|
---|
1443 | if(nside<=1) return false;
|
---|
1444 | while(nside>1) {if(nside%2!=0) return false; else nside/=2;}
|
---|
1445 | return true;
|
---|
1446 | }
|
---|
1447 |
|
---|
1448 | /* --Methode-- */
|
---|
1449 | void skymapmoduleExecutor::GetNSideGood(int_4 &nside)
|
---|
1450 | // get the nearest good nside for HealPix
|
---|
1451 | {
|
---|
1452 | double n=1.e50; int_4 ns=nside;
|
---|
1453 | for(int_4 i=1;i<66000;i*=2) {
|
---|
1454 | double v = fabs((double)(nside-i));
|
---|
1455 | if(v>n) continue;
|
---|
1456 | n=v; ns=i;
|
---|
1457 | }
|
---|
1458 | nside = ns;
|
---|
1459 | }
|
---|
1460 |
|
---|
1461 | /* --Methode-- */
|
---|
1462 | uint_2 skymapmoduleExecutor::typemap(AnyDataObj* obj)
|
---|
1463 | {
|
---|
1464 | if(obj==NULL) return UnKnown;
|
---|
1465 | if(dynamic_cast<SphereHEALPix<r_4> *>(obj)) return HealPix4;
|
---|
1466 | if(dynamic_cast<SphereHEALPix<r_8> *>(obj)) return HealPix8;
|
---|
1467 | if(dynamic_cast<SphereThetaPhi<r_4>*>(obj)) return ThetaPhi4;
|
---|
1468 | if(dynamic_cast<SphereThetaPhi<r_8>*>(obj)) return ThetaPhi8;
|
---|
1469 | if(dynamic_cast<SphericalMap<r_4> *>(obj)) return Spherical4;
|
---|
1470 | if(dynamic_cast<SphericalMap<r_8> *>(obj)) return Spherical8;
|
---|
1471 | return UnKnown;
|
---|
1472 | }
|
---|
1473 |
|
---|
1474 | ////////////////////////////////////////////////////////////
|
---|
1475 | static skymapmoduleExecutor * piaskmex = NULL;
|
---|
1476 |
|
---|
1477 | void skymapmodule_init()
|
---|
1478 | {
|
---|
1479 | if (piaskmex) delete piaskmex;
|
---|
1480 | piaskmex = new skymapmoduleExecutor;
|
---|
1481 | }
|
---|
1482 |
|
---|
1483 | void skymapmodule_end()
|
---|
1484 | {
|
---|
1485 | // Desactivation du module
|
---|
1486 | if (piaskmex) delete piaskmex;
|
---|
1487 | piaskmex = NULL;
|
---|
1488 | }
|
---|