#include "TFITS.h"
#include "TROOT.h"
#include "TImage.h"
#include "TArrayI.h"
#include "TArrayD.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TH3D.h"
#include "TVectorD.h"
#include "TMatrixD.h"
#include "TObjArray.h"
#include "TObjString.h"
#include "TCanvas.h"
#include "fitsio.h"
#include <stdlib.h>
ClassImp(TFITSHDU)
void TFITSHDU::CleanFilePath(const char *filepath_with_filter, TString &dst)
{
dst = filepath_with_filter;
Ssiz_t ndx = dst.Index("[", 1, 0, TString::kExact);
if (ndx != kNPOS) {
dst.Resize(ndx);
}
}
TFITSHDU::TFITSHDU(const char *filepath_with_filter)
{
_initialize_me();
fFilePath = filepath_with_filter;
CleanFilePath(filepath_with_filter, fBaseFilePath);
if (kFALSE == LoadHDU(fFilePath)) {
_release_resources();
throw -1;
}
}
TFITSHDU::TFITSHDU(const char *filepath, Int_t extension_number)
{
_initialize_me();
CleanFilePath(filepath, fBaseFilePath);
fFilePath.Form("%s[%d]", fBaseFilePath.Data(), extension_number);
if (kFALSE == LoadHDU(fFilePath)) {
_release_resources();
throw -1;
}
}
TFITSHDU::TFITSHDU(const char *filepath, const char *extension_name)
{
_initialize_me();
CleanFilePath(filepath, fBaseFilePath);
fFilePath.Form("%s[%s]", fBaseFilePath.Data(), extension_name);
if (kFALSE == LoadHDU(fFilePath)) {
_release_resources();
throw -1;
}
}
TFITSHDU::~TFITSHDU()
{
_release_resources();
}
void TFITSHDU::_release_resources()
{
if (fRecords) delete [] fRecords;
if (fType == kImageHDU) {
if (fSizes) delete fSizes;
if (fPixels) delete fPixels;
} else {
if (fColumnsInfo) {
if (fCells) {
for (Int_t i = 0; i < fNColumns; i++) {
if (fColumnsInfo[i].fType == kString) {
Int_t offset = i * fNRows;
for (Int_t row = 0; row < fNRows; row++) {
delete [] fCells[offset+row].fString;
}
} else if (fColumnsInfo[i].fType == kRealVector) {
Int_t offset = i * fNRows;
for (Int_t row = 0; row < fNRows; row++) {
delete [] fCells[offset+row].fRealVector;
}
}
}
delete [] fCells;
}
delete [] fColumnsInfo;
}
}
}
void TFITSHDU::_initialize_me()
{
fRecords = 0;
fPixels = 0;
fSizes = 0;
fColumnsInfo = 0;
fNColumns = fNRows = 0;
fCells = 0;
}
Bool_t TFITSHDU::LoadHDU(TString& filepath_filter)
{
fitsfile *fp=0;
int status = 0;
char errdescr[FLEN_STATUS+1];
fits_open_file(&fp, filepath_filter.Data(), READONLY, &status);
if (status) goto ERR;
int hdunum;
fits_get_hdu_num(fp, &hdunum);
fNumber = Int_t(hdunum);
int hdutype;
fits_get_hdu_type(fp, &hdutype, &status);
if (status) goto ERR;
fType = (hdutype == IMAGE_HDU) ? kImageHDU : kTableHDU;
int nkeys, morekeys;
char keyname[FLEN_KEYWORD+1];
char keyvalue[FLEN_VALUE+1];
char comment[FLEN_COMMENT+1];
fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
if (status) goto ERR;
fRecords = new struct HDURecord[nkeys];
for (int i = 1; i <= nkeys; i++) {
fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
if (status) goto ERR;
fRecords[i-1].fKeyword = keyname;
fRecords[i-1].fValue = keyvalue;
fRecords[i-1].fComment = comment;
}
fNRecords = Int_t(nkeys);
fExtensionName = "PRIMARY";
for (int i = 0; i < nkeys; i++) {
if (fRecords[i].fKeyword == "EXTNAME") {
fExtensionName = fRecords[i].fValue;
break;
}
}
if (fType == kImageHDU) {
int param_ndims=0;
long *param_dimsizes;
fits_get_img_dim(fp, ¶m_ndims, &status);
if (status) goto ERR;
if (param_ndims > 0) {
param_dimsizes = new long[param_ndims];
fits_get_img_size(fp, param_ndims, param_dimsizes, &status);
if (status) {
delete [] param_dimsizes;
goto ERR;
}
fSizes = new TArrayI(param_ndims);
fSizes = new TArrayI(param_ndims);
for (int i = 0; i < param_ndims; i++) {
fSizes->SetAt(param_dimsizes[i], i);
}
delete [] param_dimsizes;
int anynul;
long *firstpixel = new long[param_ndims];
double nulval = 0;
long npixels = 1;
for (int i = 0; i < param_ndims; i++) {
npixels *= (long) fSizes->GetAt(i);
firstpixel[i] = 1;
}
double *pixels = new double[npixels];
fits_read_pix(fp, TDOUBLE, firstpixel, npixels,
(void *) &nulval, (void *) pixels, &anynul, &status);
if (status) {
delete [] firstpixel;
delete [] pixels;
goto ERR;
}
fPixels = new TArrayD(npixels, pixels);
delete [] firstpixel;
delete [] pixels;
} else {
fSizes = new TArrayI();
fPixels = new TArrayD();
}
} else {
long table_rows;
int table_cols;
fits_get_num_rows(fp, &table_rows, &status);
if (status) goto ERR;
fNRows = Int_t(table_rows);
fits_get_num_cols(fp, &table_cols, &status);
if (status) goto ERR;
fNColumns = Int_t(table_cols);
fColumnsInfo = new struct Column[table_cols];
char colname[80];
int colnum;
fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
while (status == COL_NOT_UNIQUE)
{
fColumnsInfo[colnum-1].fName = colname;
fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
}
if (status != COL_NOT_FOUND) goto ERR;
status = 0;
fCells = new union Cell [table_rows * table_cols];
int typecode;
long repeat, width;
Int_t cellindex;
for (colnum = 0, cellindex = 0; colnum < fNColumns; colnum++) {
fits_get_coltype(fp, colnum+1, &typecode, &repeat, &width, &status);
if (status) goto ERR;
if ((typecode == TDOUBLE) || (typecode == TSHORT) || (typecode == TLONG)
|| (typecode == TFLOAT) || (typecode == TLOGICAL) || (typecode == TBIT)
|| (typecode == TBYTE) || (typecode == TSTRING)) {
fColumnsInfo[colnum].fType = (typecode == TSTRING) ? kString : kRealNumber;
if (typecode == TSTRING) {
int dispwidth=0;
fits_get_col_display_width(fp, colnum+1, &dispwidth, &status);
if (status) goto ERR;
char *nulval = (char*) "";
int anynul=0;
char **array;
if (dispwidth <= 0) {
dispwidth = 1;
}
array = new char* [table_rows];
for (long row = 0; row < table_rows; row++) {
array[row] = new char[dispwidth+1];
}
if (repeat > 0) {
fits_read_col(fp, TSTRING, colnum+1, 1, 1, table_rows, nulval, array, &anynul, &status);
if (status) {
for (long row = 0; row < table_rows; row++) {
delete [] array[row];
}
delete [] array;
goto ERR;
}
} else {
for (long row = 0; row < table_rows; row++) {
strlcpy(array[row], "-",dispwidth+1);
}
}
for (long row = 0; row < table_rows; row++) {
fCells[cellindex++].fString = array[row];
}
delete [] array;
} else {
double nulval = 0;
int anynul=0;
fColumnsInfo[colnum].fDim = (Int_t) repeat;
double *array;
if (repeat > 0) {
array = new double [table_rows * repeat];
fits_read_col(fp, TDOUBLE, colnum+1, 1, 1, table_rows * repeat, &nulval, array, &anynul, &status);
if (status) {
delete [] array;
goto ERR;
}
} else {
array = new double [table_rows];
for (long row = 0; row < table_rows; row++) {
array[row] = 0.0;
}
}
if (repeat == 1) {
for (long row = 0; row < table_rows; row++) {
fCells[cellindex++].fRealNumber = array[row];
}
} else if (repeat > 1) {
for (long row = 0; row < table_rows; row++) {
double *vec = new double [repeat];
long offset = row * repeat;
for (long component = 0; component < repeat; component++) {
vec[component] = array[offset++];
}
fCells[cellindex++].fRealVector = vec;
}
}
delete [] array;
}
} else {
Warning("LoadHDU", "error opening FITS file. Column type %d is currently not supported", typecode);
}
}
if (hdutype == ASCII_TBL) {
} else {
}
}
fits_close_file(fp, &status);
return kTRUE;
ERR:
fits_get_errstatus(status, errdescr);
Warning("LoadHDU", "error opening FITS file. Details: %s", errdescr);
status = 0;
if (fp) fits_close_file(fp, &status);
return kFALSE;
}
struct TFITSHDU::HDURecord* TFITSHDU::GetRecord(const char *keyword)
{
for (int i = 0; i < fNRecords; i++) {
if (fRecords[i].fKeyword == keyword) {
return &fRecords[i];
}
}
return 0;
}
TString& TFITSHDU::GetKeywordValue(const char *keyword)
{
HDURecord *rec = GetRecord(keyword);
if (rec) {
return rec->fValue;
} else {
return *(new TString(""));
}
}
void TFITSHDU::PrintHDUMetadata(const Option_t *) const
{
for (int i = 0; i < fNRecords; i++) {
if (fRecords[i].fComment.Length() > 0) {
printf("%-10s = %s / %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data(), fRecords[i].fComment.Data());
} else {
printf("%-10s = %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data());
}
}
}
void TFITSHDU::PrintFileMetadata(const Option_t *opt) const
{
fitsfile *fp=0;
int status = 0;
char errdescr[FLEN_STATUS+1];
int hducount, extnum;
int hdutype = IMAGE_HDU;
const char *exttype;
char extname[FLEN_CARD]="PRIMARY";
int verbose = (opt[0] ? 1 : 0);
fits_open_file(&fp, fBaseFilePath.Data(), READONLY, &status);
if (status) goto ERR;
fits_get_num_hdus(fp, &hducount, &status);
if (status) goto ERR;
printf("Total: %d HDUs\n", hducount);
extnum = 0;
while(hducount) {
fits_get_hdu_type(fp, &hdutype, &status);
if (status) goto ERR;
if (hdutype == IMAGE_HDU) {
exttype="IMAGE";
} else if (hdutype == ASCII_TBL) {
exttype="ASCII TABLE";
} else {
exttype="BINARY TABLE";
}
int nkeys, morekeys;
char keyname[FLEN_KEYWORD+1];
char keyvalue[FLEN_VALUE+1];
char comment[FLEN_COMMENT+1];
fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
if (status) goto ERR;
struct HDURecord *records = new struct HDURecord[nkeys];
for (int i = 1; i <= nkeys; i++) {
fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
if (status) {
delete [] records;
goto ERR;
}
records[i-1].fKeyword = keyname;
records[i-1].fValue = keyvalue;
records[i-1].fComment = comment;
if (strcmp(keyname, "EXTNAME") == 0) {
strlcpy(extname, keyvalue,FLEN_CARD);
}
}
printf(" [%d] %s (%s)\n", extnum, exttype, extname);
if (verbose) {
for (int i = 0; i < nkeys; i++) {
if (comment[0]) {
printf(" %-10s = %s / %s\n", records[i].fKeyword.Data(), records[i].fValue.Data(), records[i].fComment.Data());
} else {
printf(" %-10s = %s\n", records[i].fKeyword.Data(), records[i].fValue.Data());
}
}
}
printf("\n");
delete [] records;
hducount--;
extnum++;
if (hducount){
fits_movrel_hdu(fp, 1, &hdutype, &status);
if (status) goto ERR;
}
}
fits_close_file(fp, &status);
return;
ERR:
fits_get_errstatus(status, errdescr);
Warning("PrintFileMetadata", "error opening FITS file. Details: %s", errdescr);
status = 0;
if (fp) fits_close_file(fp, &status);
}
void TFITSHDU::PrintColumnInfo(const Option_t *) const
{
if (fType != kTableHDU) {
Warning("PrintColumnInfo", "this is not a table HDU.");
return;
}
for (Int_t i = 0; i < fNColumns; i++) {
printf("%-20s : %s\n", fColumnsInfo[i].fName.Data(), (fColumnsInfo[i].fType == kRealNumber) ? "REAL NUMBER" : "STRING");
}
}
void TFITSHDU::PrintFullTable(const Option_t *) const
{
int printed_chars;
if (fType != kTableHDU) {
Warning("PrintColumnInfo", "this is not a table HDU.");
return;
}
putchar('\n');
printed_chars = 0;
for (Int_t col = 0; col < fNColumns; col++) {
printed_chars += printf("%-10s| ", fColumnsInfo[col].fName.Data());
}
putchar('\n');
while(printed_chars--) {
putchar('-');
}
putchar('\n');
for (Int_t row = 0; row < fNRows; row++) {
for (Int_t col = 0; col < fNColumns; col++) {
if (fColumnsInfo[col].fType == kString) {
printf("%-10s", fCells[col * fNRows + row].fString);
} else {
printed_chars = printf("%.2lg", fCells[col * fNRows + row].fRealNumber);
printed_chars -= 10;
while (printed_chars < 0) {
putchar(' ');
printed_chars++;
}
}
if (col <= fNColumns - 1) printf("| ");
}
printf("\n");
}
}
void TFITSHDU::Print(const Option_t *opt) const
{
if ((opt[0] == 'F') || (opt[0] == 'f')) {
PrintFileMetadata((opt[1] == '+') ? "+" : "");
} else if ((opt[0] == 'T') || (opt[0] == 't')) {
if (opt[1] == '+') {
PrintFullTable("");
} else {
PrintColumnInfo("");
}
} else {
PrintHDUMetadata("");
}
}
TImage *TFITSHDU::ReadAsImage(Int_t layer, TImagePalette *pal)
{
if (fType != kImageHDU) {
Warning("ReadAsImage", "this is not an image HDU.");
return 0;
}
if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
Warning("ReadAsImage", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
return 0;
}
Int_t width, height;
UInt_t pixels_per_layer;
width = Int_t(fSizes->GetAt(0));
height = Int_t(fSizes->GetAt(1));
pixels_per_layer = UInt_t(width) * UInt_t(height);
if ( ((fSizes->GetSize() == 2) && (layer > 0))
|| (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
Warning("ReadAsImage", "layer out of bounds.");
return 0;
}
Double_t maxval = 0, minval = 0;
UInt_t i;
Double_t pixvalue;
Int_t offset = layer * pixels_per_layer;
for (i = 0; i < pixels_per_layer; i++) {
pixvalue = fPixels->GetAt(offset + i);
if (pixvalue > maxval) {
maxval = pixvalue;
}
if ((i == 0) || (pixvalue < minval)) {
minval = pixvalue;
}
}
TImage *im = TImage::Create();
if (!im) return 0;
TArrayD *layer_pixels = new TArrayD(pixels_per_layer);
if (maxval == minval) {
for (i = 0; i < pixels_per_layer; i++) {
layer_pixels->SetAt(255.0, i);
}
} else {
Double_t factor = 255.0 / (maxval-minval);
for (i = 0; i < pixels_per_layer; i++) {
pixvalue = fPixels->GetAt(offset + i);
layer_pixels->SetAt(factor * (pixvalue-minval), i) ;
}
}
if (pal == 0) {
pal = new TImagePalette(256);
for (i = 0; i < 256; i++) {
pal->fPoints[i] = ((Double_t)i)/255.0;
pal->fColorAlpha[i] = 0xffff;
pal->fColorBlue[i] = pal->fColorGreen[i] = pal->fColorRed[i] = i << 8;
}
pal->fPoints[0] = 0;
pal->fPoints[255] = 1.0;
im->SetImage(*layer_pixels, UInt_t(width), pal);
delete pal;
} else {
im->SetImage(*layer_pixels, UInt_t(width), pal);
}
delete layer_pixels;
return im;
}
void TFITSHDU::Draw(Option_t *)
{
if (fType != kImageHDU) {
Warning("Draw", "cannot draw. This is not an image HDU.");
return;
}
TImage *im = ReadAsImage(0, 0);
if (im) {
Int_t width = Int_t(fSizes->GetAt(0));
Int_t height = Int_t(fSizes->GetAt(1));
TString cname, ctitle;
cname.Form("%sHDU", this->GetName());
ctitle.Form("%d x %d", width, height);
new TCanvas(cname, ctitle, width, height);
im->Draw();
}
}
TMatrixD* TFITSHDU::ReadAsMatrix(Int_t layer, Option_t *opt)
{
if (fType != kImageHDU) {
Warning("ReadAsMatrix", "this is not an image HDU.");
return 0;
}
if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
Warning("ReadAsMatrix", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
return 0;
}
if ( ((fSizes->GetSize() == 2) && (layer > 0))
|| (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
Warning("ReadAsMatrix", "layer out of bounds.");
return 0;
}
Int_t width, height;
UInt_t pixels_per_layer;
Int_t offset;
UInt_t i;
TMatrixD *mat=0;
width = Int_t(fSizes->GetAt(0));
height = Int_t(fSizes->GetAt(1));
pixels_per_layer = UInt_t(width) * UInt_t(height);
offset = layer * pixels_per_layer;
double *layer_pixels = new double[pixels_per_layer];
if ((opt[0] == 'S') || (opt[0] == 's')) {
Double_t factor, maxval=0, minval=0;
Double_t pixvalue;
for (i = 0; i < pixels_per_layer; i++) {
pixvalue = fPixels->GetAt(offset + i);
if (pixvalue > maxval) {
maxval = pixvalue;
}
if ((i == 0) || (pixvalue < minval)) {
minval = pixvalue;
}
}
if (maxval == minval) {
for (i = 0; i < pixels_per_layer; i++) {
layer_pixels[i] = 1.0;
}
} else {
factor = 1.0 / (maxval-minval);
mat = new TMatrixD(height, width);
for (i = 0; i < pixels_per_layer; i++) {
layer_pixels[i] = factor * (fPixels->GetAt(offset + i) - minval);
}
}
} else {
mat = new TMatrixD(height, width);
for (i = 0; i < pixels_per_layer; i++) {
layer_pixels[i] = fPixels->GetAt(offset + i);
}
}
if (mat) {
memcpy(mat->GetMatrixArray(), layer_pixels, pixels_per_layer*sizeof(double));
}
delete [] layer_pixels;
return mat;
}
TH1 *TFITSHDU::ReadAsHistogram()
{
if (fType != kImageHDU) {
Warning("ReadAsHistogram", "this is not an image HDU.");
return 0;
}
TH1 *result=0;
if ((fSizes->GetSize() != 1) && (fSizes->GetSize() != 2) && (fSizes->GetSize() != 3)) {
Warning("ReadAsHistogram", "could not convert image HDU to histogram because it has %d dimensions.", fSizes->GetSize());
return 0;
}
if (fSizes->GetSize() == 1) {
UInt_t Nx = UInt_t(fSizes->GetAt(0));
UInt_t x;
TH1D *h = new TH1D("", "", Int_t(Nx), 0, Nx-1);
for (x = 0; x < Nx; x++) {
Int_t nentries = Int_t(fPixels->GetAt(x));
if (nentries < 0) nentries = 0;
h->Fill(x, nentries);
}
result = h;
} else if (fSizes->GetSize() == 2) {
UInt_t Nx = UInt_t(fSizes->GetAt(0));
UInt_t Ny = UInt_t(fSizes->GetAt(1));
UInt_t x,y;
TH2D *h = new TH2D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1);
for (y = 0; y < Ny; y++) {
UInt_t offset = y * Nx;
for (x = 0; x < Nx; x++) {
Int_t nentries = Int_t(fPixels->GetAt(offset + x));
if (nentries < 0) nentries = 0;
h->Fill(x,y, nentries);
}
}
result = h;
} else if (fSizes->GetSize() == 3) {
UInt_t Nx = UInt_t(fSizes->GetAt(0));
UInt_t Ny = UInt_t(fSizes->GetAt(1));
UInt_t Nz = UInt_t(fSizes->GetAt(2));
UInt_t x,y,z;
TH3D *h = new TH3D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1, Int_t(Nz), 0, Nz-1);
for (z = 0; z < Nz; z++) {
UInt_t offset1 = z * Nx * Ny;
for (y = 0; y < Ny; y++) {
UInt_t offset2 = y * Nx;
for (x = 0; x < Nx; x++) {
Int_t nentries = Int_t(fPixels->GetAt(offset1 + offset2 + x));
if (nentries < 0) nentries = 0;
h->Fill(x, y, z, nentries);
}
}
}
result = h;
}
return result;
}
TVectorD* TFITSHDU::GetArrayRow(UInt_t row)
{
if (fType != kImageHDU) {
Warning("GetArrayRow", "this is not an image HDU.");
return 0;
}
if (fSizes->GetSize() != 2) {
Warning("GetArrayRow", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
return 0;
}
UInt_t i, offset, W,H;
W = UInt_t(fSizes->GetAt(0));
H = UInt_t(fSizes->GetAt(1));
if (row >= H) {
Warning("GetArrayRow", "index out of bounds.");
return 0;
}
offset = W * row;
double *v = new double[W];
for (i = 0; i < W; i++) {
v[i] = fPixels->GetAt(offset+i);
}
TVectorD *vec = new TVectorD(W, v);
delete [] v;
return vec;
}
TVectorD* TFITSHDU::GetArrayColumn(UInt_t col)
{
if (fType != kImageHDU) {
Warning("GetArrayColumn", "this is not an image HDU.");
return 0;
}
if (fSizes->GetSize() != 2) {
Warning("GetArrayColumn", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
return 0;
}
UInt_t i, W,H;
W = UInt_t(fSizes->GetAt(0));
H = UInt_t(fSizes->GetAt(1));
if (col >= W) {
Warning("GetArrayColumn", "index out of bounds.");
return 0;
}
double *v = new double[H];
for (i = 0; i < H; i++) {
v[i] = fPixels->GetAt(W*i+col);
}
TVectorD *vec = new TVectorD(H, v);
delete [] v;
return vec;
}
Int_t TFITSHDU::GetColumnNumber(const char *colname)
{
Int_t colnum;
for (colnum = 0; colnum < fNColumns; colnum++) {
if (fColumnsInfo[colnum].fName == colname) {
return colnum;
}
}
return -1;
}
TObjArray* TFITSHDU::GetTabStringColumn(Int_t colnum)
{
if (fType != kTableHDU) {
Warning("GetTabStringColumn", "this is not a table HDU.");
return 0;
}
if ((colnum < 0) || (colnum >= fNColumns)) {
Warning("GetTabStringColumn", "column index out of bounds.");
return 0;
}
if (fColumnsInfo[colnum].fType != kString) {
Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
return 0;
}
Int_t offset = colnum * fNRows;
TObjArray *res = new TObjArray();
for (Int_t row = 0; row < fNRows; row++) {
res->Add(new TObjString(fCells[offset + row].fString));
}
return res;
}
TObjArray* TFITSHDU::GetTabStringColumn(const char *colname)
{
if (fType != kTableHDU) {
Warning("GetTabStringColumn", "this is not a table HDU.");
return 0;
}
Int_t colnum = GetColumnNumber(colname);
if (colnum == -1) {
Warning("GetTabStringColumn", "column not found.");
return 0;
}
if (fColumnsInfo[colnum].fType != kString) {
Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
return 0;
}
Int_t offset = colnum * fNRows;
TObjArray *res = new TObjArray();
for (Int_t row = 0; row < fNRows; row++) {
res->Add(new TObjString(fCells[offset + row].fString));
}
return res;
}
TVectorD* TFITSHDU::GetTabRealVectorColumn(Int_t colnum)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorColumn", "this is not a table HDU.");
return 0;
}
if ((colnum < 0) || (colnum >= fNColumns)) {
Warning("GetTabRealVectorColumn", "column index out of bounds.");
return 0;
}
if (fColumnsInfo[colnum].fType != kRealNumber) {
Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
return 0;
} else if (fColumnsInfo[colnum].fDim > 1) {
Warning("GetTabRealVectorColumn", "attempting to read a column whose cells have embedded vectors, not real scalars. Use GetTabRealVectorCells() instead.");
return 0;
}
Int_t offset = colnum * fNRows;
Double_t *arr = new Double_t[fNRows];
for (Int_t row = 0; row < fNRows; row++) {
arr[row] = fCells[offset + row].fRealNumber;
}
TVectorD *res = new TVectorD();
res->Use(fNRows, arr);
return res;
}
TVectorD* TFITSHDU::GetTabRealVectorColumn(const char *colname)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorColumn", "this is not a table HDU.");
return 0;
}
Int_t colnum = GetColumnNumber(colname);
if (colnum == -1) {
Warning("GetTabRealVectorColumn", "column not found.");
return 0;
}
if (fColumnsInfo[colnum].fType != kRealNumber) {
Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
return 0;
} else if (fColumnsInfo[colnum].fDim > 1) {
Warning("GetTabRealVectorColumn", "attempting to read a column whose cells have embedded vectors, not real scalars. Use GetTabRealVectorCells() instead.");
return 0;
}
Int_t offset = colnum * fNRows;
Double_t *arr = new Double_t[fNRows];
for (Int_t row = 0; row < fNRows; row++) {
arr[row] = fCells[offset + row].fRealNumber;
}
TVectorD *res = new TVectorD();
res->Use(fNRows, arr);
return res;
}
Bool_t TFITSHDU::Change(const char *filter)
{
TString tmppath;
tmppath.Form("%s%s", fBaseFilePath.Data(), filter);
_release_resources();
_initialize_me();
if (kFALSE == LoadHDU(tmppath)) {
Warning("Change", "error changing HDU. Restoring the previous one...");
_release_resources();
_initialize_me();
if (kFALSE == LoadHDU(fFilePath)) {
Warning("Change", "could not restore previous HDU. This object is no longer reliable!!");
}
return kFALSE;
}
fFilePath = tmppath;
return kTRUE;
}
Bool_t TFITSHDU::Change(Int_t extension_number)
{
TString tmppath;
tmppath.Form("[%d]", extension_number);
return Change(tmppath.Data());
}
TObjArray *TFITSHDU::GetTabRealVectorCells(Int_t colnum)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorCells", "this is not a table HDU.");
return 0;
}
if ((colnum < 0) || (colnum >= fNColumns)) {
Warning("GetTabRealVectorCells", "column index out of bounds.");
return 0;
}
if (fColumnsInfo[colnum].fType != kRealNumber) {
Warning("GetTabRealVectorCells", "attempting to read a column that is not of type 'kRealNumber'.");
return 0;
}
Int_t offset = colnum * fNRows;
TObjArray *res = new TObjArray();
Int_t dim = fColumnsInfo[colnum].fDim;
for (Int_t row = 0; row < fNRows; row++) {
TVectorD *v = new TVectorD();
v->Use(dim, fCells[offset + row].fRealVector);
res->Add(v);
}
res->SetOwner(kTRUE);
return res;
}
TObjArray *TFITSHDU::GetTabRealVectorCells(const char *colname)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorCells", "this is not a table HDU.");
return 0;
}
Int_t colnum = GetColumnNumber(colname);
if (colnum == -1) {
Warning("GetTabRealVectorCells", "column not found.");
return 0;
}
return GetTabRealVectorCells(colnum);
}
TVectorD *TFITSHDU::GetTabRealVectorCell(Int_t rownum, Int_t colnum)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorCell", "this is not a table HDU.");
return 0;
}
if ((colnum < 0) || (colnum >= fNColumns)) {
Warning("GetTabRealVectorCell", "column index out of bounds.");
return 0;
}
if ((rownum < 0) || (rownum >= fNRows)) {
Warning("GetTabRealVectorCell", "row index out of bounds.");
return 0;
}
if (fColumnsInfo[colnum].fType != kRealNumber) {
Warning("GetTabRealVectorCell", "attempting to read a column that is not of type 'kRealNumber'.");
return 0;
}
TVectorD *v = new TVectorD();
v->Use(fColumnsInfo[colnum].fDim, fCells[(colnum * fNRows) + rownum].fRealVector);
return v;
}
TVectorD *TFITSHDU::GetTabRealVectorCell(Int_t rownum, const char *colname)
{
if (fType != kTableHDU) {
Warning("GetTabRealVectorCell", "this is not a table HDU.");
return 0;
}
Int_t colnum = GetColumnNumber(colname);
if (colnum == -1) {
Warning("GetTabRealVectorCell", "column not found.");
return 0;
}
return GetTabRealVectorCell(rownum, colnum);
}
const TString& TFITSHDU::GetColumnName(Int_t colnum)
{
static TString noName;
if (fType != kTableHDU) {
Error("GetColumnName", "this is not a table HDU.");
return noName;
}
if ((colnum < 0) || (colnum >= fNColumns)) {
Error("GetColumnName", "column index out of bounds.");
return noName;
}
return fColumnsInfo[colnum].fName;
}