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

Last change on this file since 2902 was 2615, checked in by cmv, 21 years ago

using namespace sophya enleve de machdefs.h, nouveau sopnamsp.h cmv 10/09/2004

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