source: fact/tools/rootmacros/calscope_average.C@ 14098

Last change on this file since 14098 was 12172, checked in by kraehenb, 13 years ago
File size: 4.3 KB
Line 
1#include <TROOT.h>
2#include <TCanvas.h>
3#include <TProfile.h>
4#include <TH2.h>
5#include <stdint.h>
6#include <cstdio>
7
8#include "fits.h"
9#include "FOpenDataFile.c"
10#include "FOpenCalibFile.c"
11
12int calscope_average(const char *name, const char *drsname, UInt_t maxevents)
13{
14//******************************************************************************
15//Read a datafile and plot the average over all pixels of the first DRS-calibrated events
16//ATTENTION: only works for ROI=1024
17//(array indices of the calibration wrong otherwise)
18//Example call in ROOT:
19//root [74] .x calscope_average.C++("/loc_data/rawdata/2011/09/10/20110910_020.fits","/loc_data/rawdata/2011/09/09/20110909_043.drs.fits",10)
20//T. Krähenbühl, September 2011, tpk@phys.ethz.ch
21//******************************************************************************
22
23 gROOT->SetStyle("Plain");
24
25//-------------------------------------------
26//Open the file
27//-------------------------------------------
28 fits datafile(name);
29 if (!datafile)
30 {
31 cout << "Couldn't properly open the datafile." << endl;
32 return 1;
33 }
34
35//-------------------------------------------
36//Get the data
37//-------------------------------------------
38 vector<int16_t> data;
39 vector<int16_t> data_offset;
40 unsigned int data_num;
41 size_t data_n;
42 UInt_t data_px;
43 UInt_t data_roi;
44 FOpenDataFile(datafile, data, data_offset, data_num, data_n, data_roi, data_px);
45
46//-------------------------------------------
47//Get the DRS calibration
48//-------------------------------------------
49 size_t drs_n;
50 vector<float> drs_basemean;
51 vector<float> drs_gainmean;
52 vector<float> drs_triggeroffsetmean;
53 FOpenCalibFile(drsname, drs_basemean, drs_gainmean, drs_triggeroffsetmean, drs_n);
54
55//-------------------------------------------
56//Check the sizes of the data columns
57//-------------------------------------------
58 if(drs_n!=data_n)
59 {
60 cout << "Data and DRS file incompatible (Px*ROI disagree)" << endl;
61 return 1;
62 }
63
64//-------------------------------------------
65//Create the canvas, plot and title
66//-------------------------------------------
67 char title[500];
68 std::sprintf(title,"Data: %s, DRS: %s",name,drsname);
69 TCanvas *canv = new TCanvas( "canv", "Mean values of the first events", 100, 10, 700, 500 );
70 TProfile *pix = new TProfile("pix", title, 1024, -0.5, 1023.5);
71 pix->GetXaxis()->SetTitle("DRS bin (@2 GHz)");
72 pix->GetYaxis()->SetTitle("Amplitude (mV)");
73
74// char ptitle[500];
75// std::sprintf(ptitle,"Data: %s, DRS: %s, Px %i",name,drsname);
76// TCanvas *pcanv = new TCanvas( "pcanv", title, 800, 10, 700, 500 );
77// TH2F *pers = new TH2F("pers",ptitle,1024,-0.5,1023.5,1550,-49.5,1500.5);
78// pers->GetXaxis()->SetTitle("DRS bin (@2 GHz)");
79// pers->GetYaxis()->SetTitle("Amplitude (mV)");
80
81//-------------------------------------------
82//Iterate over the events
83//-------------------------------------------
84if((maxevents==0)||(maxevents>datafile.GetNumRows())) maxevents = datafile.GetNumRows();
85 for (size_t i=0; i<maxevents; i++)
86 {
87 datafile.GetRow(i);
88 cout << "Event number: " << data_num << endl;
89
90//-------------------------------------------
91//Iterate over the pixels
92//-------------------------------------------
93 for (int j=0; j<data_px; j++)
94 {
95
96//-------------------------------------------
97//Iterate over the slices
98//-------------------------------------------
99 for (UInt_t k=0; k<data_roi; k++)
100 {
101 UInt_t drs_calib_offset = (k+data_offset[j])%data_roi;
102 //Example values for the various datasets (20110803_015 / 018?9)
103 //data = -1856 +- 30.76
104 //drs_basemean = -906.6 +- 14.73
105 //drs_gainmean = 1696 +- 4
106 //drs_triggeroffsetmean = -0.01 +- 1.532
107
108 //Gain calibration to mV: 50000(DAC) / drs_gainmean(ADC) * 2500(mV) / 65536(DAC)
109 //Target value: 1.907 V für DRS-Calib files beim DAC-Wert 50000
110 float sample = (data[j*data_roi+k]*2000/4096.-drs_basemean[j*data_roi+drs_calib_offset]-drs_triggeroffsetmean[j*data_roi+k])/drs_gainmean[j*data_roi+drs_calib_offset]*1907.35;
111 pix->Fill(k,sample);
112// pers->Fill(k,sample);
113 }
114 }
115 }
116
117//-------------------------------------------
118//Draw the data
119//-------------------------------------------
120 canv->cd();
121 pix->Draw();
122 canv->Modified();
123 canv->Update();
124// pcanv->cd();
125// pers->Draw("col");
126// pcanv->Modified();
127// pcanv->Update();
128 return 0;
129}
Note: See TracBrowser for help on using the repository browser.