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