// Open a FITS file whose primary array represents // a spectrum (flux vs wavelength) void FITS_tutorial5() { TVectorD *v; printf("\n\n--------------------------------\n"); printf("WELCOME TO FITS tutorial #5 !!!!\n"); printf("--------------------------------\n"); printf("We're gonna open a FITS file that contains a\n"); printf("table with 9 rows and 8 columns. Column 4 has name\n"); printf("'mag' and contains a vector of 6 numeric components.\n"); printf("The values of vectors in rows 1 and 2 (column 4) are:\n"); printf("Row1: (99.0, 24.768, 23.215, 21.68, 21.076, 20.857)\n"); printf("Row2: (99.0, 21.689, 20.206, 18.86, 18.32 , 18.128 )\n"); printf("WARNING: when coding, row and column indices start from 0\n"); if (!gROOT->IsBatch()) { printf("Press ENTER to start..."); getchar(); printf("\n"); } //Open the table TFITSHDU *hdu = new TFITSHDU("sample4.fits[1]"); if (hdu == 0) { printf("ERROR: could not access the HDU\n"); return; } //Read vectors at rows 1 and 2 (indices 0 and 1) TVectorD *vecs[2]; vecs[0] = hdu->GetTabRealVectorCell(0, "mag"); vecs[1] = hdu->GetTabRealVectorCell(1, "mag"); for (int iVec=0; iVec < 2; iVec++) { printf("Vector %d = (", iVec+1); v = vecs[iVec]; for(int i=0; i < v->GetNoElements(); i++) { if (i>0) printf(", "); printf("%lg", (*v)[i]); //NOTE: the asterisk is for using the overloaded [] operator of the TVectorD object } printf(")\n"); } printf("\nBONUS EXAMPLE: we're gonna dump all rows using\n"); printf("the function GetTabRealVectorCells()\n"); printf("Press ENTER to continue..."); getchar(); TObjArray *vectorCollection = hdu->GetTabRealVectorCells("mag"); for (int iVec=0; iVec < vectorCollection->GetEntriesFast(); iVec++) { printf("Vector %d = (", iVec+1); v = (TVectorD *) (*vectorCollection)[iVec]; //NOTE: the asterisk is for using the overloaded [] operator of the TObjArray object for(int i=0; i < v->GetNoElements(); i++) { if (i>0) printf(", "); printf("%lg", (*v)[i]); //NOTE: the asterisk is for using the overloaded [] operator of the TVectorD object } printf(")\n"); } //Clean up delete vecs[0]; delete vecs[1]; delete vectorCollection; delete hdu; } FITS_tutorial5.C:1 FITS_tutorial5.C:2 FITS_tutorial5.C:3 FITS_tutorial5.C:4 FITS_tutorial5.C:5 FITS_tutorial5.C:6 FITS_tutorial5.C:7 FITS_tutorial5.C:8 FITS_tutorial5.C:9 FITS_tutorial5.C:10 FITS_tutorial5.C:11 FITS_tutorial5.C:12 FITS_tutorial5.C:13 FITS_tutorial5.C:14 FITS_tutorial5.C:15 FITS_tutorial5.C:16 FITS_tutorial5.C:17 FITS_tutorial5.C:18 FITS_tutorial5.C:19 FITS_tutorial5.C:20 FITS_tutorial5.C:21 FITS_tutorial5.C:22 FITS_tutorial5.C:23 FITS_tutorial5.C:24 FITS_tutorial5.C:25 FITS_tutorial5.C:26 FITS_tutorial5.C:27 FITS_tutorial5.C:28 FITS_tutorial5.C:29 FITS_tutorial5.C:30 FITS_tutorial5.C:31 FITS_tutorial5.C:32 FITS_tutorial5.C:33 FITS_tutorial5.C:34 FITS_tutorial5.C:35 FITS_tutorial5.C:36 FITS_tutorial5.C:37 FITS_tutorial5.C:38 FITS_tutorial5.C:39 FITS_tutorial5.C:40 FITS_tutorial5.C:41 FITS_tutorial5.C:42 FITS_tutorial5.C:43 FITS_tutorial5.C:44 FITS_tutorial5.C:45 FITS_tutorial5.C:46 FITS_tutorial5.C:47 FITS_tutorial5.C:48 FITS_tutorial5.C:49 FITS_tutorial5.C:50 FITS_tutorial5.C:51 FITS_tutorial5.C:52 FITS_tutorial5.C:53 FITS_tutorial5.C:54 FITS_tutorial5.C:55 FITS_tutorial5.C:56 FITS_tutorial5.C:57 FITS_tutorial5.C:58 FITS_tutorial5.C:59 FITS_tutorial5.C:60 FITS_tutorial5.C:61 FITS_tutorial5.C:62 FITS_tutorial5.C:63 FITS_tutorial5.C:64 FITS_tutorial5.C:65 FITS_tutorial5.C:66 FITS_tutorial5.C:67 FITS_tutorial5.C:68 FITS_tutorial5.C:69 |
|