NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
NcFITSIO.h
Go to the documentation of this file.
1#ifndef NcFITSIO_h
2#define NcFITSIO_h
3// Copyright(c) 2019 NCFS/IIHE, All Rights Reserved.
4// See cxx source for full Copyright notice.
5
6#include <stdlib.h>
7
8#include "fitsio.h"
9
10#include "TSystem.h"
11#include "TNamed.h"
12#include "TArrayI.h"
13#include "TArrayD.h"
14#include "TH1D.h"
15#include "TH2D.h"
16#include "TH3D.h"
17#include "TASImage.h"
18#include "TVectorD.h"
19#include "TMatrixD.h"
20#include "TObjArray.h"
21#include "TString.h"
22#include "TObjString.h"
23
25
26class NcFITSIO : public TNamed
27{
28 public:
29 // Markers of the various HDU types
31
32 // Markers of the various column data types
34
35 NcFITSIO(const char* name="NcFITSIO",const char* title="FITS data I/O interface"); // Default constructor
36 virtual ~NcFITSIO(); // Default destructor
37 NcFITSIO(const NcFITSIO &q); // Copy constructor
38 virtual TObject* Clone(const char* name="") const; // Make a deep copy and provide its pointer
39
40 // Input file handling
41 Bool_t OpenInputFile(TString specs);
42 Bool_t SelectHDU(TString extname="[0]");
43 Bool_t SelectHDU(Int_t extnumber);
44 void ListHDUHeader() const;
45 void ListFileHeader(Int_t mode=1) const;
46 TString GetKeywordValue(TString keyname,Int_t mode=0);
47 Bool_t IsTable() const;
48 Int_t GetHDUCount() const;
49
50 // Table access methods
51 Int_t GetTableNrows() const;
52 Int_t GetTableNcolumns() const;
53 Int_t GetColumnNumber(TString colname,Int_t mode=0) const;
54 TString GetColumnName(Int_t colnum) const;
55 Int_t GetTableCell(Double_t &val,Int_t row,Int_t col,Int_t layer=1);
56 Int_t GetTableCell(Double_t &val,Int_t row,TString colname,Int_t layer=1,Int_t mode=0);
57 Int_t GetTableCell(TArrayD &arr,Int_t row,Int_t col);
58 Int_t GetTableCell(TArrayD &arr,Int_t row,TString colname,Int_t mode=0);
59 Int_t GetTableCell(TString &str,Int_t row,Int_t col,Int_t layer=1);
60 Int_t GetTableCell(TString &str,Int_t row,TString colname,Int_t layer=1,Int_t mode=0);
61 Int_t GetTableCell(TString* &arr,Int_t row,Int_t col);
62 Int_t GetTableCell(TString* &arr,Int_t row,TString colname,Int_t mode=0);
63 Int_t GetTableCell(TObjArray &arr,Int_t row,Int_t col);
64 Int_t GetTableCell(TObjArray &arr,Int_t row,TString colname,Int_t mode=0);
65 Int_t GetTableColumn(TArrayD &arr,Int_t col,Int_t rstart=1,Int_t rend=0,Int_t layer=1);
66 Int_t GetTableColumn(TArrayD &arr,TString colname,Int_t rstart=1,Int_t rend=0,Int_t layer=1,Int_t mode=0);
67 Int_t GetTableColumn(TString* &arr,Int_t col,Int_t rstart=1,Int_t rend=0,Int_t layer=1);
68 Int_t GetTableColumn(TString* &arr,TString colname,Int_t rstart=1,Int_t rend=0,Int_t layer=1,Int_t mode=0);
69 Int_t GetTableColumn(TObjArray &arr,Int_t col,Int_t rstart=1,Int_t rend=0,Int_t layer=1);
70 Int_t GetTableColumn(TObjArray &arr,TString colname,Int_t rstart=1,Int_t rend=0,Int_t layer=1,Int_t mode=0);
71 void ListTable(Int_t width=-10,Int_t rstart=1,Int_t rend=0,Int_t cstart=1,Int_t cend=0);
72
73 // Image access methods
74 Int_t GetImageDimension(Int_t i=0) const;
75 Int_t GetImageLayer(TASImage &im,Int_t layer=1,Double_t* thres=0,Double_t max=-1);
76 Int_t GetImageLayer(TMatrixD &m,Int_t layer=1,Double_t* thres=0,Double_t max=-1);
77 Int_t GetImageLayer(TH2D &his,Int_t layer=1,Double_t* thres=0,Double_t max=-1);
78 UInt_t GetImageArray(TArrayD &arr,TArrayI ifirst,TArrayI ilast,TArrayI incr);
79 UInt_t GetImageArray(TArrayD &arr,TArrayI ifirst,UInt_t npix);
80
81 protected:
82 TString fFilename; // The (full path) name of the FITS file on the computer system
83 TString fFilenameFilter; // The FITS filename with the HDU selection filter
84 fitsfile* fInput; // Pointer to the FITS input file
85 fitsfile* fOutput; // Pointer to the FITS output file
86 eHDUTypes fType; // The HDU type
87 TString fExtensionName; // The HDU extension Name
88 Int_t fExtensionNumber; // The HDU extension number (0=PRIMARY)
89 Int_t fNkeys; // The number of HDU keywords
90 TString* fKeyNames; // The HDU key names
91 TString* fKeyValues; // The HDU key values
92 TString* fComments; // The HDU (key) comments
93 Int_t fNrows; // The number of table rows
94 Int_t fNcolumns; // The number of table columns
95 TString* fColumnNames; // The names of the table columns
96 eColumnTypes* fColumnTypes; // The types of the table columns
97 Int_t* fColumnLayers; // The number of layers of the table column
98 TArrayI* fSizes; // Image sizes in each dimension
99
100 void Initialize();
101 void Reset();
102 TString StripFilter(TString filename) const;
103 Bool_t LoadHeaderInfo();
104 Int_t LoadLayer(TArrayD &arr,Int_t layer);
105 void ApplyPixelThreshold(TArrayD &arr,Double_t thres);
106 void RescalePixels(TArrayD &arr,Double_t max);
107
108 ClassDef(NcFITSIO,0) // I/O interface for FITS files
109};
110#endif
I/O interface for FITS files.
Definition NcFITSIO.h:27
@ kComplexArray
Definition NcFITSIO.h:33
@ kComplexNumber
Definition NcFITSIO.h:33
@ kRealNumber
Definition NcFITSIO.h:33
@ kLogical
Definition NcFITSIO.h:33
@ kStringArray
Definition NcFITSIO.h:33
@ kString
Definition NcFITSIO.h:33
@ kLogicalArray
Definition NcFITSIO.h:33
@ kRealArray
Definition NcFITSIO.h:33
TString * fComments
Definition NcFITSIO.h:92
void Initialize()
Definition NcFITSIO.cxx:338
Bool_t LoadHeaderInfo()
Definition NcFITSIO.cxx:485
Int_t * fColumnLayers
Definition NcFITSIO.h:97
void ListTable(Int_t width=-10, Int_t rstart=1, Int_t rend=0, Int_t cstart=1, Int_t cend=0)
Int_t GetImageLayer(TASImage &im, Int_t layer=1, Double_t *thres=0, Double_t max=-1)
Int_t fNkeys
Definition NcFITSIO.h:89
Int_t GetTableCell(Double_t &val, Int_t row, Int_t col, Int_t layer=1)
NcFITSIO(const char *name="NcFITSIO", const char *title="FITS data I/O interface")
Definition NcFITSIO.cxx:297
TString * fKeyNames
Definition NcFITSIO.h:90
UInt_t GetImageArray(TArrayD &arr, TArrayI ifirst, TArrayI ilast, TArrayI incr)
Int_t GetHDUCount() const
Definition NcFITSIO.cxx:900
TString * fColumnNames
Definition NcFITSIO.h:95
Int_t GetTableNcolumns() const
Definition NcFITSIO.cxx:951
fitsfile * fOutput
Definition NcFITSIO.h:85
TString * fKeyValues
Definition NcFITSIO.h:91
void ListFileHeader(Int_t mode=1) const
TString GetKeywordValue(TString keyname, Int_t mode=0)
Definition NcFITSIO.cxx:845
void ListHDUHeader() const
eColumnTypes * fColumnTypes
Definition NcFITSIO.h:96
Int_t fNrows
Definition NcFITSIO.h:93
virtual TObject * Clone(const char *name="") const
Int_t fNcolumns
Definition NcFITSIO.h:94
TString fExtensionName
Definition NcFITSIO.h:87
Int_t GetColumnNumber(TString colname, Int_t mode=0) const
Definition NcFITSIO.cxx:962
void RescalePixels(TArrayD &arr, Double_t max)
Int_t LoadLayer(TArrayD &arr, Int_t layer)
TString fFilename
Definition NcFITSIO.h:82
TArrayI * fSizes
Definition NcFITSIO.h:98
Bool_t OpenInputFile(TString specs)
Definition NcFITSIO.cxx:413
Int_t fExtensionNumber
Definition NcFITSIO.h:88
void Reset()
Definition NcFITSIO.cxx:365
Int_t GetTableNrows() const
Definition NcFITSIO.cxx:940
virtual ~NcFITSIO()
Definition NcFITSIO.cxx:308
void ApplyPixelThreshold(TArrayD &arr, Double_t thres)
Int_t GetTableColumn(TArrayD &arr, Int_t col, Int_t rstart=1, Int_t rend=0, Int_t layer=1)
fitsfile * fInput
Definition NcFITSIO.h:84
TString fFilenameFilter
Definition NcFITSIO.h:83
eHDUTypes fType
Definition NcFITSIO.h:86
Bool_t SelectHDU(TString extname="[0]")
Definition NcFITSIO.cxx:777
Bool_t IsTable() const
Definition NcFITSIO.cxx:885
Int_t GetImageDimension(Int_t i=0) const
TString StripFilter(TString filename) const
Definition NcFITSIO.cxx:470
@ kTableHDU
Definition NcFITSIO.h:30
@ kImageHDU
Definition NcFITSIO.h:30
TString GetColumnName(Int_t colnum) const
Definition NcFITSIO.cxx:993