Logo ROOT   6.14/05
Reference Guide
TFITS.cxx
Go to the documentation of this file.
1 // @(#)root/graf2d:$Id$
2 // Author: Claudi Martinez, July 19th 2010
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2010, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /// \defgroup fitsio FITS file
13 /// \brief Interface to FITS file.
14 /// \ingroup Graphics2D
15 ///
16 /// TFITS is an interface that lets you reading Flexible Image Transport System
17 /// (FITS) files, which are generally used in astronomy. This file format
18 /// was standardized 1981 and today is still widely used among professional
19 /// and amateur astronomers. FITS is not only an image file, but also
20 /// it can contain spectrums, data tables, histograms, and multidimensional
21 /// data. Furthermore, FITS data can be described itself by containing
22 /// human-readable information that let us to interpret the data within
23 /// the FITS file. For example, a FITS could contain a 3D data cube,
24 /// but an additional description would tell us that we must read it, for
25 /// example, as a 3-layer image.
26 ///
27 /// TFITS requires CFITSIO library to be installed on your system. It
28 /// is currently maintained by NASA/GSFC and can be downloaded from
29 /// [NASA/GSFC web site](http://fits.gsfc.nasa.gov), as well as documentation.
30 ///
31 /// Using this interface is easy and straightforward. There is only 1 class
32 /// called "TFITSHDU" which has several methods to extract data from a
33 /// FITS file, more specifically, from an HDU within the file. An HDU, or
34 /// Header Data Unit, is a chunk of data with a header containing several
35 /// "keyword = value" tokens. The header describes the structure of data
36 /// within the HDU. An HDU can be of two types: an "image HDU" or a "table
37 /// HDU". The former can be any kind of multidimensional array of real numbers,
38 /// by which the name "image" may be confusing: you can store an image, but
39 /// you can also store a N-dimensional data cube. On the other hand, table
40 /// HDUs are sets of several rows and columns (a.k.a fields) which contain
41 /// generic data, as strings, real or complex numbers and even arrays.
42 ///
43 /// Please have a look to the tutorials ($ROOTSYS/tutorials/fitsio/) to see
44 /// some examples. IMPORTANT: to run tutorials it is required that
45 /// you change the current working directory of ROOT (CINT) shell to the
46 /// tutorials directory. Example:
47 /// ~~~ {.cpp}
48 /// root [1] gSystem->ChangeDirectory("tutorials/fitsio")
49 /// root [1] .x FITS_tutorial1.C
50 /// ~~~
51 /// LIST OF TODO
52 /// - Support for complex values within data tables
53 /// - Support for reading arrays from table cells
54 /// - Support for grouping
55 ///
56 /// IMPLEMENTATION NOTES:
57 ///
58 /// CFITSIO library uses standard C types ('int', 'long', ...). To avoid
59 /// confusion, the same types are used internally by the access methods.
60 /// However, class's fields are ROOT-defined types.
61 
62 /** \class TFITSHDU
63 \ingroup fitsio
64 
65 FITS file interface class
66 
67 TFITSHDU is a class that allows extracting images and data from FITS files and contains
68 several methods to manage them.
69 */
70 
71 #include "TFITS.h"
72 #include "TROOT.h"
73 #include "TImage.h"
74 #include "TArrayI.h"
75 #include "TArrayD.h"
76 #include "TH1D.h"
77 #include "TH2D.h"
78 #include "TH3D.h"
79 #include "TVectorD.h"
80 #include "TMatrixD.h"
81 #include "TObjArray.h"
82 #include "TObjString.h"
83 #include "TCanvas.h"
84 
85 #include "fitsio.h"
86 #include <stdlib.h>
87 
88 
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Clean path from possible filter and put the result in 'dst'.
93 
94 void TFITSHDU::CleanFilePath(const char *filepath_with_filter, TString &dst)
95 {
96  dst = filepath_with_filter;
97 
98  Ssiz_t ndx = dst.Index("[", 1, 0, TString::kExact);
99  if (ndx != kNPOS) {
100  dst.Resize(ndx);
101  }
102 }
103 
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// TFITSHDU constructor from file path with HDU selection filter.
107 /// Please refer to CFITSIO manual for more information about
108 /// HDU selection filters.
109 ///
110 /// Examples:
111 /// - `TFITSHDU("/path/to/myfile.fits")`: just open the PRIMARY HDU
112 /// - `TFITSHDU("/path/to/myfile.fits[1]")`: open HDU #1
113 /// - `TFITSHDU("/path/to/myfile.fits[PICS]")`: open HDU called 'PICS'
114 /// - `TFITSHDU("/path/to/myfile.fits[ACQ][EXPOSURE > 5]")`: open the (table) HDU called 'ACQ' and
115 /// selects the rows that have column 'EXPOSURE'
116 /// greater than 5.
117 
118 TFITSHDU::TFITSHDU(const char *filepath_with_filter)
119 {
120  _initialize_me();
121 
122  fFilePath = filepath_with_filter;
123  CleanFilePath(filepath_with_filter, fBaseFilePath);
124 
125  if (kFALSE == LoadHDU(fFilePath)) {
127  throw -1;
128  }
129 }
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// TFITSHDU constructor from filepath and extension number.
133 
134 TFITSHDU::TFITSHDU(const char *filepath, Int_t extension_number)
135 {
136  _initialize_me();
137  CleanFilePath(filepath, fBaseFilePath);
138 
139  //Add "by extension number" filter
140  fFilePath.Form("%s[%d]", fBaseFilePath.Data(), extension_number);
141 
142  if (kFALSE == LoadHDU(fFilePath)) {
144  throw -1;
145  }
146 }
147 
148 ////////////////////////////////////////////////////////////////////////////////
149 /// TFITSHDU constructor from filepath and extension name.
150 
151 TFITSHDU::TFITSHDU(const char *filepath, const char *extension_name)
152 {
153  _initialize_me();
154  CleanFilePath(filepath, fBaseFilePath);
155 
156  //Add "by extension number" filter
157  fFilePath.Form("%s[%s]", fBaseFilePath.Data(), extension_name);
158 
159 
160  if (kFALSE == LoadHDU(fFilePath)) {
162  throw -1;
163  }
164 }
165 
166 ////////////////////////////////////////////////////////////////////////////////
167 /// TFITSHDU destructor.
168 
170 {
172 }
173 
174 ////////////////////////////////////////////////////////////////////////////////
175 /// Release internal resources.
176 
178 {
179  if (fRecords) delete [] fRecords;
180 
181  if (fType == kImageHDU) {
182  if (fSizes) delete fSizes;
183  if (fPixels) delete fPixels;
184  } else {
185  if (fColumnsInfo) {
186  if (fCells) {
187  for (Int_t i = 0; i < fNColumns; i++) {
188  if (fColumnsInfo[i].fType == kString) {
189  //Deallocate character arrays allocated for kString columns
190  Int_t offset = i * fNRows;
191  for (Int_t row = 0; row < fNRows; row++) {
192  delete [] fCells[offset+row].fString;
193  }
194  } else if (fColumnsInfo[i].fType == kRealVector) {
195  //Deallocate character arrays allocated for kString columns
196  Int_t offset = i * fNRows;
197  for (Int_t row = 0; row < fNRows; row++) {
198  delete [] fCells[offset+row].fRealVector;
199  }
200  }
201  }
202 
203  delete [] fCells;
204  }
205 
206  delete [] fColumnsInfo;
207  }
208 
209 
210  }
211 }
212 
213 ////////////////////////////////////////////////////////////////////////////////
214 /// Do some initializations.
215 
217 {
218  fRecords = 0;
219  fPixels = 0;
220  fSizes = 0;
221  fColumnsInfo = 0;
222  fNColumns = fNRows = 0;
223  fCells = 0;
224 }
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// Load HDU from fits file satisfying the specified filter.
228 /// Returns kTRUE if success. Otherwise kFALSE.
229 /// If filter == "" then the primary array is selected
230 
232 {
233  fitsfile *fp=0;
234  int status = 0;
235  char errdescr[FLEN_STATUS+1];
236 
237  // Open file with filter
238  fits_open_file(&fp, filepath_filter.Data(), READONLY, &status);
239  if (status) goto ERR;
240 
241  // Read HDU number
242  int hdunum;
243  fits_get_hdu_num(fp, &hdunum);
244  fNumber = Int_t(hdunum);
245 
246  // Read HDU type
247  int hdutype;
248  fits_get_hdu_type(fp, &hdutype, &status);
249  if (status) goto ERR;
250  fType = (hdutype == IMAGE_HDU) ? kImageHDU : kTableHDU;
251 
252  //Read HDU header records
253  int nkeys, morekeys;
254  char keyname[FLEN_KEYWORD+1];
255  char keyvalue[FLEN_VALUE+1];
256  char comment[FLEN_COMMENT+1];
257 
258  fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
259  if (status) goto ERR;
260 
261  fRecords = new struct HDURecord[nkeys];
262 
263  for (int i = 1; i <= nkeys; i++) {
264  fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
265  if (status) goto ERR;
266  fRecords[i-1].fKeyword = keyname;
267  fRecords[i-1].fValue = keyvalue;
268  fRecords[i-1].fComment = comment;
269  }
270 
271  fNRecords = Int_t(nkeys);
272 
273  //Set extension name
274  fExtensionName = "PRIMARY"; //Default
275  for (int i = 0; i < nkeys; i++) {
276  if (fRecords[i].fKeyword == "EXTNAME") {
278  break;
279  }
280  }
281 
282  //Read HDU's data
283  if (fType == kImageHDU) {
284  //Image
285  int param_ndims=0;
286  long *param_dimsizes;
287 
288  //Read image number of dimensions
289  fits_get_img_dim(fp, &param_ndims, &status);
290  if (status) goto ERR;
291  if (param_ndims > 0) {
292  //Read image sizes in each dimension
293  param_dimsizes = new long[param_ndims];
294  fits_get_img_size(fp, param_ndims, param_dimsizes, &status);
295  if (status) {
296  delete [] param_dimsizes;
297  goto ERR;
298  }
299 
300  fSizes = new TArrayI(param_ndims);
301  fSizes = new TArrayI(param_ndims);
302  for (int i = 0; i < param_ndims; i++) { //Use for loop to copy values instead of passing array to constructor, since 'Int_t' size may differ from 'long' size
303  fSizes->SetAt(param_dimsizes[i], i);
304  }
305 
306  delete [] param_dimsizes;
307 
308  //Read pixels
309  int anynul;
310  long *firstpixel = new long[param_ndims];
311  double nulval = 0;
312  long npixels = 1;
313 
314  for (int i = 0; i < param_ndims; i++) {
315  npixels *= (long) fSizes->GetAt(i); //Compute total number of pixels
316  firstpixel[i] = 1; //Set first pixel to read from.
317  }
318 
319  double *pixels = new double[npixels];
320 
321  fits_read_pix(fp, TDOUBLE, firstpixel, npixels,
322  (void *) &nulval, (void *) pixels, &anynul, &status);
323 
324  if (status) {
325  delete [] firstpixel;
326  delete [] pixels;
327  goto ERR;
328  }
329 
330  fPixels = new TArrayD(npixels, pixels);
331 
332  delete [] firstpixel;
333  delete [] pixels;
334 
335  } else {
336  //Null array
337  fSizes = new TArrayI();
338  fPixels = new TArrayD();
339  }
340  } else {
341  // Table
342 
343  // Get table's number of rows and columns
344  long table_rows;
345  int table_cols;
346 
347  fits_get_num_rows(fp, &table_rows, &status);
348  if (status) goto ERR;
349 
350  fNRows = Int_t(table_rows);
351 
352  fits_get_num_cols(fp, &table_cols, &status);
353  if (status) goto ERR;
354 
355  fNColumns = Int_t(table_cols);
356 
357  // Allocate column info array
358  fColumnsInfo = new struct Column[table_cols];
359 
360  // Read column names
361  char colname[80];
362  int colnum;
363 
364  fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
365  while (status == COL_NOT_UNIQUE)
366  {
367  fColumnsInfo[colnum-1].fName = colname;
368  fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
369  }
370  if (status != COL_NOT_FOUND) goto ERR;
371  status = 0;
372 
373  //Allocate cells
374  fCells = new union Cell [table_rows * table_cols];
375 
376  // Read columns
377  int typecode;
378  long repeat, width;
379  Int_t cellindex;
380 
381 
382  for (colnum = 0, cellindex = 0; colnum < fNColumns; colnum++) {
383  fits_get_coltype(fp, colnum+1, &typecode, &repeat, &width, &status);
384 
385  if (status) goto ERR;
386 
387  if ((typecode == TDOUBLE) || (typecode == TSHORT) || (typecode == TLONG)
388  || (typecode == TFLOAT) || (typecode == TLOGICAL) || (typecode == TBIT)
389  || (typecode == TBYTE) || (typecode == TSTRING)) {
390 
391  fColumnsInfo[colnum].fType = (typecode == TSTRING) ? kString : kRealNumber;
392 
393  if (typecode == TSTRING) {
394  // String column
395  int dispwidth=0;
396  fits_get_col_display_width(fp, colnum+1, &dispwidth, &status);
397  if (status) goto ERR;
398 
399 
400  char *nulval = (char*) "";
401  int anynul=0;
402  char **array;
403 
404  if (dispwidth <= 0) {
405  dispwidth = 1;
406  }
407 
408  array = new char* [table_rows];
409  for (long row = 0; row < table_rows; row++) {
410  array[row] = new char[dispwidth+1]; //also room for end null!
411  }
412 
413  if (repeat > 0) {
414  fits_read_col(fp, TSTRING, colnum+1, 1, 1, table_rows, nulval, array, &anynul, &status);
415  if (status) {
416  for (long row = 0; row < table_rows; row++) {
417  delete [] array[row];
418  }
419  delete [] array;
420  goto ERR;
421  }
422 
423  } else {
424  //No elements: set dummy
425  for (long row = 0; row < table_rows; row++) {
426  strlcpy(array[row], "-",dispwidth+1);
427  }
428  }
429 
430 
431  //Save values
432  for (long row = 0; row < table_rows; row++) {
433  fCells[cellindex++].fString = array[row];
434  }
435 
436  delete [] array; //Delete temporal array holding pointer to strings, but not delete strings themselves!
437 
438 
439  } else {
440  //Numeric or vector column
441  double nulval = 0;
442  int anynul=0;
443 
444  fColumnsInfo[colnum].fDim = (Int_t) repeat;
445 
446  double *array = 0;
447  char *arrayl = 0;
448 
449  if (repeat > 0) {
450 
451  if (typecode == TLOGICAL) {
452  arrayl = new char[table_rows * repeat];
453  fits_read_col(fp, TLOGICAL, colnum + 1, 1, 1, table_rows * repeat, &nulval, arrayl, &anynul,
454  &status);
455  if (status) {
456  delete[] arrayl;
457  goto ERR;
458  }
459  } else {
460  array = new double[table_rows * repeat]; // Hope you got a big machine! Ask China otherwise :-)
461  fits_read_col(fp, TDOUBLE, colnum + 1, 1, 1, table_rows * repeat, &nulval, array, &anynul,
462  &status);
463  if (status) {
464  delete[] array;
465  goto ERR;
466  }
467  }
468 
469  } else {
470  // No elements: set dummy
471  array = new double[table_rows];
472  for (long row = 0; row < table_rows; row++) {
473  array[row] = 0.0;
474  }
475  }
476 
477  // Save values
478  if (repeat == 1) {
479  // Scalar
480  if (typecode == TLOGICAL) {
481  for (long row = 0; row < table_rows; row++) {
482  int temp = (signed char)arrayl[row];
483  fCells[cellindex++].fRealNumber = (double)temp;
484  }
485  delete[] arrayl;
486  } else {
487  for (long row = 0; row < table_rows; row++) {
488  fCells[cellindex++].fRealNumber = array[row];
489  }
490  delete[] array;
491  }
492  } else if (repeat > 1) {
493  // Vector
494  if (typecode == TLOGICAL) {
495  for (long row = 0; row < table_rows; row++) {
496  double *vec = new double[repeat];
497  long offset = row * repeat;
498  for (long component = 0; component < repeat; component++) {
499  int temp = (signed char)arrayl[offset++];
500  vec[component] = (double)temp;
501  }
502  fCells[cellindex++].fRealVector = vec;
503  }
504  delete[] arrayl;
505  } else {
506  for (long row = 0; row < table_rows; row++) {
507  double *vec = new double[repeat];
508  long offset = row * repeat;
509  for (long component = 0; component < repeat; component++) {
510  vec[component] = array[offset++];
511  }
512  fCells[cellindex++].fRealVector = vec;
513  }
514  delete[] array;
515  }
516  }
517  }
518  } else {
519  Warning("LoadHDU", "error opening FITS file. Column type %d is currently not supported", typecode);
520  }
521  }
522 
523  if (hdutype == ASCII_TBL) {
524  // ASCII table
525 
526  } else {
527  // Binary table
528  }
529  }
530 
531  // Close file
532  fits_close_file(fp, &status);
533  return kTRUE;
534 
535 ERR:
536  fits_get_errstatus(status, errdescr);
537  Warning("LoadHDU", "error opening FITS file. Details: %s", errdescr);
538  status = 0;
539  if (fp) fits_close_file(fp, &status);
540  return kFALSE;
541 }
542 
543 ////////////////////////////////////////////////////////////////////////////////
544 /// Get record by keyword.
545 
546 struct TFITSHDU::HDURecord* TFITSHDU::GetRecord(const char *keyword)
547 {
548  for (int i = 0; i < fNRecords; i++) {
549  if (fRecords[i].fKeyword == keyword) {
550  return &fRecords[i];
551  }
552  }
553  return 0;
554 }
555 
556 ////////////////////////////////////////////////////////////////////////////////
557 /// Get the value of a given keyword. Return "" if not found.
558 
559 TString& TFITSHDU::GetKeywordValue(const char *keyword)
560 {
561  HDURecord *rec = GetRecord(keyword);
562  if (rec) {
563  return rec->fValue;
564  } else {
565  return *(new TString(""));
566  }
567 }
568 
569 ////////////////////////////////////////////////////////////////////////////////
570 /// Print records.
571 
573 {
574  for (int i = 0; i < fNRecords; i++) {
575  if (fRecords[i].fComment.Length() > 0) {
576  printf("%-10s = %s / %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data(), fRecords[i].fComment.Data());
577  } else {
578  printf("%-10s = %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data());
579  }
580  }
581 }
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// Print HDU's parent file's metadata.
585 
586 void TFITSHDU::PrintFileMetadata(const Option_t *opt) const
587 {
588  fitsfile *fp=0;
589  int status = 0;
590  char errdescr[FLEN_STATUS+1];
591  int hducount, extnum;
592  int hdutype = IMAGE_HDU;
593  const char *exttype;
594  char extname[FLEN_CARD]="PRIMARY"; //room enough
595  int verbose = (opt[0] ? 1 : 0);
596 
597  // Open file with no filters: current HDU will be the primary one.
598  fits_open_file(&fp, fBaseFilePath.Data(), READONLY, &status);
599  if (status) goto ERR;
600 
601  // Read HDU count
602  fits_get_num_hdus(fp, &hducount, &status);
603  if (status) goto ERR;
604  printf("Total: %d HDUs\n", hducount);
605 
606  extnum = 0;
607  while(hducount) {
608  // Read HDU type
609  fits_get_hdu_type(fp, &hdutype, &status);
610  if (status) goto ERR;
611 
612  if (hdutype == IMAGE_HDU) {
613  exttype="IMAGE";
614  } else if (hdutype == ASCII_TBL) {
615  exttype="ASCII TABLE";
616  } else {
617  exttype="BINARY TABLE";
618  }
619 
620  //Read HDU header records
621  int nkeys, morekeys;
622  char keyname[FLEN_KEYWORD+1];
623  char keyvalue[FLEN_VALUE+1];
624  char comment[FLEN_COMMENT+1];
625 
626  fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
627  if (status) goto ERR;
628 
629  struct HDURecord *records = new struct HDURecord[nkeys];
630 
631  for (int i = 1; i <= nkeys; i++) {
632  fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
633  if (status) {
634  delete [] records;
635  goto ERR;
636  }
637 
638  records[i-1].fKeyword = keyname;
639  records[i-1].fValue = keyvalue;
640  records[i-1].fComment = comment;
641 
642  if (strcmp(keyname, "EXTNAME") == 0) {
643  //Extension name
644  strlcpy(extname, keyvalue,FLEN_CARD);
645  }
646  }
647 
648  //HDU info
649  printf(" [%d] %s (%s)\n", extnum, exttype, extname);
650 
651  //HDU records
652  if (verbose) {
653  for (int i = 0; i < nkeys; i++) {
654  if (comment[0]) {
655  printf(" %-10s = %s / %s\n", records[i].fKeyword.Data(), records[i].fValue.Data(), records[i].fComment.Data());
656  } else {
657  printf(" %-10s = %s\n", records[i].fKeyword.Data(), records[i].fValue.Data());
658  }
659  }
660  }
661  printf("\n");
662 
663  delete [] records;
664 
665  hducount--;
666  extnum++;
667  if (hducount){
668  fits_movrel_hdu(fp, 1, &hdutype, &status);
669  if (status) goto ERR;
670  }
671  }
672 
673  // Close file
674  fits_close_file(fp, &status);
675  return;
676 
677 ERR:
678  fits_get_errstatus(status, errdescr);
679  Warning("PrintFileMetadata", "error opening FITS file. Details: %s", errdescr);
680  status = 0;
681  if (fp) fits_close_file(fp, &status);
682 }
683 
684 ////////////////////////////////////////////////////////////////////////////////
685 /// Print column information
686 
688 {
689  if (fType != kTableHDU) {
690  Warning("PrintColumnInfo", "this is not a table HDU.");
691  return;
692  }
693 
694  for (Int_t i = 0; i < fNColumns; i++) {
695  printf("%-20s : %s\n", fColumnsInfo[i].fName.Data(), (fColumnsInfo[i].fType == kRealNumber) ? "REAL NUMBER" : "STRING");
696  }
697 }
698 
699 ////////////////////////////////////////////////////////////////////////////////
700 /// Print full table contents
701 
703 {
704  int printed_chars;
705 
706  if (fType != kTableHDU) {
707  Warning("PrintColumnInfo", "this is not a table HDU.");
708  return;
709  }
710 
711  // Dump header
712  putchar('\n');
713  printed_chars = 0;
714  for (Int_t col = 0; col < fNColumns; col++) {
715  printed_chars += printf("%-10s| ", fColumnsInfo[col].fName.Data());
716  }
717  putchar('\n');
718  while(printed_chars--) {
719  putchar('-');
720  }
721  putchar('\n');
722 
723  // Dump rows
724  for (Int_t row = 0; row < fNRows; row++) {
725  for (Int_t col = 0; col < fNColumns; col++) {
726  if (fColumnsInfo[col].fType == kString) {
727  printf("%-10s", fCells[col * fNRows + row].fString);
728  } else {
729  printed_chars = printf("%.2lg", fCells[col * fNRows + row].fRealNumber);
730  printed_chars -= 10;
731  while (printed_chars < 0) {
732  putchar(' ');
733  printed_chars++;
734  }
735  }
736 
737  if (col <= fNColumns - 1) printf("| ");
738  }
739  printf("\n");
740  }
741 }
742 
743 ////////////////////////////////////////////////////////////////////////////////
744 /// Print metadata.
745 /// Currently supported options:
746 ///
747 /// - "" : print HDU record data
748 /// - "F" : print FITS file's extension names, numbers and types
749 /// - "F+": print FITS file's extension names and types and their record data
750 /// - "T" : print column information when HDU is a table
751 /// - "T+" : print full table (columns header and rows)
752 
753 void TFITSHDU::Print(const Option_t *opt) const
754 {
755  if ((opt[0] == 'F') || (opt[0] == 'f')) {
756  PrintFileMetadata((opt[1] == '+') ? "+" : "");
757  } else if ((opt[0] == 'T') || (opt[0] == 't')) {
758  if (opt[1] == '+') {
759  PrintFullTable("");
760  } else {
761  PrintColumnInfo("");
762  }
763 
764  } else {
765  PrintHDUMetadata("");
766  }
767 }
768 
769 
770 ////////////////////////////////////////////////////////////////////////////////
771 /// Read image HDU as a displayable image. Return 0 if conversion cannot be done.
772 /// If the HDU seems to be a multilayer image, 'layer' parameter can be used
773 /// to retrieve the specified layer (starting from 0)
774 
776 {
777  if (fType != kImageHDU) {
778  Warning("ReadAsImage", "this is not an image HDU.");
779  return 0;
780  }
781 
782 
783  if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
784  Warning("ReadAsImage", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
785  return 0;
786  }
787 
788  Int_t width, height;
789  UInt_t pixels_per_layer;
790 
791  width = Int_t(fSizes->GetAt(0));
792  height = Int_t(fSizes->GetAt(1));
793 
794  pixels_per_layer = UInt_t(width) * UInt_t(height);
795 
796  if ( ((fSizes->GetSize() == 2) && (layer > 0))
797  || (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
798 
799  Warning("ReadAsImage", "layer out of bounds.");
800  return 0;
801  }
802 
803  // Get the maximum and minimum pixel values in the layer to auto-stretch pixels
804  Double_t maxval = 0, minval = 0;
805  UInt_t i;
806  Double_t pixvalue;
807  Int_t offset = layer * pixels_per_layer;
808 
809  for (i = 0; i < pixels_per_layer; i++) {
810  pixvalue = fPixels->GetAt(offset + i);
811 
812  if (pixvalue > maxval) {
813  maxval = pixvalue;
814  }
815 
816  if ((i == 0) || (pixvalue < minval)) {
817  minval = pixvalue;
818  }
819  }
820 
821  //Build the image stretching pixels into a range from 0.0 to 255.0
822  //TImage *im = new TImage(width, height);
823  TImage *im = TImage::Create();
824  if (!im) return 0;
825  TArrayD *layer_pixels = new TArrayD(pixels_per_layer);
826 
827 
828  if (maxval == minval) {
829  //plain image
830  for (i = 0; i < pixels_per_layer; i++) {
831  layer_pixels->SetAt(255.0, i);
832  }
833  } else {
834  Double_t factor = 255.0 / (maxval-minval);
835  for (i = 0; i < pixels_per_layer; i++) {
836  pixvalue = fPixels->GetAt(offset + i);
837  layer_pixels->SetAt(factor * (pixvalue-minval), i) ;
838  }
839  }
840 
841  if (pal == 0) {
842  // Default palette: grayscale palette
843  pal = new TImagePalette(256);
844  for (i = 0; i < 256; i++) {
845  pal->fPoints[i] = ((Double_t)i)/255.0;
846  pal->fColorAlpha[i] = 0xffff;
847  pal->fColorBlue[i] = pal->fColorGreen[i] = pal->fColorRed[i] = i << 8;
848  }
849  pal->fPoints[0] = 0;
850  pal->fPoints[255] = 1.0;
851 
852  im->SetImage(*layer_pixels, UInt_t(width), pal);
853 
854  delete pal;
855 
856  } else {
857  im->SetImage(*layer_pixels, UInt_t(width), pal);
858  }
859 
860  delete layer_pixels;
861 
862  return im;
863 }
864 
865 ////////////////////////////////////////////////////////////////////////////////
866 /// If the HDU is an image, draw the first layer of the primary array
867 /// To set a title to the canvas, pass it in "opt"
868 
870 {
871  if (fType != kImageHDU) {
872  Warning("Draw", "cannot draw. This is not an image HDU.");
873  return;
874  }
875 
876  TImage *im = ReadAsImage(0, 0);
877  if (im) {
878  Int_t width = Int_t(fSizes->GetAt(0));
879  Int_t height = Int_t(fSizes->GetAt(1));
880  TString cname, ctitle;
881  cname.Form("%sHDU", this->GetName());
882  ctitle.Form("%d x %d", width, height);
883  new TCanvas(cname, ctitle, width, height);
884  im->Draw();
885  }
886 }
887 
888 
889 ////////////////////////////////////////////////////////////////////////////////
890 /// Read image HDU as a matrix. Return 0 if conversion cannot be done
891 /// If the HDU seems to be a multilayer image, 'layer' parameter can be used
892 /// to retrieve the specified layer (starting from 0) in matrix form.
893 /// Options (value of 'opt'):
894 /// "S": stretch pixel values to a range from 0.0 to 1.0
895 
897 {
898  if (fType != kImageHDU) {
899  Warning("ReadAsMatrix", "this is not an image HDU.");
900  return 0;
901  }
902 
903  if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
904  Warning("ReadAsMatrix", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
905  return 0;
906  }
907 
908 
909  if ( ((fSizes->GetSize() == 2) && (layer > 0))
910  || (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
911  Warning("ReadAsMatrix", "layer out of bounds.");
912  return 0;
913  }
914 
915  Int_t width, height;
916  UInt_t pixels_per_layer;
917  Int_t offset;
918  UInt_t i;
919  TMatrixD *mat=0;
920 
921  width = Int_t(fSizes->GetAt(0));
922  height = Int_t(fSizes->GetAt(1));
923 
924  pixels_per_layer = UInt_t(width) * UInt_t(height);
925  offset = layer * pixels_per_layer;
926 
927  double *layer_pixels = new double[pixels_per_layer];
928 
929  if ((opt[0] == 'S') || (opt[0] == 's')) {
930  //Stretch
931  // Get the maximum and minimum pixel values in the layer to auto-stretch pixels
932  Double_t factor, maxval=0, minval=0;
933  Double_t pixvalue;
934  for (i = 0; i < pixels_per_layer; i++) {
935  pixvalue = fPixels->GetAt(offset + i);
936 
937  if (pixvalue > maxval) {
938  maxval = pixvalue;
939  }
940 
941  if ((i == 0) || (pixvalue < minval)) {
942  minval = pixvalue;
943  }
944  }
945 
946  if (maxval == minval) {
947  //plain image
948  for (i = 0; i < pixels_per_layer; i++) {
949  layer_pixels[i] = 1.0;
950  }
951  } else {
952  factor = 1.0 / (maxval-minval);
953  mat = new TMatrixD(height, width);
954 
955  for (i = 0; i < pixels_per_layer; i++) {
956  layer_pixels[i] = factor * (fPixels->GetAt(offset + i) - minval);
957  }
958  }
959 
960  } else {
961  //No stretching
962  mat = new TMatrixD(height, width);
963 
964  for (i = 0; i < pixels_per_layer; i++) {
965  layer_pixels[i] = fPixels->GetAt(offset + i);
966  }
967  }
968 
969  if (mat) {
970  // mat->Use(height, width, layer_pixels);
971  memcpy(mat->GetMatrixArray(), layer_pixels, pixels_per_layer*sizeof(double));
972  }
973 
974  delete [] layer_pixels;
975  return mat;
976 }
977 
978 
979 ////////////////////////////////////////////////////////////////////////////////
980 /// Read image HDU as a histogram. Return 0 if conversion cannot be done.
981 /// The returned object can be TH1D, TH2D or TH3D depending on data dimensionality.
982 /// Please, check condition (returnedValue->IsA() == TH*D::Class()) to
983 /// determine the object class.
984 ///
985 /// NOTE: do not confuse with image histogram! This function interprets
986 /// the array as a histogram. It does not compute the histogram of pixel
987 /// values of an image! Here "pixels" are interpreted as number of entries.
988 
990 {
991  if (fType != kImageHDU) {
992  Warning("ReadAsHistogram", "this is not an image HDU.");
993  return 0;
994  }
995 
996  TH1 *result=0;
997 
998  if ((fSizes->GetSize() != 1) && (fSizes->GetSize() != 2) && (fSizes->GetSize() != 3)) {
999  Warning("ReadAsHistogram", "could not convert image HDU to histogram because it has %d dimensions.", fSizes->GetSize());
1000  return 0;
1001  }
1002 
1003  if (fSizes->GetSize() == 1) {
1004  //1D
1005  UInt_t Nx = UInt_t(fSizes->GetAt(0));
1006  UInt_t x;
1007 
1008  TH1D *h = new TH1D("", "", Int_t(Nx), 0, Nx-1);
1009 
1010  for (x = 0; x < Nx; x++) {
1011  Int_t nentries = Int_t(fPixels->GetAt(x));
1012  if (nentries < 0) nentries = 0; //Crop negative values
1013  h->Fill(x, nentries);
1014  }
1015 
1016  result = h;
1017 
1018  } else if (fSizes->GetSize() == 2) {
1019  //2D
1020  UInt_t Nx = UInt_t(fSizes->GetAt(0));
1021  UInt_t Ny = UInt_t(fSizes->GetAt(1));
1022  UInt_t x,y;
1023 
1024  TH2D *h = new TH2D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1);
1025 
1026  for (y = 0; y < Ny; y++) {
1027  UInt_t offset = y * Nx;
1028  for (x = 0; x < Nx; x++) {
1029  Int_t nentries = Int_t(fPixels->GetAt(offset + x));
1030  if (nentries < 0) nentries = 0; //Crop negative values
1031  h->Fill(x,y, nentries);
1032  }
1033  }
1034 
1035  result = h;
1036 
1037  } else if (fSizes->GetSize() == 3) {
1038  //3D
1039  UInt_t Nx = UInt_t(fSizes->GetAt(0));
1040  UInt_t Ny = UInt_t(fSizes->GetAt(1));
1041  UInt_t Nz = UInt_t(fSizes->GetAt(2));
1042  UInt_t x,y,z;
1043 
1044  TH3D *h = new TH3D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1, Int_t(Nz), 0, Nz-1);
1045 
1046 
1047  for (z = 0; z < Nz; z++) {
1048  UInt_t offset1 = z * Nx * Ny;
1049  for (y = 0; y < Ny; y++) {
1050  UInt_t offset2 = y * Nx;
1051  for (x = 0; x < Nx; x++) {
1052  Int_t nentries = Int_t(fPixels->GetAt(offset1 + offset2 + x));
1053  if (nentries < 0) nentries = 0; //Crop negative values
1054  h->Fill(x, y, z, nentries);
1055  }
1056  }
1057  }
1058 
1059  result = h;
1060  }
1061 
1062  return result;
1063 }
1064 
1065 ////////////////////////////////////////////////////////////////////////////////
1066 /// Get a row from the image HDU when it's a 2D array.
1067 
1069 {
1070  if (fType != kImageHDU) {
1071  Warning("GetArrayRow", "this is not an image HDU.");
1072  return 0;
1073  }
1074 
1075  if (fSizes->GetSize() != 2) {
1076  Warning("GetArrayRow", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
1077  return 0;
1078  }
1079 
1080  UInt_t i, offset, W,H;
1081 
1082  W = UInt_t(fSizes->GetAt(0));
1083  H = UInt_t(fSizes->GetAt(1));
1084 
1085 
1086  if (row >= H) {
1087  Warning("GetArrayRow", "index out of bounds.");
1088  return 0;
1089  }
1090 
1091  offset = W * row;
1092  double *v = new double[W];
1093 
1094  for (i = 0; i < W; i++) {
1095  v[i] = fPixels->GetAt(offset+i);
1096  }
1097 
1098  TVectorD *vec = new TVectorD(W, v);
1099 
1100  delete [] v;
1101 
1102  return vec;
1103 }
1104 
1105 ////////////////////////////////////////////////////////////////////////////////
1106 /// Get a column from the image HDU when it's a 2D array.
1107 
1109 {
1110  if (fType != kImageHDU) {
1111  Warning("GetArrayColumn", "this is not an image HDU.");
1112  return 0;
1113  }
1114 
1115  if (fSizes->GetSize() != 2) {
1116  Warning("GetArrayColumn", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
1117  return 0;
1118  }
1119 
1120  UInt_t i, W,H;
1121 
1122  W = UInt_t(fSizes->GetAt(0));
1123  H = UInt_t(fSizes->GetAt(1));
1124 
1125 
1126  if (col >= W) {
1127  Warning("GetArrayColumn", "index out of bounds.");
1128  return 0;
1129  }
1130 
1131  double *v = new double[H];
1132 
1133  for (i = 0; i < H; i++) {
1134  v[i] = fPixels->GetAt(W*i+col);
1135  }
1136 
1137  TVectorD *vec = new TVectorD(H, v);
1138 
1139  delete [] v;
1140 
1141  return vec;
1142 }
1143 
1144 
1145 ////////////////////////////////////////////////////////////////////////////////
1146 ///Get column number given its name
1147 
1148 Int_t TFITSHDU::GetColumnNumber(const char *colname)
1149 {
1150  Int_t colnum;
1151  for (colnum = 0; colnum < fNColumns; colnum++) {
1152  if (fColumnsInfo[colnum].fName == colname) {
1153  return colnum;
1154  }
1155  }
1156  return -1;
1157 }
1158 
1159 ////////////////////////////////////////////////////////////////////////////////
1160 /// Get a string-typed column from a table HDU given its column index (>=0).
1161 
1163 {
1164  if (fType != kTableHDU) {
1165  Warning("GetTabStringColumn", "this is not a table HDU.");
1166  return 0;
1167  }
1168 
1169  if ((colnum < 0) || (colnum >= fNColumns)) {
1170  Warning("GetTabStringColumn", "column index out of bounds.");
1171  return 0;
1172  }
1173 
1174  if (fColumnsInfo[colnum].fType != kString) {
1175  Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
1176  return 0;
1177  }
1178 
1179  Int_t offset = colnum * fNRows;
1180 
1181  TObjArray *res = new TObjArray();
1182  for (Int_t row = 0; row < fNRows; row++) {
1183  res->Add(new TObjString(fCells[offset + row].fString));
1184  }
1185 
1186  return res;
1187 }
1188 
1189 ////////////////////////////////////////////////////////////////////////////////
1190 /// Get a string-typed column from a table HDU given its name
1191 
1193 {
1194  if (fType != kTableHDU) {
1195  Warning("GetTabStringColumn", "this is not a table HDU.");
1196  return 0;
1197  }
1198 
1199 
1200  Int_t colnum = GetColumnNumber(colname);
1201 
1202  if (colnum == -1) {
1203  Warning("GetTabStringColumn", "column not found.");
1204  return 0;
1205  }
1206 
1207  if (fColumnsInfo[colnum].fType != kString) {
1208  Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
1209  return 0;
1210  }
1211 
1212  Int_t offset = colnum * fNRows;
1213 
1214  TObjArray *res = new TObjArray();
1215  for (Int_t row = 0; row < fNRows; row++) {
1216  res->Add(new TObjString(fCells[offset + row].fString));
1217  }
1218 
1219  return res;
1220 }
1221 
1222 ////////////////////////////////////////////////////////////////////////////////
1223 /// Get a real number-typed column from a table HDU given its column index (>=0).
1224 
1226 {
1227  if (fType != kTableHDU) {
1228  Warning("GetTabRealVectorColumn", "this is not a table HDU.");
1229  return 0;
1230  }
1231 
1232  if ((colnum < 0) || (colnum >= fNColumns)) {
1233  Warning("GetTabRealVectorColumn", "column index out of bounds.");
1234  return 0;
1235  }
1236 
1237  if (fColumnsInfo[colnum].fType != kRealNumber) {
1238  Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
1239  return 0;
1240  } else if (fColumnsInfo[colnum].fDim > 1) {
1241  Warning("GetTabRealVectorColumn", "attempting to read a column whose cells have embedded vectors, not real scalars. Use GetTabRealVectorCells() instead.");
1242  return 0;
1243  }
1244 
1245  Int_t offset = colnum * fNRows;
1246 
1247  Double_t *arr = new Double_t[fNRows];
1248 
1249  for (Int_t row = 0; row < fNRows; row++) {
1250  arr[row] = fCells[offset + row].fRealNumber;
1251  }
1252 
1253  TVectorD *res = new TVectorD();
1254  res->Use(fNRows, arr);
1255 
1256  return res;
1257 }
1258 
1259 ////////////////////////////////////////////////////////////////////////////////
1260 /// Get a real number-typed column from a table HDU given its name
1261 
1263 {
1264  if (fType != kTableHDU) {
1265  Warning("GetTabRealVectorColumn", "this is not a table HDU.");
1266  return 0;
1267  }
1268 
1269  Int_t colnum = GetColumnNumber(colname);
1270 
1271  if (colnum == -1) {
1272  Warning("GetTabRealVectorColumn", "column not found.");
1273  return 0;
1274  }
1275 
1276  if (fColumnsInfo[colnum].fType != kRealNumber) {
1277  Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
1278  return 0;
1279  } else if (fColumnsInfo[colnum].fDim > 1) {
1280  Warning("GetTabRealVectorColumn", "attempting to read a column whose cells have embedded vectors, not real scalars. Use GetTabRealVectorCells() instead.");
1281  return 0;
1282  }
1283 
1284  Int_t offset = colnum * fNRows;
1285 
1286  Double_t *arr = new Double_t[fNRows];
1287 
1288  for (Int_t row = 0; row < fNRows; row++) {
1289  arr[row] = fCells[offset + row].fRealNumber;
1290  }
1291 
1292  TVectorD *res = new TVectorD();
1293  res->Use(fNRows, arr);
1294 
1295  return res;
1296 }
1297 
1298 ////////////////////////////////////////////////////////////////////////////////
1299 /// Change to another HDU given by "filter".
1300 /// The parameter "filter" will be appended to the
1301 /// FITS file's base path. For example:
1302 /// hduObject.Change("[EVENTS][TIME > 5]");
1303 /// Please, see documentation of TFITSHDU(const char *filepath_with_filter) constructor
1304 /// for further information.
1305 
1306 Bool_t TFITSHDU::Change(const char *filter)
1307 {
1308  TString tmppath;
1309  tmppath.Form("%s%s", fBaseFilePath.Data(), filter);
1310 
1312  _initialize_me();
1313 
1314  if (kFALSE == LoadHDU(tmppath)) {
1315  //Failed! Restore previous hdu
1316  Warning("Change", "error changing HDU. Restoring the previous one...");
1317 
1319  _initialize_me();
1320 
1321  if (kFALSE == LoadHDU(fFilePath)) {
1322  Warning("Change", "could not restore previous HDU. This object is no longer reliable!!");
1323  }
1324  return kFALSE;
1325  }
1326 
1327  //Set new full path
1328  fFilePath = tmppath;
1329  return kTRUE;
1330 }
1331 
1332 ////////////////////////////////////////////////////////////////////////////////
1333 /// Change to another HDU given by extension_number
1334 
1335 Bool_t TFITSHDU::Change(Int_t extension_number)
1336 {
1337  TString tmppath;
1338  tmppath.Form("[%d]", extension_number);
1339 
1340  return Change(tmppath.Data());
1341 }
1342 
1343 ////////////////////////////////////////////////////////////////////////////////
1344 /// Get a collection of real vectors embedded in cells along a given column from a table HDU. colnum >= 0.
1345 
1347 {
1348  if (fType != kTableHDU) {
1349  Warning("GetTabRealVectorCells", "this is not a table HDU.");
1350  return 0;
1351  }
1352 
1353  if ((colnum < 0) || (colnum >= fNColumns)) {
1354  Warning("GetTabRealVectorCells", "column index out of bounds.");
1355  return 0;
1356  }
1357 
1358  if (fColumnsInfo[colnum].fType != kRealNumber) {
1359  Warning("GetTabRealVectorCells", "attempting to read a column that is not of type 'kRealNumber'.");
1360  return 0;
1361  }
1362 
1363  Int_t offset = colnum * fNRows;
1364 
1365  TObjArray *res = new TObjArray();
1366  Int_t dim = fColumnsInfo[colnum].fDim;
1367 
1368  for (Int_t row = 0; row < fNRows; row++) {
1369  TVectorD *v = new TVectorD();
1370  v->Use(dim, fCells[offset + row].fRealVector);
1371  res->Add(v);
1372  }
1373 
1374  //Make the collection to own the allocated TVectorD objects, so when
1375  //destroying the collection, the vectors will be destroyed too.
1376  res->SetOwner(kTRUE);
1377 
1378  return res;
1379 }
1380 
1381 ////////////////////////////////////////////////////////////////////////////////
1382 /// Get a collection of real vectors embedded in cells along a given column from a table HDU by name
1383 
1385 {
1386  if (fType != kTableHDU) {
1387  Warning("GetTabRealVectorCells", "this is not a table HDU.");
1388  return 0;
1389  }
1390 
1391  Int_t colnum = GetColumnNumber(colname);
1392 
1393  if (colnum == -1) {
1394  Warning("GetTabRealVectorCells", "column not found.");
1395  return 0;
1396  }
1397 
1398  return GetTabRealVectorCells(colnum);
1399 }
1400 
1401 ////////////////////////////////////////////////////////////////////////////////
1402 /// Get a real vector embedded in a cell given by (row>=0, column>=0)
1403 
1405 {
1406  if (fType != kTableHDU) {
1407  Warning("GetTabRealVectorCell", "this is not a table HDU.");
1408  return 0;
1409  }
1410 
1411  if ((colnum < 0) || (colnum >= fNColumns)) {
1412  Warning("GetTabRealVectorCell", "column index out of bounds.");
1413  return 0;
1414  }
1415 
1416  if ((rownum < 0) || (rownum >= fNRows)) {
1417  Warning("GetTabRealVectorCell", "row index out of bounds.");
1418  return 0;
1419  }
1420 
1421  if (fColumnsInfo[colnum].fType != kRealNumber) {
1422  Warning("GetTabRealVectorCell", "attempting to read a column that is not of type 'kRealNumber'.");
1423  return 0;
1424  }
1425 
1426  TVectorD *v = new TVectorD();
1427  v->Use(fColumnsInfo[colnum].fDim, fCells[(colnum * fNRows) + rownum].fRealVector);
1428 
1429  return v;
1430 }
1431 
1432 ////////////////////////////////////////////////////////////////////////////////
1433 /// Get a real vector embedded in a cell given by (row>=0, column name)
1434 
1435 TVectorD *TFITSHDU::GetTabRealVectorCell(Int_t rownum, const char *colname)
1436 {
1437  if (fType != kTableHDU) {
1438  Warning("GetTabRealVectorCell", "this is not a table HDU.");
1439  return 0;
1440  }
1441 
1442  Int_t colnum = GetColumnNumber(colname);
1443 
1444  if (colnum == -1) {
1445  Warning("GetTabRealVectorCell", "column not found.");
1446  return 0;
1447  }
1448 
1449  return GetTabRealVectorCell(rownum, colnum);
1450 }
1451 
1452 
1453 ////////////////////////////////////////////////////////////////////////////////
1454 /// Get the name of a column given its index (column>=0).
1455 /// In case of error the column name is "".
1456 
1458 {
1459  static TString noName;
1460  if (fType != kTableHDU) {
1461  Error("GetColumnName", "this is not a table HDU.");
1462  return noName;
1463  }
1464 
1465  if ((colnum < 0) || (colnum >= fNColumns)) {
1466  Error("GetColumnName", "column index out of bounds.");
1467  return noName;
1468  }
1469  return fColumnsInfo[colnum].fName;
1470 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
TString fKeyword
Definition: TFITS.h:53
Double_t * fPoints
[fNumPoints] value of each anchor point [0..1]
Definition: TAttImage.h:37
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3251
Double_t * fRealVector
Definition: TFITS.h:68
An array of TObjects.
Definition: TObjArray.h:37
TString fValue
Definition: TFITS.h:54
void SetAt(Double_t v, Int_t i)
Definition: TArrayD.h:51
Int_t fNRows
Number of rows (when fType == kTableHDU)
Definition: TFITS.h:83
TString fExtensionName
Extension Name.
Definition: TFITS.h:77
Collectable string class.
Definition: TObjString.h:28
TFITSHDU(const char *filepath_with_filter)
TFITSHDU constructor from file path with HDU selection filter.
Definition: TFITS.cxx:118
const char Option_t
Definition: RtypesCore.h:62
void Draw(Option_t *opt="")
If the HDU is an image, draw the first layer of the primary array To set a title to the canvas...
Definition: TFITS.cxx:869
Int_t fNRecords
Number of records.
Definition: TFITS.h:75
virtual const Element * GetMatrixArray() const
Definition: TMatrixT.h:222
TArrayI * fSizes
Image sizes in each dimension (when fType == kImageHDU)
Definition: TFITS.h:79
const Ssiz_t kNPOS
Definition: RtypesCore.h:111
image html pict1_TGaxis_012 png width
Define new text attributes for the label number "labNum".
Definition: TGaxis.cxx:2551
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TMatrixD * ReadAsMatrix(Int_t layer=0, Option_t *opt="")
Read image HDU as a matrix.
Definition: TFITS.cxx:896
#define H(x, y, z)
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:634
Basic string class.
Definition: TString.h:131
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Int_t fDim
Definition: TFITS.h:61
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
void _release_resources()
Release internal resources.
Definition: TFITS.cxx:177
Bool_t LoadHDU(TString &filepath_filter)
Load HDU from fits file satisfying the specified filter.
Definition: TFITS.cxx:231
virtual void SetImage(const Double_t *, UInt_t, UInt_t, TImagePalette *=0)
Definition: TImage.h:116
Array of integers (32 bits per element).
Definition: TArrayI.h:27
UShort_t * fColorRed
[fNumPoints] red color at each anchor point
Definition: TAttImage.h:38
void Print(const Option_t *opt="") const
Print metadata.
Definition: TFITS.cxx:753
~TFITSHDU()
TFITSHDU destructor.
Definition: TFITS.cxx:169
struct HDURecord * GetRecord(const char *keyword)
Get record by keyword.
Definition: TFITS.cxx:546
Double_t fRealNumber
Definition: TFITS.h:67
TVectorT< Element > & Use(Int_t lwb, Int_t upb, Element *data)
Use the array data to fill the vector lwb..upb].
Definition: TVectorT.cxx:347
Double_t x[n]
Definition: legend1.C:17
FITS file interface class.
Definition: TFITS.h:34
void _initialize_me()
Do some initializations.
Definition: TFITS.cxx:216
void PrintFileMetadata(const Option_t *opt="") const
Print HDU&#39;s parent file&#39;s metadata.
Definition: TFITS.cxx:586
TVectorD * GetArrayRow(UInt_t row)
Get a row from the image HDU when it&#39;s a 2D array.
Definition: TFITS.cxx:1068
enum EHDUTypes fType
HDU type.
Definition: TFITS.h:76
Double_t GetAt(Int_t i) const
Definition: TArrayI.h:45
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:22
TObjArray * GetTabStringColumn(Int_t colnum)
Get a string-typed column from a table HDU given its column index (>=0).
Definition: TFITS.cxx:1162
TH1 * ReadAsHistogram()
Read image HDU as a histogram.
Definition: TFITS.cxx:989
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:22
THist< 3, double, THistStatContent, THistStatUncertainty > TH3D
Definition: THist.hxx:296
Bool_t Change(const char *filter)
Change to another HDU given by "filter".
Definition: TFITS.cxx:1306
enum EColumnTypes fType
Definition: TFITS.h:60
void PrintColumnInfo(const Option_t *) const
Print column information.
Definition: TFITS.cxx:687
3-D histogram with a double per channel (see TH1 documentation)}
Definition: TH3.h:304
Int_t GetSize() const
Definition: TArray.h:47
SVector< double, 2 > v
Definition: Dict.h:5
UShort_t * fColorAlpha
[fNumPoints] alpha at each anchor point
Definition: TAttImage.h:41
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2264
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
Ssiz_t Length() const
Definition: TString.h:405
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH3.cxx:286
void PrintHDUMetadata(const Option_t *opt="") const
Print records.
Definition: TFITS.cxx:572
UShort_t * fColorBlue
[fNumPoints] blue color at each anchor point
Definition: TAttImage.h:40
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:610
TString fName
Definition: TNamed.h:32
#define h(i)
Definition: RSha256.hxx:106
const Bool_t kFALSE
Definition: RtypesCore.h:88
int Ssiz_t
Definition: RtypesCore.h:63
The Canvas class.
Definition: TCanvas.h:31
TArrayD * fPixels
Image pixels (when fType == kImageHDU)
Definition: TFITS.h:80
#define ClassImp(name)
Definition: Rtypes.h:359
void SetAt(Double_t v, Int_t i)
Definition: TArrayI.h:51
double Double_t
Definition: RtypesCore.h:55
TString fName
Definition: TFITS.h:59
TVectorD * GetTabRealVectorColumn(Int_t colnum)
Get a real number-typed column from a table HDU given its column index (>=0).
Definition: TFITS.cxx:1225
Char_t * fString
Definition: TFITS.h:66
Int_t fNumber
HDU number (1=PRIMARY)
Definition: TFITS.h:78
TString fComment
Definition: TFITS.h:55
int nentries
Definition: THbookFile.cxx:89
Double_t y[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:56
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:290
A class to define a conversion from pixel values to pixel color.
Definition: TAttImage.h:33
Array of doubles (64 bits per element).
Definition: TArrayD.h:27
TString fBaseFilePath
Path to HDU&#39;s file excluding filter.
Definition: TFITS.h:73
TVectorD * GetTabRealVectorCell(Int_t rownum, Int_t colnum)
Get a real vector embedded in a cell given by (row>=0, column>=0)
Definition: TFITS.cxx:1404
TString fFilePath
Path to HDU&#39;s file including filter.
Definition: TFITS.h:72
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
UShort_t * fColorGreen
[fNumPoints] green color at each anchor point
Definition: TAttImage.h:39
TVectorD * GetArrayColumn(UInt_t col)
Get a column from the image HDU when it&#39;s a 2D array.
Definition: TFITS.cxx:1108
union Cell * fCells
Table cells (when fType == kTableHDU).
Definition: TFITS.h:84
TString & GetKeywordValue(const char *keyword)
Get the value of a given keyword. Return "" if not found.
Definition: TFITS.cxx:559
static void CleanFilePath(const char *filepath_with_filter, TString &dst)
Clean path from possible filter and put the result in &#39;dst&#39;.
Definition: TFITS.cxx:94
struct HDURecord * fRecords
HDU metadata records.
Definition: TFITS.h:74
Double_t GetAt(Int_t i) const
Definition: TArrayD.h:45
void PrintFullTable(const Option_t *) const
Print full table contents.
Definition: TFITS.cxx:702
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:284
void Add(TObject *obj)
Definition: TObjArray.h:73
const TString & GetColumnName(Int_t colnum)
Get the name of a column given its index (column>=0).
Definition: TFITS.cxx:1457
static TImage * Create()
Create an image.
Definition: TImage.cxx:36
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
struct Column * fColumnsInfo
Information about columns (when fType == kTableHDU)
Definition: TFITS.h:81
const Bool_t kTRUE
Definition: RtypesCore.h:87
Int_t fNColumns
Number of columns (when fType == kTableHDU)
Definition: TFITS.h:82
TObjArray * GetTabRealVectorCells(Int_t colnum)
Get a collection of real vectors embedded in cells along a given column from a table HDU...
Definition: TFITS.cxx:1346
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
void Resize(Ssiz_t n)
Resize the string. Truncate or add blanks as necessary.
Definition: TString.cxx:1070
TImage * ReadAsImage(Int_t layer=0, TImagePalette *pal=0)
Read image HDU as a displayable image.
Definition: TFITS.cxx:775
Int_t GetColumnNumber(const char *colname)
Get column number given its name.
Definition: TFITS.cxx:1148
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:291
const char * Data() const
Definition: TString.h:364