source: Sophya/trunk/SophyaPI/ProgPI/skymapmodule.cc@ 2468

Last change on this file since 2468 was 2468, checked in by ansari, 22 years ago

Ajout commentaire groupe de help - Reza 2/12/2003

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