Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
FITS_tutorial1.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_FITS
3 /// \notebook -draw
4 ///
5 /// Open a FITS file and retrieve the first plane of the image array
6 /// as a TImage object
7 ///
8 /// \macro_image
9 /// \macro_code
10 /// \macro_output
11 ///
12 /// \author Claudi Martinez
13 
14 void FITS_tutorial1()
15 {
16  // Here we open a FITS file that contains only the primary HDU, consisting on an image.
17  // The object you will see is a snapshot of the NGC7662 nebula,
18  // which was taken by the author on November 2009 in Barcelona.
19 
20  TString dir = gROOT->GetTutorialDir();
21 
22  // Open primary HDU from file
23  TFITSHDU hdu(dir + "/fitsio/sample1.fits");
24 
25  // Dump the HDUs within the FITS file
26  // and also their metadata
27  // printf("Press ENTER to see summary of all data stored in the file:"); getchar();
28 
29  hdu.Print("F+");
30 
31  // Here we get the exposure time.
32  printf("Exposure time = %s\n", hdu.GetKeywordValue("EXPTIME").Data());
33 
34  // Read the primary array as a matrix, selecting only layer 0.
35  // This function may be useful to do image processing, e.g. custom filtering
36 
37  std::unique_ptr<TMatrixD> mat(hdu.ReadAsMatrix(0));
38  mat->Print();
39 
40  // Read the primary array as an image, selecting only layer 0.
41 
42  std::unique_ptr<TImage> im(hdu.ReadAsImage(0));
43 
44  // Read the primary array as a histogram. Depending on array dimensions, the returned
45  // histogram will be 1D, 2D or 3D.
46  std::unique_ptr<TH1> hist(hdu.ReadAsHistogram());
47 
48  auto c = new TCanvas("c1", "FITS tutorial #1", 1400, 800);
49  c->Divide(2, 1);
50  c->cd(1);
51  im->DrawClone();
52  c->cd(2);
53  hist->DrawClone("COL");
54 }