ROOT
6.07/01
Reference Guide
ROOT Home Page
Main Page
Tutorials
User's Classes
Namespaces
All Classes
Files
Release Notes
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Groups
Pages
tutorials
fitsio
FITS_tutorial1.C
Go to the documentation of this file.
1
// Open a FITS file and retrieve the first plane of the image array
2
// as a TImage object
3
void
FITS_tutorial1
()
4
{
5
printf
(
"\n\n--------------------------------\n"
);
6
printf
(
"WELCOME TO FITS tutorial #1 !!!!\n"
);
7
printf
(
"--------------------------------\n"
);
8
printf
(
"We're gonna open a FITS file that contains only the\n"
);
9
printf
(
"primary HDU, consisting on an image.\n"
);
10
printf
(
"The object you will see is a snapshot of the NGC7662 nebula,\n"
);
11
printf
(
"which was taken by the author on November 2009 in Barcelona (CATALONIA).\n\n"
);
12
13
if
(!
gROOT
->IsBatch()) {
14
//printf("Press ENTER to start..."); getchar();
15
}
16
TString
dir
=
gSystem
->
DirName
(__FILE__);
17
18
// Open primary HDU from file
19
TFITSHDU
*hdu =
new
TFITSHDU
(dir+
"/sample1.fits"
);
20
if
(hdu == 0) {
21
printf
(
"ERROR: could not access the HDU\n"
);
return
;
22
}
23
printf
(
"File successfully open!\n"
);
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
printf
(
"....................................\n"
);
32
// Here we get the exposure time.
33
//printf("Press ENTER to retrieve the exposure time from the HDU metadata..."); getchar();
34
printf
(
"Exposure time = %s\n"
, hdu->
GetKeywordValue
(
"EXPTIME"
).
Data
());
35
36
37
// Read the primary array as a matrix,
38
// selecting only layer 0.
39
// This function may be useful to
40
// do image processing.
41
printf
(
"....................................\n"
);
42
printf
(
"We can read the image as a matrix of values.\n"
);
43
printf
(
"This feature is useful to do image processing, e.g:\n"
);
44
printf
(
"histogram equalization, custom filtering, ...\n"
);
45
//printf("Press ENTER to continue..."); getchar();
46
47
TMatrixD
*mat = hdu->
ReadAsMatrix
(0);
48
mat->
Print
();
49
delete
mat;
50
51
// Read the primary array as an image,
52
// selecting only layer 0.
53
printf
(
"....................................\n"
);
54
printf
(
"Now the primary array will be read both as an image and as a histogram,\n"
);
55
printf
(
"and they will be shown in a canvas.\n"
);
56
//printf("Press ENTER to continue..."); getchar();
57
58
TImage
*im = hdu->
ReadAsImage
(0);
59
60
// Read the primary array as a histogram.
61
// Depending on array dimensions, returned
62
// histogram will be 1D, 2D or 3D
63
TH1
*hist = hdu->
ReadAsHistogram
();
64
65
66
TCanvas
*
c
=
new
TCanvas
(
"c1"
,
"FITS tutorial #1"
, 800, 300);
67
c->
Divide
(2,1);
68
c->
cd
(1);
69
im->
Draw
();
70
c->
cd
(2);
71
hist->
Draw
(
"COL"
);
72
73
// Clean up
74
delete
hdu;
75
}
76
77
c
return c
Definition:
entrylist_figure1.C:47
TFITSHDU::ReadAsMatrix
TMatrixD * ReadAsMatrix(Int_t layer=0, Option_t *opt="")
Read image HDU as a matrix.
Definition:
TFITS.cxx:865
TCanvas::cd
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition:
TCanvas.cxx:659
gROOT
#define gROOT
Definition:
TROOT.h:344
TString
Basic string class.
Definition:
TString.h:137
TSystem::DirName
virtual const char * DirName(const char *pathname)
Return the directory name in pathname.
Definition:
TSystem.cxx:980
TObject::Draw
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition:
TObject.cxx:254
TImage
An abstract interface to image processing library.
Definition:
TImage.h:45
TMatrixT< Double_t >
TString::Data
const char * Data() const
Definition:
TString.h:349
TFITSHDU
FITS file interface class.
Definition:
TFITS.h:40
TFITSHDU::Print
void Print(const Option_t *opt="") const
Print metadata.
Definition:
TFITS.cxx:722
TFITSHDU::ReadAsHistogram
TH1 * ReadAsHistogram()
Read image HDU as a histogram.
Definition:
TFITS.cxx:958
TMatrixTBase::Print
void Print(Option_t *name="") const
Print the matrix as a table of elements.
Definition:
TMatrixTBase.cxx:815
gSystem
R__EXTERN TSystem * gSystem
Definition:
TSystem.h:545
TH1::Draw
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition:
TH1.cxx:2878
TCanvas
The Canvas class.
Definition:
TCanvas.h:48
printf
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
dir
void dir(char *path=0)
Definition:
rootalias.C:30
TH1
The TH1 histogram class.
Definition:
TH1.h:80
FITS_tutorial1
void FITS_tutorial1()
Definition:
FITS_tutorial1.C:3
TPad::Divide
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:1073
TFITSHDU::GetKeywordValue
TString & GetKeywordValue(const char *keyword)
Get the value of a given keyword. Return "" if not found.
Definition:
TFITS.cxx:528
TFITSHDU::ReadAsImage
TImage * ReadAsImage(Int_t layer=0, TImagePalette *pal=0)
Read image HDU as a displayable image.
Definition:
TFITS.cxx:744