00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065 #include "TFITS.h"
00066 #include "TROOT.h"
00067 #include "TImage.h"
00068 #include "TArrayI.h"
00069 #include "TArrayD.h"
00070 #include "TH1D.h"
00071 #include "TH2D.h"
00072 #include "TH3D.h"
00073 #include "TVectorD.h"
00074 #include "TMatrixD.h"
00075 #include "TObjArray.h"
00076 #include "TObjString.h"
00077 #include "TCanvas.h"
00078
00079 #include "fitsio.h"
00080 #include <stdlib.h>
00081
00082
00083 ClassImp(TFITSHDU)
00084
00085
00086
00087
00088
00089
00090 void TFITSHDU::CleanFilePath(const char *filepath_with_filter, TString &dst)
00091 {
00092
00093
00094 dst = filepath_with_filter;
00095
00096 Ssiz_t ndx = dst.Index("[", 1, 0, TString::kExact);
00097 if (ndx != kNPOS) {
00098 dst.Resize(ndx);
00099 }
00100 }
00101
00102
00103
00104 TFITSHDU::TFITSHDU(const char *filepath_with_filter)
00105 {
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117 _initialize_me();
00118
00119 fFilePath = filepath_with_filter;
00120 CleanFilePath(filepath_with_filter, fBaseFilePath);
00121
00122 if (kFALSE == LoadHDU(fFilePath)) {
00123 _release_resources();
00124 throw -1;
00125 }
00126 }
00127
00128
00129 TFITSHDU::TFITSHDU(const char *filepath, Int_t extension_number)
00130 {
00131
00132
00133 _initialize_me();
00134 CleanFilePath(filepath, fBaseFilePath);
00135
00136
00137 fFilePath.Form("%s[%d]", fBaseFilePath.Data(), extension_number);
00138
00139 if (kFALSE == LoadHDU(fFilePath)) {
00140 _release_resources();
00141 throw -1;
00142 }
00143 }
00144
00145
00146 TFITSHDU::TFITSHDU(const char *filepath, const char *extension_name)
00147 {
00148
00149
00150 _initialize_me();
00151 CleanFilePath(filepath, fBaseFilePath);
00152
00153
00154 fFilePath.Form("%s[%s]", fBaseFilePath.Data(), extension_name);
00155
00156
00157 if (kFALSE == LoadHDU(fFilePath)) {
00158 _release_resources();
00159 throw -1;
00160 }
00161 }
00162
00163
00164 TFITSHDU::~TFITSHDU()
00165 {
00166
00167
00168 _release_resources();
00169 }
00170
00171
00172 void TFITSHDU::_release_resources()
00173 {
00174
00175
00176 if (fRecords) delete [] fRecords;
00177
00178 if (fType == kImageHDU) {
00179 if (fSizes) delete fSizes;
00180 if (fPixels) delete fPixels;
00181 } else {
00182 if (fColumnsInfo) {
00183 if (fCells) {
00184 for (Int_t i = 0; i < fNColumns; i++) {
00185 if (fColumnsInfo[i].fType == kString) {
00186
00187 Int_t offset = i * fNRows;
00188 for (Int_t row = 0; row < fNRows; row++) {
00189 delete [] fCells[offset+row].fString;
00190 }
00191 }
00192 }
00193
00194 delete [] fCells;
00195 }
00196
00197 delete [] fColumnsInfo;
00198 }
00199
00200
00201 }
00202 }
00203
00204
00205 void TFITSHDU::_initialize_me()
00206 {
00207
00208
00209 fRecords = 0;
00210 fPixels = 0;
00211 fSizes = 0;
00212 fColumnsInfo = 0;
00213 fNColumns = fNRows = 0;
00214 fCells = 0;
00215 }
00216
00217
00218 Bool_t TFITSHDU::LoadHDU(TString& filepath_filter)
00219 {
00220
00221
00222
00223
00224 fitsfile *fp=0;
00225 int status = 0;
00226 char errdescr[FLEN_STATUS+1];
00227
00228
00229 fits_open_file(&fp, filepath_filter.Data(), READONLY, &status);
00230 if (status) goto ERR;
00231
00232
00233 int hdunum;
00234 fits_get_hdu_num(fp, &hdunum);
00235 fNumber = Int_t(hdunum);
00236
00237
00238 int hdutype;
00239 fits_get_hdu_type(fp, &hdutype, &status);
00240 if (status) goto ERR;
00241 fType = (hdutype == IMAGE_HDU) ? kImageHDU : kTableHDU;
00242
00243
00244 int nkeys, morekeys;
00245 char keyname[FLEN_KEYWORD+1];
00246 char keyvalue[FLEN_VALUE+1];
00247 char comment[FLEN_COMMENT+1];
00248
00249 fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
00250 if (status) goto ERR;
00251
00252 fRecords = new struct HDURecord[nkeys];
00253
00254 for (int i = 1; i <= nkeys; i++) {
00255 fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
00256 if (status) goto ERR;
00257 fRecords[i-1].fKeyword = keyname;
00258 fRecords[i-1].fValue = keyvalue;
00259 fRecords[i-1].fComment = comment;
00260 }
00261
00262 fNRecords = Int_t(nkeys);
00263
00264
00265 fExtensionName = "PRIMARY";
00266 for (int i = 0; i < nkeys; i++) {
00267 if (fRecords[i].fKeyword == "EXTNAME") {
00268 fExtensionName = fRecords[i].fValue;
00269 break;
00270 }
00271 }
00272
00273
00274 if (fType == kImageHDU) {
00275
00276 int param_ndims=0;
00277 long *param_dimsizes;
00278
00279
00280 fits_get_img_dim(fp, ¶m_ndims, &status);
00281 if (status) goto ERR;
00282 if (param_ndims > 0) {
00283
00284 param_dimsizes = new long[param_ndims];
00285 fits_get_img_size(fp, param_ndims, param_dimsizes, &status);
00286 if (status) goto ERR;
00287
00288 fSizes = new TArrayI(param_ndims);
00289 fSizes = new TArrayI(param_ndims);
00290 for (int i = 0; i < param_ndims; i++) {
00291 fSizes->SetAt(param_dimsizes[i], i);
00292 }
00293
00294 delete [] param_dimsizes;
00295
00296
00297 int anynul;
00298 long *firstpixel = new long[param_ndims];
00299 double nulval = 0;
00300 long npixels = 1;
00301
00302 for (int i = 0; i < param_ndims; i++) {
00303 npixels *= (long) fSizes->GetAt(i);
00304 firstpixel[i] = 1;
00305 }
00306
00307 double *pixels = new double[npixels];
00308
00309 fits_read_pix(fp, TDOUBLE, firstpixel, npixels,
00310 (void *) &nulval, (void *) pixels, &anynul, &status);
00311
00312 if (status) {
00313 delete [] firstpixel;
00314 delete [] pixels;
00315 goto ERR;
00316 }
00317
00318 fPixels = new TArrayD(npixels, pixels);
00319
00320 delete [] firstpixel;
00321 delete [] pixels;
00322
00323 } else {
00324
00325 fSizes = new TArrayI();
00326 fPixels = new TArrayD();
00327 }
00328 } else {
00329
00330
00331
00332 long table_rows;
00333 int table_cols;
00334
00335 fits_get_num_rows(fp, &table_rows, &status);
00336 if (status) goto ERR;
00337
00338 fNRows = Int_t(table_rows);
00339
00340 fits_get_num_cols(fp, &table_cols, &status);
00341 if (status) goto ERR;
00342
00343 fNColumns = Int_t(table_cols);
00344
00345
00346 fColumnsInfo = new struct Column[table_cols];
00347
00348
00349 char colname[80];
00350 int colnum;
00351
00352 fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
00353 while (status == COL_NOT_UNIQUE)
00354 {
00355 fColumnsInfo[colnum-1].fName = colname;
00356 fits_get_colname(fp, CASEINSEN, (char*) "*", colname, &colnum, &status);
00357 }
00358 if (status != COL_NOT_FOUND) goto ERR;
00359 status = 0;
00360
00361
00362 fCells = new union Cell [table_rows * table_cols];
00363
00364
00365 int typecode;
00366 long repeat, width;
00367 Int_t cellindex;
00368
00369 fColumnTypes = new enum EColumnTypes [table_cols];
00370
00371 for (colnum = 0, cellindex = 0; colnum < fNColumns; colnum++) {
00372 fits_get_coltype(fp, colnum+1, &typecode, &repeat, &width, &status);
00373 if (status) goto ERR;
00374
00375 if ((typecode == TDOUBLE) || (typecode == TSHORT) || (typecode == TLONG)
00376 || (typecode == TFLOAT) || (typecode == TLOGICAL) || (typecode == TBIT)
00377 || (typecode == TBYTE) || (typecode == TSTRING)) {
00378
00379 fColumnsInfo[colnum].fType = (typecode == TSTRING) ? kString : kRealNumber;
00380
00381 if (typecode == TSTRING) {
00382
00383 int dispwidth=0;
00384 fits_get_col_display_width(fp, colnum+1, &dispwidth, &status);
00385 if (status) goto ERR;
00386
00387
00388 char *nulval = (char*) "";
00389 int anynul=0;
00390 char **array;
00391
00392 if (dispwidth <= 0) {
00393 dispwidth = 1;
00394 }
00395
00396 array = new char* [table_rows];
00397 for (long row = 0; row < table_rows; row++) {
00398 array[row] = new char[dispwidth+1];
00399 }
00400
00401 if (repeat > 0) {
00402 fits_read_col(fp, TSTRING, colnum+1, 1, 1, table_rows, nulval, array, &anynul, &status);
00403 if (status) {
00404 for (long row = 0; row < table_rows; row++) {
00405 delete [] array[row];
00406 }
00407 delete [] array;
00408 goto ERR;
00409 }
00410
00411 } else {
00412
00413 for (long row = 0; row < table_rows; row++) {
00414 strlcpy(array[row], "-",dispwidth+1);
00415 }
00416 }
00417
00418
00419
00420 for (long row = 0; row < table_rows; row++) {
00421 fCells[cellindex++].fString = array[row];
00422 }
00423
00424 delete [] array;
00425
00426
00427 } else {
00428
00429 double nulval = 0;
00430 int anynul=0;
00431 double *array = new double [table_rows];
00432
00433 if (repeat > 0) {
00434 fits_read_col(fp, TDOUBLE, colnum+1, 1, 1, table_rows, &nulval, array, &anynul, &status);
00435
00436 if (status) {
00437 delete [] array;
00438 goto ERR;
00439 }
00440 } else {
00441
00442 for (long row = 0; row < table_rows; row++) {
00443 array[row] = 0.0;
00444 }
00445 }
00446
00447
00448 for (long row = 0; row < table_rows; row++) {
00449 fCells[cellindex++].fRealNumber = array[row];
00450 }
00451
00452 delete [] array;
00453
00454 }
00455
00456 } else {
00457 Warning("LoadHDU", "error opening FITS file. Column type %d is currently not supported", typecode);
00458 }
00459 }
00460
00461
00462
00463 if (hdutype == ASCII_TBL) {
00464
00465
00466 } else {
00467
00468 }
00469
00470 }
00471
00472
00473 fits_close_file(fp, &status);
00474 return kTRUE;
00475
00476 ERR:
00477 fits_get_errstatus(status, errdescr);
00478 Warning("LoadHDU", "error opening FITS file. Details: %s", errdescr);
00479 status = 0;
00480 if (fp) fits_close_file(fp, &status);
00481 return kFALSE;
00482 }
00483
00484
00485 struct TFITSHDU::HDURecord* TFITSHDU::GetRecord(const char *keyword)
00486 {
00487
00488
00489 for (int i = 0; i < fNRecords; i++) {
00490 if (fRecords[i].fKeyword == keyword) {
00491 return &fRecords[i];
00492 }
00493 }
00494 return 0;
00495 }
00496
00497
00498 TString& TFITSHDU::GetKeywordValue(const char *keyword)
00499 {
00500
00501
00502 HDURecord *rec = GetRecord(keyword);
00503 if (rec) {
00504 return rec->fValue;
00505 } else {
00506 return *(new TString(""));
00507 }
00508 }
00509
00510
00511 void TFITSHDU::PrintHDUMetadata(const Option_t *) const
00512 {
00513
00514
00515 for (int i = 0; i < fNRecords; i++) {
00516 if (fRecords[i].fComment.Length() > 0) {
00517 printf("%-10s = %s / %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data(), fRecords[i].fComment.Data());
00518 } else {
00519 printf("%-10s = %s\n", fRecords[i].fKeyword.Data(), fRecords[i].fValue.Data());
00520 }
00521 }
00522 }
00523
00524
00525 void TFITSHDU::PrintFileMetadata(const Option_t *opt) const
00526 {
00527
00528
00529 fitsfile *fp=0;
00530 int status = 0;
00531 char errdescr[FLEN_STATUS+1];
00532 int hducount, extnum;
00533 int hdutype = IMAGE_HDU;
00534 const char *exttype;
00535 char extname[FLEN_CARD]="PRIMARY";
00536 int verbose = (opt[0] ? 1 : 0);
00537
00538
00539 fits_open_file(&fp, fBaseFilePath.Data(), READONLY, &status);
00540 if (status) goto ERR;
00541
00542
00543 fits_get_num_hdus(fp, &hducount, &status);
00544 if (status) goto ERR;
00545 printf("Total: %d HDUs\n", hducount);
00546
00547 extnum = 0;
00548 while(hducount) {
00549
00550 fits_get_hdu_type(fp, &hdutype, &status);
00551 if (status) goto ERR;
00552
00553 if (hdutype == IMAGE_HDU) {
00554 exttype="IMAGE";
00555 } else if (hdutype == ASCII_TBL) {
00556 exttype="ASCII TABLE";
00557 } else {
00558 exttype="BINARY TABLE";
00559 }
00560
00561
00562 int nkeys, morekeys;
00563 char keyname[FLEN_KEYWORD+1];
00564 char keyvalue[FLEN_VALUE+1];
00565 char comment[FLEN_COMMENT+1];
00566
00567 fits_get_hdrspace(fp, &nkeys, &morekeys, &status);
00568 if (status) goto ERR;
00569
00570 struct HDURecord *records = new struct HDURecord[nkeys];
00571
00572 for (int i = 1; i <= nkeys; i++) {
00573 fits_read_keyn(fp, i, keyname, keyvalue, comment, &status);
00574 if (status) {
00575 delete [] records;
00576 goto ERR;
00577 }
00578
00579 records[i-1].fKeyword = keyname;
00580 records[i-1].fValue = keyvalue;
00581 records[i-1].fComment = comment;
00582
00583 if (strcmp(keyname, "EXTNAME") == 0) {
00584
00585 strlcpy(extname, keyvalue,FLEN_CARD);
00586 }
00587 }
00588
00589
00590 printf(" [%d] %s (%s)\n", extnum, exttype, extname);
00591
00592
00593 if (verbose) {
00594 for (int i = 0; i < nkeys; i++) {
00595 if (comment[0]) {
00596 printf(" %-10s = %s / %s\n", records[i].fKeyword.Data(), records[i].fValue.Data(), records[i].fComment.Data());
00597 } else {
00598 printf(" %-10s = %s\n", records[i].fKeyword.Data(), records[i].fValue.Data());
00599 }
00600 }
00601 }
00602 printf("\n");
00603
00604 delete [] records;
00605
00606 hducount--;
00607 extnum++;
00608 if (hducount){
00609 fits_movrel_hdu(fp, 1, &hdutype, &status);
00610 if (status) goto ERR;
00611 }
00612 }
00613
00614
00615 fits_close_file(fp, &status);
00616 return;
00617
00618 ERR:
00619 fits_get_errstatus(status, errdescr);
00620 Warning("PrintFileMetadata", "error opening FITS file. Details: %s", errdescr);
00621 status = 0;
00622 if (fp) fits_close_file(fp, &status);
00623 }
00624
00625
00626 void TFITSHDU::PrintColumnInfo(const Option_t *) const
00627 {
00628
00629
00630 if (fType != kTableHDU) {
00631 Warning("PrintColumnInfo", "this is not a table HDU.");
00632 return;
00633 }
00634
00635 for (Int_t i = 0; i < fNColumns; i++) {
00636 printf("%-20s : %s\n", fColumnsInfo[i].fName.Data(), (fColumnsInfo[i].fType == kRealNumber) ? "REAL NUMBER" : "STRING");
00637 }
00638 }
00639
00640
00641 void TFITSHDU::PrintFullTable(const Option_t *) const
00642 {
00643
00644
00645 int printed_chars;
00646
00647 if (fType != kTableHDU) {
00648 Warning("PrintColumnInfo", "this is not a table HDU.");
00649 return;
00650 }
00651
00652
00653 putchar('\n');
00654 printed_chars = 0;
00655 for (Int_t col = 0; col < fNColumns; col++) {
00656 printed_chars += printf("%-10s| ", fColumnsInfo[col].fName.Data());
00657 }
00658 putchar('\n');
00659 while(printed_chars--) {
00660 putchar('-');
00661 }
00662 putchar('\n');
00663
00664
00665 for (Int_t row = 0; row < fNRows; row++) {
00666 for (Int_t col = 0; col < fNColumns; col++) {
00667 if (fColumnsInfo[col].fType == kString) {
00668 printf("%-10s", fCells[col * fNRows + row].fString);
00669 } else {
00670 printed_chars = printf("%.2lg", fCells[col * fNRows + row].fRealNumber);
00671 printed_chars -= 10;
00672 while (printed_chars < 0) {
00673 putchar(' ');
00674 printed_chars++;
00675 }
00676 }
00677
00678 if (col <= fNColumns - 1) printf("| ");
00679 }
00680 printf("\n");
00681 }
00682 }
00683
00684
00685 void TFITSHDU::Print(const Option_t *opt) const
00686 {
00687
00688
00689
00690
00691
00692
00693
00694
00695 if ((opt[0] == 'F') || (opt[0] == 'f')) {
00696 PrintFileMetadata((opt[1] == '+') ? "+" : "");
00697 } else if ((opt[0] == 'T') || (opt[0] == 't')) {
00698 if (opt[1] == '+') {
00699 PrintFullTable("");
00700 } else {
00701 PrintColumnInfo("");
00702 }
00703
00704 } else {
00705 PrintHDUMetadata("");
00706 }
00707 }
00708
00709
00710
00711 TImage *TFITSHDU::ReadAsImage(Int_t layer, TImagePalette *pal)
00712 {
00713
00714
00715
00716
00717 if (fType != kImageHDU) {
00718 Warning("ReadAsImage", "this is not an image HDU.");
00719 return 0;
00720 }
00721
00722
00723 if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
00724 Warning("ReadAsImage", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
00725 return 0;
00726 }
00727
00728 Int_t width, height;
00729 UInt_t pixels_per_layer;
00730
00731 width = Int_t(fSizes->GetAt(0));
00732 height = Int_t(fSizes->GetAt(1));
00733
00734 pixels_per_layer = UInt_t(width) * UInt_t(height);
00735
00736 if ( ((fSizes->GetSize() == 2) && (layer > 0))
00737 || (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
00738
00739 Warning("ReadAsImage", "layer out of bounds.");
00740 return 0;
00741 }
00742
00743
00744 Double_t maxval = 0, minval = 0;
00745 register UInt_t i;
00746 Double_t pixvalue;
00747 Int_t offset = layer * pixels_per_layer;
00748
00749 for (i = 0; i < pixels_per_layer; i++) {
00750 pixvalue = fPixels->GetAt(offset + i);
00751
00752 if (pixvalue > maxval) {
00753 maxval = pixvalue;
00754 }
00755
00756 if ((i == 0) || (pixvalue < minval)) {
00757 minval = pixvalue;
00758 }
00759 }
00760
00761
00762
00763 TImage *im = TImage::Create();
00764 TArrayD *layer_pixels = new TArrayD(pixels_per_layer);
00765
00766
00767 if (maxval == minval) {
00768
00769 for (i = 0; i < pixels_per_layer; i++) {
00770 layer_pixels->SetAt(255.0, i);
00771 }
00772 } else {
00773 Double_t factor = 255.0 / (maxval-minval);
00774 for (i = 0; i < pixels_per_layer; i++) {
00775 pixvalue = fPixels->GetAt(offset + i);
00776 layer_pixels->SetAt(factor * (pixvalue-minval), i) ;
00777 }
00778 }
00779
00780 if (pal == 0) {
00781
00782 pal = new TImagePalette(256);
00783 for (i = 0; i < 256; i++) {
00784 pal->fPoints[i] = ((Double_t)i)/255.0;
00785 pal->fColorAlpha[i] = 0xffff;
00786 pal->fColorBlue[i] = pal->fColorGreen[i] = pal->fColorRed[i] = i << 8;
00787 }
00788 pal->fPoints[0] = 0;
00789 pal->fPoints[255] = 1.0;
00790
00791 im->SetImage(*layer_pixels, UInt_t(width), pal);
00792
00793 delete pal;
00794
00795 } else {
00796 im->SetImage(*layer_pixels, UInt_t(width), pal);
00797 }
00798
00799 delete layer_pixels;
00800
00801 return im;
00802 }
00803
00804
00805 void TFITSHDU::Draw(Option_t *)
00806 {
00807
00808
00809
00810 if (fType != kImageHDU) {
00811 Warning("Draw", "cannot draw. This is not an image HDU.");
00812 return;
00813 }
00814
00815 TImage *im = ReadAsImage(0, 0);
00816 if (im) {
00817 Int_t width = Int_t(fSizes->GetAt(0));
00818 Int_t height = Int_t(fSizes->GetAt(1));
00819 TString cname, ctitle;
00820 cname.Form("%sHDU", this->GetName());
00821 ctitle.Form("%d x %d", width, height);
00822 new TCanvas(cname, ctitle, width, height);
00823 im->Draw();
00824 }
00825 }
00826
00827
00828
00829 TMatrixD* TFITSHDU::ReadAsMatrix(Int_t layer, Option_t *opt)
00830 {
00831
00832
00833
00834
00835
00836
00837 if (fType != kImageHDU) {
00838 Warning("ReadAsMatrix", "this is not an image HDU.");
00839 return 0;
00840 }
00841
00842 if (((fSizes->GetSize() != 2) && (fSizes->GetSize() != 3) && (fSizes->GetSize() != 4)) || ((fSizes->GetSize() == 4) && (fSizes->GetAt(3) > 1))) {
00843 Warning("ReadAsMatrix", "could not convert image HDU to image because it has %d dimensions.", fSizes->GetSize());
00844 return 0;
00845 }
00846
00847
00848 if ( ((fSizes->GetSize() == 2) && (layer > 0))
00849 || (((fSizes->GetSize() == 3) || (fSizes->GetSize() == 4)) && (layer >= fSizes->GetAt(2)))) {
00850 Warning("ReadAsMatrix", "layer out of bounds.");
00851 return 0;
00852 }
00853
00854 Int_t width, height;
00855 UInt_t pixels_per_layer;
00856 Int_t offset;
00857 register UInt_t i;
00858 TMatrixD *mat=0;
00859 double *layer_pixels=0;
00860
00861 width = Int_t(fSizes->GetAt(0));
00862 height = Int_t(fSizes->GetAt(1));
00863
00864 pixels_per_layer = UInt_t(width) * UInt_t(height);
00865 offset = layer * pixels_per_layer;
00866
00867 if ((opt[0] == 'S') || (opt[0] == 's')) {
00868
00869
00870 Double_t factor, maxval=0, minval=0;
00871 Double_t pixvalue;
00872 for (i = 0; i < pixels_per_layer; i++) {
00873 pixvalue = fPixels->GetAt(offset + i);
00874
00875 if (pixvalue > maxval) {
00876 maxval = pixvalue;
00877 }
00878
00879 if ((i == 0) || (pixvalue < minval)) {
00880 minval = pixvalue;
00881 }
00882 }
00883
00884 if (maxval == minval) {
00885
00886 for (i = 0; i < pixels_per_layer; i++) {
00887 layer_pixels[i] = 1.0;
00888 }
00889 } else {
00890 factor = 1.0 / (maxval-minval);
00891 mat = new TMatrixD(height, width);
00892 layer_pixels = new double[pixels_per_layer];
00893
00894 for (i = 0; i < pixels_per_layer; i++) {
00895 layer_pixels[i] = factor * (fPixels->GetAt(offset + i) - minval);
00896 }
00897 }
00898
00899 } else {
00900
00901 layer_pixels = new double[pixels_per_layer];
00902 mat = new TMatrixD(height, width);
00903
00904 for (i = 0; i < pixels_per_layer; i++) {
00905 layer_pixels[i] = fPixels->GetAt(offset + i);
00906 }
00907 }
00908
00909 mat->Use(height, width, layer_pixels);
00910
00911 return mat;
00912 }
00913
00914
00915
00916 TH1 *TFITSHDU::ReadAsHistogram()
00917 {
00918
00919
00920
00921
00922
00923
00924
00925
00926 if (fType != kImageHDU) {
00927 Warning("ReadAsHistogram", "this is not an image HDU.");
00928 return 0;
00929 }
00930
00931 TH1 *result=0;
00932
00933 if ((fSizes->GetSize() != 1) && (fSizes->GetSize() != 2) && (fSizes->GetSize() != 3)) {
00934 Warning("ReadAsHistogram", "could not convert image HDU to histogram because it has %d dimensions.", fSizes->GetSize());
00935 return 0;
00936 }
00937
00938 if (fSizes->GetSize() == 1) {
00939
00940 UInt_t Nx = UInt_t(fSizes->GetAt(0));
00941 UInt_t x;
00942
00943 TH1D *h = new TH1D("", "", Int_t(Nx), 0, Nx-1);
00944
00945 for (x = 0; x < Nx; x++) {
00946 Int_t nentries = Int_t(fPixels->GetAt(x));
00947 if (nentries < 0) nentries = 0;
00948 h->Fill(x, nentries);
00949 }
00950
00951 result = h;
00952
00953 } else if (fSizes->GetSize() == 2) {
00954
00955 UInt_t Nx = UInt_t(fSizes->GetAt(0));
00956 UInt_t Ny = UInt_t(fSizes->GetAt(1));
00957 UInt_t x,y;
00958
00959 TH2D *h = new TH2D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1);
00960
00961 for (y = 0; y < Ny; y++) {
00962 UInt_t offset = y * Nx;
00963 for (x = 0; x < Nx; x++) {
00964 Int_t nentries = Int_t(fPixels->GetAt(offset + x));
00965 if (nentries < 0) nentries = 0;
00966 h->Fill(x,y, nentries);
00967 }
00968 }
00969
00970 result = h;
00971
00972 } else if (fSizes->GetSize() == 3) {
00973
00974 UInt_t Nx = UInt_t(fSizes->GetAt(0));
00975 UInt_t Ny = UInt_t(fSizes->GetAt(1));
00976 UInt_t Nz = UInt_t(fSizes->GetAt(2));
00977 UInt_t x,y,z;
00978
00979 TH3D *h = new TH3D("", "", Int_t(Nx), 0, Nx-1, Int_t(Ny), 0, Ny-1, Int_t(Nz), 0, Nz-1);
00980
00981
00982 for (z = 0; z < Nz; z++) {
00983 UInt_t offset1 = z * Nx * Ny;
00984 for (y = 0; y < Ny; y++) {
00985 UInt_t offset2 = y * Nx;
00986 for (x = 0; x < Nx; x++) {
00987 Int_t nentries = Int_t(fPixels->GetAt(offset1 + offset2 + x));
00988 if (nentries < 0) nentries = 0;
00989 h->Fill(x, y, z, nentries);
00990 }
00991 }
00992 }
00993
00994 result = h;
00995 }
00996
00997 return result;
00998 }
00999
01000
01001 TVectorD* TFITSHDU::GetArrayRow(UInt_t row)
01002 {
01003
01004
01005 if (fType != kImageHDU) {
01006 Warning("GetArrayRow", "this is not an image HDU.");
01007 return 0;
01008 }
01009
01010 if (fSizes->GetSize() != 2) {
01011 Warning("GetArrayRow", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
01012 return 0;
01013 }
01014
01015 UInt_t i, offset, W,H;
01016
01017 W = UInt_t(fSizes->GetAt(0));
01018 H = UInt_t(fSizes->GetAt(1));
01019
01020
01021 if (row >= H) {
01022 Warning("GetArrayRow", "index out of bounds.");
01023 return 0;
01024 }
01025
01026 offset = W * row;
01027 double *v = new double[W];
01028
01029 for (i = 0; i < W; i++) {
01030 v[i] = fPixels->GetAt(offset+i);
01031 }
01032
01033 TVectorD *vec = new TVectorD(W, v);
01034
01035 delete [] v;
01036
01037 return vec;
01038 }
01039
01040
01041 TVectorD* TFITSHDU::GetArrayColumn(UInt_t col)
01042 {
01043
01044
01045 if (fType != kImageHDU) {
01046 Warning("GetArrayColumn", "this is not an image HDU.");
01047 return 0;
01048 }
01049
01050 if (fSizes->GetSize() != 2) {
01051 Warning("GetArrayColumn", "could not get row from HDU because it has %d dimensions.", fSizes->GetSize());
01052 return 0;
01053 }
01054
01055 UInt_t i, W,H;
01056
01057 W = UInt_t(fSizes->GetAt(0));
01058 H = UInt_t(fSizes->GetAt(1));
01059
01060
01061 if (col >= W) {
01062 Warning("GetArrayColumn", "index out of bounds.");
01063 return 0;
01064 }
01065
01066 double *v = new double[H];
01067
01068 for (i = 0; i < H; i++) {
01069 v[i] = fPixels->GetAt(W*i+col);
01070 }
01071
01072 TVectorD *vec = new TVectorD(H, v);
01073
01074 delete [] v;
01075
01076 return vec;
01077 }
01078
01079
01080
01081 Int_t TFITSHDU::GetColumnNumber(const char *colname)
01082 {
01083
01084
01085 Int_t colnum;
01086 for (colnum = 0; colnum < fNColumns; colnum++) {
01087 if (fColumnsInfo[colnum].fName == colname) {
01088 return colnum;
01089 }
01090 }
01091 return -1;
01092 }
01093
01094
01095 TObjArray* TFITSHDU::GetTabStringColumn(Int_t colnum)
01096 {
01097
01098
01099 if (fType != kTableHDU) {
01100 Warning("GetTabStringColumn", "this is not a table HDU.");
01101 return 0;
01102 }
01103
01104 if ((colnum < 0) || (colnum >= fNColumns)) {
01105 Warning("GetTabStringColumn", "column index out of bounds.");
01106 return 0;
01107 }
01108
01109 if (fColumnsInfo[colnum].fType != kString) {
01110 Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
01111 return 0;
01112 }
01113
01114 Int_t offset = colnum * fNRows;
01115
01116 TObjArray *res = new TObjArray();
01117 for (Int_t row = 0; row < fNRows; row++) {
01118 res->Add(new TObjString(fCells[offset + row].fString));
01119 }
01120
01121 return res;
01122 }
01123
01124
01125 TObjArray* TFITSHDU::GetTabStringColumn(const char *colname)
01126 {
01127
01128
01129 if (fType != kTableHDU) {
01130 Warning("GetTabStringColumn", "this is not a table HDU.");
01131 return 0;
01132 }
01133
01134
01135 Int_t colnum = GetColumnNumber(colname);
01136
01137 if (colnum == -1) {
01138 Warning("GetTabStringColumn", "column not found.");
01139 return 0;
01140 }
01141
01142 if (fColumnsInfo[colnum].fType != kString) {
01143 Warning("GetTabStringColumn", "attempting to read a column that is not of type 'kString'.");
01144 return 0;
01145 }
01146
01147 Int_t offset = colnum * fNRows;
01148
01149 TObjArray *res = new TObjArray();
01150 for (Int_t row = 0; row < fNRows; row++) {
01151 res->Add(new TObjString(fCells[offset + row].fString));
01152 }
01153
01154 return res;
01155 }
01156
01157
01158 TVectorD* TFITSHDU::GetTabRealVectorColumn(Int_t colnum)
01159 {
01160
01161
01162 if (fType != kTableHDU) {
01163 Warning("GetTabRealVectorColumn", "this is not a table HDU.");
01164 return 0;
01165 }
01166
01167 if ((colnum < 0) || (colnum >= fNColumns)) {
01168 Warning("GetTabRealVectorColumn", "column index out of bounds.");
01169 return 0;
01170 }
01171
01172 if (fColumnsInfo[colnum].fType != kRealNumber) {
01173 Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
01174 return 0;
01175 }
01176
01177 Int_t offset = colnum * fNRows;
01178
01179 Double_t *arr = new Double_t[fNRows];
01180
01181 for (Int_t row = 0; row < fNRows; row++) {
01182 arr[row] = fCells[offset + row].fRealNumber;
01183 }
01184
01185 TVectorD *res = new TVectorD();
01186 res->Use(fNRows, arr);
01187
01188 return res;
01189 }
01190
01191
01192 TVectorD* TFITSHDU::GetTabRealVectorColumn(const char *colname)
01193 {
01194
01195
01196 if (fType != kTableHDU) {
01197 Warning("GetTabRealVectorColumn", "this is not a table HDU.");
01198 return 0;
01199 }
01200
01201 Int_t colnum = GetColumnNumber(colname);
01202
01203 if (colnum == -1) {
01204 Warning("GetTabRealVectorColumn", "column not found.");
01205 return 0;
01206 }
01207
01208 if (fColumnsInfo[colnum].fType != kRealNumber) {
01209 Warning("GetTabRealVectorColumn", "attempting to read a column that is not of type 'kRealNumber'.");
01210 return 0;
01211 }
01212
01213 Int_t offset = colnum * fNRows;
01214
01215 Double_t *arr = new Double_t[fNRows];
01216
01217 for (Int_t row = 0; row < fNRows; row++) {
01218 arr[row] = fCells[offset + row].fRealNumber;
01219 }
01220
01221 TVectorD *res = new TVectorD();
01222 res->Use(fNRows, arr);
01223
01224 return res;
01225 }
01226
01227
01228 Bool_t TFITSHDU::Change(const char *filter)
01229 {
01230
01231
01232
01233
01234
01235
01236
01237 TString tmppath;
01238 tmppath.Form("%s%s", fBaseFilePath.Data(), filter);
01239
01240 _release_resources();
01241 _initialize_me();
01242
01243 if (kFALSE == LoadHDU(tmppath)) {
01244
01245 Warning("Change", "error changing HDU. Restoring the previous one...");
01246
01247 _release_resources();
01248 _initialize_me();
01249
01250 if (kFALSE == LoadHDU(fFilePath)) {
01251 Warning("Change", "could not restore previous HDU. This object is no longer reliable!!");
01252 }
01253 return kFALSE;
01254 }
01255
01256
01257 fFilePath = tmppath;
01258 return kTRUE;
01259 }
01260
01261
01262 Bool_t TFITSHDU::Change(Int_t extension_number)
01263 {
01264
01265
01266 TString tmppath;
01267 tmppath.Form("[%d]", extension_number);
01268
01269 return Change(tmppath.Data());
01270 }
01271
01272
01273
01274
01275
01276