source: fact/tools/rootmacros/calscope.C@ 12814

Last change on this file since 12814 was 12362, checked in by kraehenb, 13 years ago
New version of the calscope-script which is easier to adapt. The included file TPKplotevent was dropped for the FCalibrateEvent.C, which takes an event and returns an array with the DRS-calibrated pipeline.
File size: 3.0 KB
Line 
1#include <TROOT.h>
2#include <TCanvas.h>
3#include <TProfile.h>
4
5#include <stdint.h>
6#include <cstdio>
7
8#define HAVE_ZLIB
9#include "fits.h"
10#include "FOpenDataFile.c"
11#include "FOpenCalibFile.c"
12#include "FCalibrateEvent.c"
13
14int calscope(const char *name, const char *drsname, size_t eventnr, size_t pixelnr)
15{
16//******************************************************************************
17//Read a datafile and plot the DRS-calibrated data
18//ATTENTION: only works for ROI=1024
19// (array indices of the calibration wrong otherwise)
20//Example call in ROOT:
21//root [74] .x calscope.C++("20110804_024.fits","20110804_023.drs.fits",10,1348)
22//T. Kr�henb�hl, August 2011, tpk@phys.ethz.ch
23//******************************************************************************
24
25 gROOT->SetStyle("Plain");
26
27//-------------------------------------------
28//Open the file
29//-------------------------------------------
30 fits datafile(name);
31 if (!datafile)
32 {
33 cout << "Couldn't properly open the datafile." << endl;
34 return 1;
35 }
36
37//-------------------------------------------
38//Get the data
39//-------------------------------------------
40 vector<int16_t> data;
41 vector<int16_t> data_offset;
42 unsigned int data_num;
43 size_t data_n;
44 UInt_t data_px;
45 UInt_t data_roi;
46 FOpenDataFile(datafile, data, data_offset, data_num, data_n, data_roi, data_px);
47
48//-------------------------------------------
49//Get the DRS calibration
50//-------------------------------------------
51 size_t drs_n;
52 vector<float> drs_basemean;
53 vector<float> drs_gainmean;
54 vector<float> drs_triggeroffsetmean;
55 FOpenCalibFile(drsname, drs_basemean, drs_gainmean, drs_triggeroffsetmean, drs_n);
56
57//-------------------------------------------
58//Check the sizes of the data columns
59//-------------------------------------------
60 if(drs_n!=data_n)
61 {
62 cout << "Data and DRS file incompatible (Px*ROI disagree)" << endl;
63 return 1;
64 }
65
66//-------------------------------------------
67//Create the title
68//-------------------------------------------
69 char title[500];
70 std::sprintf(title,"Data: %s, DRS: %s, Px %i Ev %i",name,drsname,pixelnr,eventnr);
71
72//-------------------------------------------
73//Get the calibrated event
74//-------------------------------------------
75 vector<float> calevent(data_px*data_roi); //Vector for the calibrated event
76
77 cout << "--------------------- Data --------------------" << endl;
78 datafile.GetRow(eventnr);
79 cout << "Event number: " << data_num << endl;
80
81 FCalibrateEvent(data, data_offset, drs_basemean, drs_gainmean, drs_triggeroffsetmean, calevent, data_px, data_roi);
82
83//-------------------------------------------
84//Draw the data
85//-------------------------------------------
86 TCanvas *canv = new TCanvas( "canv", "Mean values of the first event", 100, 10, 700, 500 );
87 TProfile *pix = new TProfile("pix", title, 1024, -0.5, 1023.5);
88
89 for (UInt_t k=0; k<data_roi; k++)
90 {
91 pix->Fill(k,calevent[pixelnr*data_roi+k]);
92 }
93
94 pix->Draw();
95 canv->Modified();
96 canv->Update();
97
98 return 0;
99}
Note: See TracBrowser for help on using the repository browser.