Main Page | Modules | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members | Related Pages

map2tga.cc

00001 /*
00002  *  This file is part of Healpix_cxx.
00003  *
00004  *  Healpix_cxx is free software; you can redistribute it and/or modify
00005  *  it under the terms of the GNU General Public License as published by
00006  *  the Free Software Foundation; either version 2 of the License, or
00007  *  (at your option) any later version.
00008  *
00009  *  Healpix_cxx is distributed in the hope that it will be useful,
00010  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00011  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00012  *  GNU General Public License for more details.
00013  *
00014  *  You should have received a copy of the GNU General Public License
00015  *  along with Healpix_cxx; if not, write to the Free Software
00016  *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017  *
00018  *  For more information about HEALPix, see http://healpix.jpl.nasa.gov
00019  */
00020 
00021 /*
00022  *  Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik
00023  *  and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
00024  *  (DLR).
00025  */
00026 
00027 /*
00028  *  Copyright (C) 2003, 2004 Max-Planck-Society
00029  *  Author: Martin Reinecke
00030  */
00031 
00032 #include <sstream>
00033 #include <iomanip>
00034 #include "healpix_map.h"
00035 #include "healpix_map_fitsio.h"
00036 #include "rotmatrix.h"
00037 #include "pointing.h"
00038 #include "tga_image.h"
00039 #include "paramfile.h"
00040 
00041 using namespace std;
00042 
00043 void histo_eq (arr2<float> &img, float &minv, float &maxv, arr<double> &newpos)
00044   {
00045   const int nbins=100;
00046   arr<int> bincnt (nbins);
00047   bincnt.fill(0);
00048   int pixels=0;
00049 
00050   double fact = 1./(maxv-minv);
00051   for (int i=0; i<img.size1(); ++i)
00052     for (int j=0; j<img.size2(); ++j)
00053       if (img[i][j]>-1e30)
00054         {
00055         img[i][j] = (img[i][j]-minv)*fact; 
00056         int idx = int(img[i][j]*nbins);
00057         idx=max(0,min(idx,nbins-1));
00058         ++bincnt[idx];
00059         ++pixels;
00060         }
00061 
00062   newpos.alloc(nbins+1);
00063   int accu=0;
00064   for (int m=0; m<nbins; ++m)
00065     {
00066     newpos[m] = double(accu)/pixels;
00067     accu += bincnt[m];
00068     }
00069   newpos[nbins]=1.;
00070 
00071   for (int i=0; i<img.size1(); ++i)
00072     for (int j=0; j<img.size2(); ++j)
00073       if (img[i][j]>-1e30)
00074         {
00075         int idx = int(img[i][j]*nbins);
00076         idx=max(0,min(idx,nbins-1));
00077         double frac = nbins*img[i][j] - idx;
00078         img[i][j] = (1-frac)*newpos[idx] + frac*newpos[idx+1];
00079         img[i][j] = minv+(maxv-minv)*img[i][j];
00080         }
00081   }
00082 
00083 void pro_mollw (const Healpix_Map<float> &map, double lon0, double lat0,
00084   int xsize, arr2<float> &img, float &minv, float &maxv, int smooth)
00085   {
00086   int ysize=xsize/2;
00087   img.alloc(xsize,ysize);
00088   img.fill(-1e35);
00089   double xc=(xsize-1)/2., yc=(ysize-1)/2.;
00090   double lon0rad = lon0*degr2rad;
00091   double lat0rad = lat0*degr2rad;
00092 
00093   rotmatrix rot;
00094   rot.Make_CPAC_Euler_Matrix(0,-lat0rad,-lon0rad);
00095 
00096   minv=1e30;
00097   maxv=-1e30;
00098   for (int i=0; i<img.size1(); ++i)
00099     for (int j=0; j<img.size2(); ++j)
00100       {
00101       double u = 2*(i-xc)/(xc/1.02);
00102       double v = (j-yc)/(yc/1.02);
00103       bool mask = ((u*u/4 + v*v) <= 1);
00104       if (mask)
00105         {
00106         pointing ptg (halfpi-(asin(2/pi*(asin(v) + v*sqrt((1-v)*(1+v))))),
00107                       -halfpi*u/max(sqrt((1-v)*(1+v)),1e-6));
00108         vec3 pnt = rot.Transform(ptg.to_vec3());
00109         if (smooth==1)
00110           img[i][j] = map.interpolated_value(pnt);
00111         else if (smooth==2)
00112           img[i][j] = map.interpolated_value2(pnt);
00113         else
00114           img[i][j]=map[map.ang2pix(pnt)];
00115         if (!approx(img[i][j],Healpix_undef))
00116           {
00117           if (img[i][j]<minv) minv=img[i][j];
00118           if (img[i][j]>maxv) maxv=img[i][j];
00119           }
00120         }
00121       }
00122   }
00123 
00124 void pro_gno (const Healpix_Map<float> &map, double lon0, double lat0,
00125   int xsize, double resgrid, arr2<float> &img, float &minv, float &maxv,
00126   int smooth)
00127   {
00128   double lon0rad = lon0*degr2rad;
00129   double lat0rad = lat0*degr2rad;
00130 
00131   rotmatrix rot;
00132   rot.Make_CPAC_Euler_Matrix(lon0rad,-lat0rad,0);
00133 
00134   double delta=resgrid*degr2rad/60.;
00135   double start=-(xsize/2.)*delta;
00136   img.alloc(xsize,xsize);
00137   minv=1e30;
00138   maxv=-1e30;
00139   for (int i=0; i<img.size1(); ++i)
00140     for (int j=0; j<img.size2(); ++j)
00141       {
00142       vec3 pnt (1,-(start+i*delta), start+j*delta);
00143       pnt = rot.Transform(pnt);
00144         if (smooth==1)
00145           img[i][j] = map.interpolated_value(pnt);
00146         else if (smooth==2)
00147           img[i][j] = map.interpolated_value2(pnt);
00148         else
00149           img[i][j]=map[map.ang2pix(pnt)];
00150       if (!approx(img[i][j],Healpix_undef))
00151         {
00152         if (img[i][j]<minv) minv=img[i][j];
00153         if (img[i][j]>maxv) maxv=img[i][j];
00154         }
00155       }
00156   }
00157 
00158 void colorbar (TGA_Image &img, const Palette &pal, int xmin, int xmax,
00159   int ymin, int ymax, bool flippal, const arr<double> &newpos)
00160   {
00161   int nbins = newpos.size()-1;
00162   for (int i=xmin; i<=xmax; ++i)
00163     {
00164     double val = (double(i)-xmin)/(xmax-xmin);
00165     if (nbins>0)
00166       {
00167       int idx = int(val*nbins);
00168       idx=max(0,min(idx,nbins-1));
00169       double frac = nbins*val - idx;
00170       val = (1-frac)*newpos[idx] + frac*newpos[idx+1];
00171       }  
00172     if (flippal) val=1-val;
00173     Colour c = pal.Get_Colour(val);
00174     for (int j=ymin; j<=ymax; ++j)
00175       img.put_pixel(i,j,c);
00176     }
00177   }
00178 
00179 string conv (double val)
00180   {
00181   ostringstream os;
00182   if (abs(val)>100 || abs(val)<0.01)
00183     {
00184     os << setw(10) << setprecision(3) << scientific << val;
00185     return os.str();
00186     }
00187   os << setw(10) << setprecision(6) << fixed << val;
00188   return trim(os.str());
00189   }
00190 
00191 void usage()
00192   {
00193   cout <<
00194     "\nUsage:\n"
00195     "  map2tga <parameter file>\n\n"
00196     "or:\n"
00197     "  map2tga <input file> <output file> [-sig <int>] [-pal <int>]\n"
00198     "    [-xsz <int>] [-bar] [-log] [-asinh] [-lon <float>] [-lat <float>]\n"
00199     "    [-mul <float>] [-add <float>] [-min <float>] [-max <float>]\n"
00200     "    [-res <float>] [-title <string>] [-flippal] [-gnomonic]\n"
00201     "    [-interpol <0|1|2>] [-equalize] \n\n";
00202   exit(1);
00203   }
00204 
00205 int main (int argc, const char **argv)
00206   {
00207 PLANCK_DIAGNOSIS_BEGIN
00208   announce ("map2tga");
00209   if (argc<2) usage();
00210 
00211   string infile = "";
00212   string outfile = "";
00213   string title = "";
00214   int xres = 1024;
00215   bool bar = false, mollpro=true;
00216   double lat0=0, lon0=0, res=1;
00217   double usermin=-1e30, usermax=1e30, offset=0, factor=1;
00218   bool min_supplied=false, max_supplied=false;
00219   bool logflag=false, asinhflag=false, eqflag=false;
00220   int palnr = 4;
00221   int colnum=1;
00222   bool flippal = false;
00223   int interpol=0;
00224 
00225   if (argc>2)
00226     {
00227     infile = argv[1];
00228     outfile = argv[2];
00229     int m=3;
00230     while (m<argc)
00231       {
00232       int mstart = m;
00233       string arg = argv[m];
00234       if (arg == "-sig") { stringToData(argv[m+1],colnum); m+=2; }
00235       if (arg == "-pal") { stringToData(argv[m+1],palnr); m+=2; }
00236       if (arg == "-xsz") { stringToData(argv[m+1],xres); m+=2; }
00237       if (arg == "-bar") { bar=true; ++m; }
00238       if (arg == "-log") { logflag=true; ++m; }
00239       if (arg == "-equalize") { eqflag=true; ++m; }
00240       if (arg == "-asinh") { asinhflag=true; ++m; }
00241       if (arg == "-lon") { stringToData(argv[m+1],lon0); m+=2; }
00242       if (arg == "-lat") { stringToData(argv[m+1],lat0); m+=2; }
00243       if (arg == "-mul") { stringToData(argv[m+1],factor); m+=2; }
00244       if (arg == "-add") { stringToData(argv[m+1],offset); m+=2; }
00245       if (arg == "-min")
00246         {
00247         stringToData(argv[m+1],usermin);
00248         min_supplied=true;
00249         m+=2;
00250         }
00251       if (arg == "-max")
00252         {
00253         stringToData(argv[m+1],usermax);
00254         max_supplied=true;
00255         m+=2;
00256         }
00257       if (arg == "-res") { stringToData(argv[m+1],res); m+=2; }
00258       if (arg == "-title") { title = argv[m+1]; m+=2; }
00259       if (arg == "-flippal") { flippal=true; ++m; }
00260       if (arg == "-gnomonic") { mollpro=false; ++m; }
00261       if (arg == "-interpol") { stringToData(argv[m+1],interpol); m+=2; }
00262 
00263       if (mstart==m)
00264         {
00265         cout << "unrecognized option: " + arg << endl;
00266         usage();
00267         }
00268       }
00269     }
00270   else
00271     {
00272     paramfile params (argv[1]);
00273     infile = params.find<string>("infile");
00274     outfile = params.find<string>("outfile");
00275     colnum = params.find<int>("sig",colnum);
00276     palnr = params.find<int>("pal",palnr);
00277     xres = params.find<int>("xsz",xres);
00278     bar = params.find<bool>("bar",bar);
00279     logflag = params.find<bool>("log",logflag);
00280     eqflag = params.find<bool>("equalize",eqflag);
00281     asinhflag = params.find<bool>("asinh",asinhflag);
00282     lon0 = params.find<double>("lon",lon0);
00283     lat0 = params.find<double>("lat",lat0);
00284     factor = params.find<double>("mul",factor);
00285     offset = params.find<double>("add",offset);
00286     usermin = params.find<double>("min", usermin);
00287     if (usermin>-1e29) min_supplied = true;
00288     usermax = params.find<double>("max", usermax);
00289     if (usermax<1e29) max_supplied = true;
00290     res = params.find<double>("res",res);
00291     title = params.find<string>("title","");
00292     flippal = params.find<bool>("flippal",flippal);
00293     interpol = params.find<int>("interpol",interpol);
00294     string tmp = params.find<string>("pro","");
00295     if (tmp == "gno") mollpro=false;
00296     }
00297 
00298   Healpix_Map<float> map(0,RING);
00299   read_Healpix_map_from_fits(infile,map,colnum,2);
00300   for (int m=0; m<map.Npix(); ++m)
00301     {
00302     if (!approx(map[m],Healpix_undef))
00303       {
00304       map[m] = (map[m]+offset)*factor;
00305       if (logflag)
00306         {
00307         if (map[m]<=0)
00308           map[m] = Healpix_undef;
00309         else
00310           map[m] = log(double(map[m]))/ln10;
00311         }
00312       if (asinhflag)
00313         {
00314         if (map[m]>=0)
00315           map[m] = log(double(map[m]+sqrt(map[m]*map[m]+1)));
00316         else
00317           map[m] = -log(double(-map[m]+sqrt(map[m]*map[m]+1)));
00318         }
00319       if (min_supplied) if (map[m] < usermin) map[m] = usermin;
00320       if (max_supplied) if (map[m] > usermax) map[m] = usermax;
00321       }
00322     }
00323 
00324   arr2<float> imgarr;
00325   float minv, maxv;
00326   if (mollpro)
00327     pro_mollw (map,lon0,lat0,xres,imgarr,minv,maxv,interpol);
00328   else
00329     pro_gno (map,lon0,lat0,xres,res,imgarr,minv,maxv,interpol);
00330 
00331   arr<double> newpos;
00332   if (eqflag) histo_eq(imgarr,minv,maxv,newpos);
00333 
00334   if (min_supplied) minv = usermin;
00335   if (max_supplied) maxv = usermax;
00336   if (maxv==minv) maxv=minv+1e-10;
00337 
00338   int xsz=imgarr.size1();
00339   int ysz=imgarr.size2();
00340   int yofs=ysz;
00341   int scale = max(1,xsz/800);
00342   if (bar) ysz+=60*scale;
00343   if (title != "") { ysz+=50*scale; yofs+=50*scale; }
00344   TGA_Image img(xsz,ysz);
00345   img.fill(Colour(1,1,1));
00346   img.set_font (giant_font);
00347   Palette pal;
00348   pal.setPredefined(palnr);
00349   if (title != "")
00350     img.annotate_centered(xsz/2,25*scale,Colour(0,0,0),title,scale);
00351   for (int i=0; i<imgarr.size1(); ++i)
00352     for (int j=0; j<imgarr.size2(); ++j)
00353       {
00354       if (imgarr[i][j]>-1e32)
00355         {
00356         Colour c(0.5,0.5,0.5);
00357         if (!approx (imgarr[i][j],Healpix_undef))
00358           {
00359           int col = int((imgarr[i][j]-minv)/(maxv-minv)*256);
00360           col = min(255,max(col,0));
00361           float colfrac = (imgarr[i][j]-minv)/(maxv-minv);
00362           if (flippal) colfrac = 1-colfrac;
00363           c = pal.Get_Colour(colfrac);
00364           }
00365         img.put_pixel(i,yofs-j-1,c);
00366         }
00367       }
00368 
00369   if (bar)
00370     {
00371     colorbar (img,pal,xsz/10,(xsz*9)/10,ysz-40*scale,ysz-20*scale,flippal,
00372       newpos);
00373     img.set_font (medium_bold_font);
00374     img.annotate_centered (xsz/20,ysz-30*scale,Colour(0,0,0),conv(minv),scale);
00375     img.annotate_centered ((xsz*19)/20,ysz-30*scale,Colour(0,0,0),conv(maxv),
00376       scale);
00377     }
00378   img.write(outfile);
00379 PLANCK_DIAGNOSIS_END
00380   }

Generated on Fri Jul 8 09:37:15 2005 for Healpix C++