Logo ROOT   6.14/05
Reference Guide
FITS_tutorial1.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_FITS
3 /// Open a FITS file and retrieve the first plane of the image array
4 /// as a TImage object
5 ///
6 /// \macro_code
7 ///
8 /// \author Claudi Martinez
9 
10 void FITS_tutorial1()
11 {
12  printf("\n\n--------------------------------\n");
13  printf("WELCOME TO FITS tutorial #1 !!!!\n");
14  printf("--------------------------------\n");
15  printf("We're gonna open a FITS file that contains only the\n");
16  printf("primary HDU, consisting on an image.\n");
17  printf("The object you will see is a snapshot of the NGC7662 nebula,\n");
18  printf("which was taken by the author on November 2009 in Barcelona (CATALONIA).\n\n");
19 
20  if (!gROOT->IsBatch()) {
21  //printf("Press ENTER to start..."); getchar();
22  }
23  TString dir = gSystem->DirName(__FILE__);
24 
25  // Open primary HDU from file
26  TFITSHDU *hdu = new TFITSHDU(dir+"/sample1.fits");
27  if (hdu == 0) {
28  printf("ERROR: could not access the HDU\n"); return;
29  }
30  printf("File successfully open!\n");
31 
32  // Dump the HDUs within the FITS file
33  // and also their metadata
34  //printf("Press ENTER to see summary of all data stored in the file:"); getchar();
35 
36  hdu->Print("F+");
37 
38  printf("....................................\n");
39  // Here we get the exposure time.
40  //printf("Press ENTER to retrieve the exposure time from the HDU metadata..."); getchar();
41  printf("Exposure time = %s\n", hdu->GetKeywordValue("EXPTIME").Data());
42 
43 
44  // Read the primary array as a matrix,
45  // selecting only layer 0.
46  // This function may be useful to
47  // do image processing.
48  printf("....................................\n");
49  printf("We can read the image as a matrix of values.\n");
50  printf("This feature is useful to do image processing, e.g:\n");
51  printf("histogram equalization, custom filtering, ...\n");
52  //printf("Press ENTER to continue..."); getchar();
53 
54  TMatrixD *mat = hdu->ReadAsMatrix(0);
55  mat->Print();
56  delete mat;
57 
58  // Read the primary array as an image,
59  // selecting only layer 0.
60  printf("....................................\n");
61  printf("Now the primary array will be read both as an image and as a histogram,\n");
62  printf("and they will be shown in a canvas.\n");
63  //printf("Press ENTER to continue..."); getchar();
64 
65  TImage *im = hdu->ReadAsImage(0);
66 
67  // Read the primary array as a histogram.
68  // Depending on array dimensions, returned
69  // histogram will be 1D, 2D or 3D
70  TH1 *hist = hdu->ReadAsHistogram();
71 
72 
73  TCanvas *c = new TCanvas("c1", "FITS tutorial #1", 800, 300);
74  c->Divide(2,1);
75  c->cd(1);
76  im->Draw();
77  c->cd(2);
78  hist->Draw("COL");
79 
80  // Clean up
81  delete hdu;
82 }
TMatrixD * ReadAsMatrix(Int_t layer=0, Option_t *opt="")
Read image HDU as a matrix.
Definition: TFITS.cxx:896
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
#define gROOT
Definition: TROOT.h:410
Basic string class.
Definition: TString.h:131
virtual const char * DirName(const char *pathname)
Return the directory name in pathname.
Definition: TSystem.cxx:1004
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:195
An abstract interface to image processing library.
Definition: TImage.h:29
TMatrixT.
Definition: TMatrixDfwd.h:22
void Print(const Option_t *opt="") const
Print metadata.
Definition: TFITS.cxx:753
FITS file interface class.
Definition: TFITS.h:34
TH1 * ReadAsHistogram()
Read image HDU as a histogram.
Definition: TFITS.cxx:989
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
The Canvas class.
Definition: TCanvas.h:31
The TH1 histogram class.
Definition: TH1.h:56
void Print(Option_t *name="") const
Print the matrix as a table of elements.
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
TString & GetKeywordValue(const char *keyword)
Get the value of a given keyword. Return "" if not found.
Definition: TFITS.cxx:559
#define c(i)
Definition: RSha256.hxx:101
TImage * ReadAsImage(Int_t layer=0, TImagePalette *pal=0)
Read image HDU as a displayable image.
Definition: TFITS.cxx:775
const char * Data() const
Definition: TString.h:364