94 void TFITSHDU::CleanFilePath(
const char *filepath_with_filter, TString &dst)
96 dst = filepath_with_filter;
98 Ssiz_t ndx = dst.Index(
"[", 1, 0, TString::kExact);
118 TFITSHDU::TFITSHDU(
const char *filepath_with_filter)
122 fFilePath = filepath_with_filter;
123 CleanFilePath(filepath_with_filter, fBaseFilePath);
125 if (kFALSE == LoadHDU(fFilePath)) {
126 _release_resources();
134 TFITSHDU::TFITSHDU(
const char *filepath, Int_t extension_number)
137 CleanFilePath(filepath, fBaseFilePath);
140 fFilePath.Form(
"%s[%d]", fBaseFilePath.Data(), extension_number);
142 if (kFALSE == LoadHDU(fFilePath)) {
143 _release_resources();
151 TFITSHDU::TFITSHDU(
const char *filepath,
const char *extension_name)
154 CleanFilePath(filepath, fBaseFilePath);
157 fFilePath.Form(
"%s[%s]", fBaseFilePath.Data(), extension_name);
160 if (kFALSE == LoadHDU(fFilePath)) {
161 _release_resources();
169 TFITSHDU::~TFITSHDU()
171 _release_resources();
177 void TFITSHDU::_release_resources()
179 if (fRecords)
delete [] fRecords;
181 if (fType == kImageHDU) {
182 if (fSizes)
delete fSizes;
183 if (fPixels)
delete fPixels;
187 for (Int_t i = 0; i < fNColumns; i++) {
188 if (fColumnsInfo[i].fType == kString) {
190 Int_t offset = i * fNRows;
191 for (Int_t row = 0; row < fNRows; row++) {
192 delete [] fCells[offset+row].fString;
194 }
else if (fColumnsInfo[i].fType == kRealVector) {
196 Int_t offset = i * fNRows;
197 for (Int_t row = 0; row < fNRows; row++) {
198 delete [] fCells[offset+row].fRealVector;
206 delete [] fColumnsInfo;
216 void TFITSHDU::_initialize_me()
222 fNColumns = fNRows = 0;
231 Bool_t TFITSHDU::LoadHDU(TString& filepath_filter)
235 char errdescr[FLEN_STATUS+1];
238 fits_open_file(&fp, filepath_filter.Data(), READONLY, &status);
239 if (status)
goto ERR;
243 fits_get_hdu_num(fp, &hdunum);
244 fNumber = Int_t(hdunum);
248 fits_get_hdu_type(fp, &hdutype, &status);
249 if (status)
goto ERR;
250 fType = (hdutype == IMAGE_HDU) ? kImageHDU : kTableHDU;
254 char keyname[FLEN_KEYWORD+1];
255 char keyvalue[FLEN_VALUE+1];
256 char comment[FLEN_COMMENT+1];
258 fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
259 if (status)
goto ERR;
261 fRecords =
new struct HDURecord[nkeys];
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;
271 fNRecords = Int_t(nkeys);
274 fExtensionName =
"PRIMARY";
275 for (
int i = 0; i < nkeys; i++) {
276 if (fRecords[i].fKeyword ==
"EXTNAME") {
277 fExtensionName = fRecords[i].fValue;
283 if (fType == kImageHDU) {
286 long *param_dimsizes;
289 fits_get_img_dim(fp, ¶m_ndims, &status);
290 if (status)
goto ERR;
291 if (param_ndims > 0) {
293 param_dimsizes =
new long[param_ndims];
294 fits_get_img_size(fp, param_ndims, param_dimsizes, &status);
296 delete [] param_dimsizes;
300 fSizes =
new TArrayI(param_ndims);
301 fSizes =
new TArrayI(param_ndims);
302 for (
int i = 0; i < param_ndims; i++) {
303 fSizes->SetAt(param_dimsizes[i], i);
306 delete [] param_dimsizes;
310 long *firstpixel =
new long[param_ndims];
314 for (
int i = 0; i < param_ndims; i++) {
315 npixels *= (long) fSizes->GetAt(i);
319 double *pixels =
new double[npixels];
321 fits_read_pix(fp, TDOUBLE, firstpixel, npixels,
322 (
void *) &nulval, (
void *) pixels, &anynul, &status);
325 delete [] firstpixel;
330 fPixels =
new TArrayD(npixels, pixels);
332 delete [] firstpixel;
337 fSizes =
new TArrayI();
338 fPixels =
new TArrayD();
347 fits_get_num_rows(fp, &table_rows, &status);
348 if (status)
goto ERR;
350 fNRows = Int_t(table_rows);
352 fits_get_num_cols(fp, &table_cols, &status);
353 if (status)
goto ERR;
355 fNColumns = Int_t(table_cols);
358 fColumnsInfo =
new struct Column[table_cols];
364 fits_get_colname(fp, CASEINSEN, (
char*)
"*", colname, &colnum, &status);
365 while (status == COL_NOT_UNIQUE)
367 fColumnsInfo[colnum-1].fName = colname;
368 fits_get_colname(fp, CASEINSEN, (
char*)
"*", colname, &colnum, &status);
370 if (status != COL_NOT_FOUND)
goto ERR;
374 fCells =
new union Cell [table_rows * table_cols];
382 for (colnum = 0, cellindex = 0; colnum < fNColumns; colnum++) {
383 fits_get_coltype(fp, colnum+1, &typecode, &repeat, &width, &status);
385 if (status)
goto ERR;
387 if ((typecode == TDOUBLE) || (typecode == TSHORT) || (typecode == TLONG)
388 || (typecode == TFLOAT) || (typecode == TLOGICAL) || (typecode == TBIT)
389 || (typecode == TBYTE) || (typecode == TSTRING)) {
391 fColumnsInfo[colnum].fType = (typecode == TSTRING) ? kString : kRealNumber;
393 if (typecode == TSTRING) {
396 fits_get_col_display_width(fp, colnum+1, &dispwidth, &status);
397 if (status)
goto ERR;
400 char *nulval = (
char*)
"";
404 if (dispwidth <= 0) {
408 array =
new char* [table_rows];
409 for (
long row = 0; row < table_rows; row++) {
410 array[row] =
new char[dispwidth+1];
414 fits_read_col(fp, TSTRING, colnum+1, 1, 1, table_rows, nulval, array, &anynul, &status);
416 for (
long row = 0; row < table_rows; row++) {
417 delete [] array[row];
425 for (
long row = 0; row < table_rows; row++) {
426 strlcpy(array[row],
"-",dispwidth+1);
432 for (
long row = 0; row < table_rows; row++) {
433 fCells[cellindex++].fString = array[row];
444 fColumnsInfo[colnum].fDim = (Int_t) repeat;
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,
460 array =
new double[table_rows * repeat];
461 fits_read_col(fp, TDOUBLE, colnum + 1, 1, 1, table_rows * repeat, &nulval, array, &anynul,
471 array =
new double[table_rows];
472 for (
long row = 0; row < table_rows; row++) {
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;
487 for (
long row = 0; row < table_rows; row++) {
488 fCells[cellindex++].fRealNumber = array[row];
492 }
else if (repeat > 1) {
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;
502 fCells[cellindex++].fRealVector = vec;
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++];
512 fCells[cellindex++].fRealVector = vec;
519 Warning(
"LoadHDU",
"error opening FITS file. Column type %d is currently not supported", typecode);
523 if (hdutype == ASCII_TBL) {
532 fits_close_file(fp, &status);
536 fits_get_errstatus(status, errdescr);
537 Warning(
"LoadHDU",
"error opening FITS file. Details: %s", errdescr);
539 if (fp) fits_close_file(fp, &status);
546 struct TFITSHDU::HDURecord* TFITSHDU::GetRecord(
const char *keyword)
548 for (
int i = 0; i < fNRecords; i++) {
549 if (fRecords[i].fKeyword == keyword) {
559 TString& TFITSHDU::GetKeywordValue(
const char *keyword)
561 HDURecord *rec = GetRecord(keyword);
565 return *(
new TString(
""));
572 void TFITSHDU::PrintHDUMetadata(
const Option_t *)
const
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());
578 printf(
"%-10s = %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data());
586 void TFITSHDU::PrintFileMetadata(
const Option_t *opt)
const
590 char errdescr[FLEN_STATUS+1];
591 int hducount, extnum;
592 int hdutype = IMAGE_HDU;
594 char extname[FLEN_CARD]=
"PRIMARY";
595 int verbose = (opt[0] ? 1 : 0);
598 fits_open_file(&fp, fBaseFilePath.Data(), READONLY, &status);
599 if (status)
goto ERR;
602 fits_get_num_hdus(fp, &hducount, &status);
603 if (status)
goto ERR;
604 printf(
"Total: %d HDUs\n", hducount);
609 fits_get_hdu_type(fp, &hdutype, &status);
610 if (status)
goto ERR;
612 if (hdutype == IMAGE_HDU) {
614 }
else if (hdutype == ASCII_TBL) {
615 exttype=
"ASCII TABLE";
617 exttype=
"BINARY TABLE";
622 char keyname[FLEN_KEYWORD+1];
623 char keyvalue[FLEN_VALUE+1];
624 char comment[FLEN_COMMENT+1];
626 fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
627 if (status)
goto ERR;
629 struct HDURecord *records =
new struct HDURecord[nkeys];
631 for (
int i = 1; i <= nkeys; i++) {
632 fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
638 records[i-1].fKeyword = keyname;
639 records[i-1].fValue = keyvalue;
640 records[i-1].fComment = comment;
642 if (strcmp(keyname,
"EXTNAME") == 0) {
644 strlcpy(extname, keyvalue,FLEN_CARD);
649 printf(
" [%d] %s (%s)\n", extnum, exttype, extname);
653 for (
int i = 0; i < nkeys; i++) {
655 printf(
" %-10s = %s / %s\n", records[i].fKeyword.Data(), records[i].fValue.Data(), records[i].fComment.Data());
657 printf(
" %-10s = %s\n", records[i].fKeyword.Data(), records[i].fValue.Data());
668 fits_movrel_hdu(fp, 1, &hdutype, &status);
669 if (status)
goto ERR;
674 fits_close_file(fp, &status);
678 fits_get_errstatus(status, errdescr);
679 Warning(
"PrintFileMetadata",
"error opening FITS file. Details: %s", errdescr);
681 if (fp) fits_close_file(fp, &status);
687 void TFITSHDU::PrintColumnInfo(
const Option_t *)
const
689 if (fType != kTableHDU) {
690 Warning(
"PrintColumnInfo",
"this is not a table HDU.");
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");
702 void TFITSHDU::PrintFullTable(
const Option_t *)
const
706 if (fType != kTableHDU) {
707 Warning(
"PrintColumnInfo",
"this is not a table HDU.");
714 for (Int_t col = 0; col < fNColumns; col++) {
715 printed_chars += printf(
"%-10s| ", fColumnsInfo[col].fName.Data());
718 while(printed_chars--) {
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);
729 printed_chars = printf(
"%.2lg", fCells[col * fNRows + row].fRealNumber);
731 while (printed_chars < 0) {
737 if (col <= fNColumns - 1) printf(
"| ");
753 void TFITSHDU::Print(
const Option_t *opt)
const
755 if ((opt[0] ==
'F') || (opt[0] ==
'f')) {
756 PrintFileMetadata((opt[1] ==
'+') ?
"+" :
"");
757 }
else if ((opt[0] ==
'T') || (opt[0] ==
't')) {
765 PrintHDUMetadata(
"");
775 TImage *TFITSHDU::ReadAsImage(Int_t layer, TImagePalette *pal)
777 if (fType != kImageHDU) {
778 Warning(
"ReadAsImage",
"this is not an image HDU.");
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());
789 UInt_t pixels_per_layer;
791 width = Int_t(fSizes->GetAt(0));
792 height = Int_t(fSizes->GetAt(1));
794 pixels_per_layer = UInt_t(width) * UInt_t(height);
796 if ( ((fSizes->GetSize() == 2) && (layer > 0))
797 || (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
799 Warning(
"ReadAsImage",
"layer out of bounds.");
804 Double_t maxval = 0, minval = 0;
807 Int_t offset = layer * pixels_per_layer;
809 for (i = 0; i < pixels_per_layer; i++) {
810 pixvalue = fPixels->GetAt(offset + i);
812 if (pixvalue > maxval) {
816 if ((i == 0) || (pixvalue < minval)) {
823 TImage *im = TImage::Create();
825 TArrayD *layer_pixels =
new TArrayD(pixels_per_layer);
828 if (maxval == minval) {
830 for (i = 0; i < pixels_per_layer; i++) {
831 layer_pixels->SetAt(255.0, i);
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) ;
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;
850 pal->fPoints[255] = 1.0;
852 im->SetImage(*layer_pixels, UInt_t(width), pal);
857 im->SetImage(*layer_pixels, UInt_t(width), pal);
869 void TFITSHDU::Draw(Option_t *)
871 if (fType != kImageHDU) {
872 Warning(
"Draw",
"cannot draw. This is not an image HDU.");
876 TImage *im = ReadAsImage(0, 0);
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);
896 TMatrixD* TFITSHDU::ReadAsMatrix(Int_t layer, Option_t *opt)
898 if (fType != kImageHDU) {
899 Warning(
"ReadAsMatrix",
"this is not an image HDU.");
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());
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.");
916 UInt_t pixels_per_layer;
921 width = Int_t(fSizes->GetAt(0));
922 height = Int_t(fSizes->GetAt(1));
924 pixels_per_layer = UInt_t(width) * UInt_t(height);
925 offset = layer * pixels_per_layer;
927 double *layer_pixels =
new double[pixels_per_layer];
929 if ((opt[0] ==
'S') || (opt[0] ==
's')) {
932 Double_t factor, maxval=0, minval=0;
934 for (i = 0; i < pixels_per_layer; i++) {
935 pixvalue = fPixels->GetAt(offset + i);
937 if (pixvalue > maxval) {
941 if ((i == 0) || (pixvalue < minval)) {
946 if (maxval == minval) {
948 for (i = 0; i < pixels_per_layer; i++) {
949 layer_pixels[i] = 1.0;
952 factor = 1.0 / (maxval-minval);
953 mat =
new TMatrixD(height, width);
955 for (i = 0; i < pixels_per_layer; i++) {
956 layer_pixels[i] = factor * (fPixels->GetAt(offset + i) - minval);
962 mat =
new TMatrixD(height, width);
964 for (i = 0; i < pixels_per_layer; i++) {
965 layer_pixels[i] = fPixels->GetAt(offset + i);
971 memcpy(mat->GetMatrixArray(), layer_pixels, pixels_per_layer*
sizeof(double));
974 delete [] layer_pixels;
989 TH1 *TFITSHDU::ReadAsHistogram()
991 if (fType != kImageHDU) {
992 Warning(
"ReadAsHistogram",
"this is not an image HDU.");
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());
1003 if (fSizes->GetSize() == 1) {
1005 UInt_t Nx = UInt_t(fSizes->GetAt(0));
1008 TH1D *h =
new TH1D(
"",
"", Int_t(Nx), 0, Nx-1);
1010 for (x = 0; x < Nx; x++) {
1011 Int_t nentries = Int_t(fPixels->GetAt(x));
1012 if (nentries < 0) nentries = 0;
1013 h->Fill(x, nentries);
1018 }
else if (fSizes->GetSize() == 2) {
1020 UInt_t Nx = UInt_t(fSizes->GetAt(0));
1021 UInt_t Ny = UInt_t(fSizes->GetAt(1));
1024 TH2D *h =
new TH2D(
"",
"", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1);
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;
1031 h->Fill(x,y, nentries);
1037 }
else if (fSizes->GetSize() == 3) {
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));
1044 TH3D *h =
new TH3D(
"",
"", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1, Int_t(Nz), 0, Nz-1);
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;
1054 h->Fill(x, y, z, nentries);
1068 TVectorD* TFITSHDU::GetArrayRow(UInt_t row)
1070 if (fType != kImageHDU) {
1071 Warning(
"GetArrayRow",
"this is not an image HDU.");
1075 if (fSizes->GetSize() != 2) {
1076 Warning(
"GetArrayRow",
"could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
1080 UInt_t i, offset, W,H;
1082 W = UInt_t(fSizes->GetAt(0));
1083 H = UInt_t(fSizes->GetAt(1));
1087 Warning(
"GetArrayRow",
"index out of bounds.");
1092 double *v =
new double[W];
1094 for (i = 0; i < W; i++) {
1095 v[i] = fPixels->GetAt(offset+i);
1098 TVectorD *vec =
new TVectorD(W, v);
1108 TVectorD* TFITSHDU::GetArrayColumn(UInt_t col)
1110 if (fType != kImageHDU) {
1111 Warning(
"GetArrayColumn",
"this is not an image HDU.");
1115 if (fSizes->GetSize() != 2) {
1116 Warning(
"GetArrayColumn",
"could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
1122 W = UInt_t(fSizes->GetAt(0));
1123 H = UInt_t(fSizes->GetAt(1));
1127 Warning(
"GetArrayColumn",
"index out of bounds.");
1131 double *v =
new double[H];
1133 for (i = 0; i < H; i++) {
1134 v[i] = fPixels->GetAt(W*i+col);
1137 TVectorD *vec =
new TVectorD(H, v);
1148 Int_t TFITSHDU::GetColumnNumber(
const char *colname)
1151 for (colnum = 0; colnum < fNColumns; colnum++) {
1152 if (fColumnsInfo[colnum].fName == colname) {
1162 TObjArray* TFITSHDU::GetTabStringColumn(Int_t colnum)
1164 if (fType != kTableHDU) {
1165 Warning(
"GetTabStringColumn",
"this is not a table HDU.");
1169 if ((colnum < 0) || (colnum >= fNColumns)) {
1170 Warning(
"GetTabStringColumn",
"column index out of bounds.");
1174 if (fColumnsInfo[colnum].fType != kString) {
1175 Warning(
"GetTabStringColumn",
"attempting to read a column that is not of type 'kString'.");
1179 Int_t offset = colnum * fNRows;
1181 TObjArray *res =
new TObjArray();
1182 for (Int_t row = 0; row < fNRows; row++) {
1183 res->Add(
new TObjString(fCells[offset + row].fString));
1192 TObjArray* TFITSHDU::GetTabStringColumn(
const char *colname)
1194 if (fType != kTableHDU) {
1195 Warning(
"GetTabStringColumn",
"this is not a table HDU.");
1200 Int_t colnum = GetColumnNumber(colname);
1203 Warning(
"GetTabStringColumn",
"column not found.");
1207 if (fColumnsInfo[colnum].fType != kString) {
1208 Warning(
"GetTabStringColumn",
"attempting to read a column that is not of type 'kString'.");
1212 Int_t offset = colnum * fNRows;
1214 TObjArray *res =
new TObjArray();
1215 for (Int_t row = 0; row < fNRows; row++) {
1216 res->Add(
new TObjString(fCells[offset + row].fString));
1225 TVectorD* TFITSHDU::GetTabRealVectorColumn(Int_t colnum)
1227 if (fType != kTableHDU) {
1228 Warning(
"GetTabRealVectorColumn",
"this is not a table HDU.");
1232 if ((colnum < 0) || (colnum >= fNColumns)) {
1233 Warning(
"GetTabRealVectorColumn",
"column index out of bounds.");
1237 if (fColumnsInfo[colnum].fType != kRealNumber) {
1238 Warning(
"GetTabRealVectorColumn",
"attempting to read a column that is not of type 'kRealNumber'.");
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.");
1245 Int_t offset = colnum * fNRows;
1247 Double_t *arr =
new Double_t[fNRows];
1249 for (Int_t row = 0; row < fNRows; row++) {
1250 arr[row] = fCells[offset + row].fRealNumber;
1253 TVectorD *res =
new TVectorD();
1254 res->Use(fNRows, arr);
1262 TVectorD* TFITSHDU::GetTabRealVectorColumn(
const char *colname)
1264 if (fType != kTableHDU) {
1265 Warning(
"GetTabRealVectorColumn",
"this is not a table HDU.");
1269 Int_t colnum = GetColumnNumber(colname);
1272 Warning(
"GetTabRealVectorColumn",
"column not found.");
1276 if (fColumnsInfo[colnum].fType != kRealNumber) {
1277 Warning(
"GetTabRealVectorColumn",
"attempting to read a column that is not of type 'kRealNumber'.");
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.");
1284 Int_t offset = colnum * fNRows;
1286 Double_t *arr =
new Double_t[fNRows];
1288 for (Int_t row = 0; row < fNRows; row++) {
1289 arr[row] = fCells[offset + row].fRealNumber;
1292 TVectorD *res =
new TVectorD();
1293 res->Use(fNRows, arr);
1306 Bool_t TFITSHDU::Change(
const char *filter)
1309 tmppath.Form(
"%s%s", fBaseFilePath.Data(), filter);
1311 _release_resources();
1314 if (kFALSE == LoadHDU(tmppath)) {
1316 Warning(
"Change",
"error changing HDU. Restoring the previous one...");
1318 _release_resources();
1321 if (kFALSE == LoadHDU(fFilePath)) {
1322 Warning(
"Change",
"could not restore previous HDU. This object is no longer reliable!!");
1328 fFilePath = tmppath;
1335 Bool_t TFITSHDU::Change(Int_t extension_number)
1338 tmppath.Form(
"[%d]", extension_number);
1340 return Change(tmppath.Data());
1346 TObjArray *TFITSHDU::GetTabRealVectorCells(Int_t colnum)
1348 if (fType != kTableHDU) {
1349 Warning(
"GetTabRealVectorCells",
"this is not a table HDU.");
1353 if ((colnum < 0) || (colnum >= fNColumns)) {
1354 Warning(
"GetTabRealVectorCells",
"column index out of bounds.");
1358 if (fColumnsInfo[colnum].fType != kRealNumber) {
1359 Warning(
"GetTabRealVectorCells",
"attempting to read a column that is not of type 'kRealNumber'.");
1363 Int_t offset = colnum * fNRows;
1365 TObjArray *res =
new TObjArray();
1366 Int_t dim = fColumnsInfo[colnum].fDim;
1368 for (Int_t row = 0; row < fNRows; row++) {
1369 TVectorD *v =
new TVectorD();
1370 v->Use(dim, fCells[offset + row].fRealVector);
1376 res->SetOwner(kTRUE);
1384 TObjArray *TFITSHDU::GetTabRealVectorCells(
const char *colname)
1386 if (fType != kTableHDU) {
1387 Warning(
"GetTabRealVectorCells",
"this is not a table HDU.");
1391 Int_t colnum = GetColumnNumber(colname);
1394 Warning(
"GetTabRealVectorCells",
"column not found.");
1398 return GetTabRealVectorCells(colnum);
1404 TVectorD *TFITSHDU::GetTabRealVectorCell(Int_t rownum, Int_t colnum)
1406 if (fType != kTableHDU) {
1407 Warning(
"GetTabRealVectorCell",
"this is not a table HDU.");
1411 if ((colnum < 0) || (colnum >= fNColumns)) {
1412 Warning(
"GetTabRealVectorCell",
"column index out of bounds.");
1416 if ((rownum < 0) || (rownum >= fNRows)) {
1417 Warning(
"GetTabRealVectorCell",
"row index out of bounds.");
1421 if (fColumnsInfo[colnum].fType != kRealNumber) {
1422 Warning(
"GetTabRealVectorCell",
"attempting to read a column that is not of type 'kRealNumber'.");
1426 TVectorD *v =
new TVectorD();
1427 v->Use(fColumnsInfo[colnum].fDim, fCells[(colnum * fNRows) + rownum].fRealVector);
1435 TVectorD *TFITSHDU::GetTabRealVectorCell(Int_t rownum,
const char *colname)
1437 if (fType != kTableHDU) {
1438 Warning(
"GetTabRealVectorCell",
"this is not a table HDU.");
1442 Int_t colnum = GetColumnNumber(colname);
1445 Warning(
"GetTabRealVectorCell",
"column not found.");
1449 return GetTabRealVectorCell(rownum, colnum);
1457 const TString& TFITSHDU::GetColumnName(Int_t colnum)
1459 static TString noName;
1460 if (fType != kTableHDU) {
1461 Error(
"GetColumnName",
"this is not a table HDU.");
1465 if ((colnum < 0) || (colnum >= fNColumns)) {
1466 Error(
"GetColumnName",
"column index out of bounds.");
1469 return fColumnsInfo[colnum].fName;