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

Last change on this file since 2218 was 2211, checked in by ansari, 23 years ago

Compil SGI - 2 Reza - 12/10

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