source: fact/tools/rootmacros/calscope_persistent.C@ 13182

Last change on this file since 13182 was 12166, checked in by neise, 13 years ago
initial commit
File size: 5.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#include "fits.h"
8
9int calscope_persistent(const char *name, const char *drsname, size_t pixelnr, UInt_t maxevents)
10{
11//******************************************************************************
12//Read a datafile and plot the DRS-calibrated data
13//ATTENTION: only works for ROI=1024
14// (array indices of the calibration wrong otherwise)
15//Example call in ROOT:
16//root [74] .x calscope.C++("20110804_024.fits","20110804_023.drs.fits",10,1348)
17//T. Krähenbühl, August 2011, tpk@phys.ethz.ch
18//******************************************************************************
19
20 gROOT->SetStyle("Plain");
21
22//-------------------------------------------
23//Open the files
24//-------------------------------------------
25 fits datafile(name);
26 fits drsfile(drsname);
27
28 if (!datafile)
29 {
30 cout << "Couldn't properly open the datafile." << endl;
31 return 1;
32 }
33
34 if (!drsfile)
35 {
36 cout << "Couldn't properly open the drsfile." << endl;
37 return 1;
38 }
39
40//-------------------------------------------
41//Print the headers
42//-------------------------------------------
43 cout << "-------------------- Data Header -------------------" << endl;
44 datafile.PrintKeys();
45 cout << "----------------- DRS Calib Header -----------------" << endl;
46 drsfile.PrintKeys();
47 cout << "------------------- Data Columns -------------------" << endl;
48 datafile.PrintColumns();
49 cout << "---------------- DRS Calib Columns -----------------" << endl;
50 drsfile.PrintColumns();
51 cout << "--------------------- Data --------------------" << endl;
52
53//-------------------------------------------
54//Get the sizes of the data columns
55//-------------------------------------------
56 const size_t data_n = datafile.GetN("Data"); //Size of column "Data" = #Pixel x ROI
57 const size_t drs_n = drsfile.GetN("BaselineMean");
58 if(drs_n!=data_n)
59 {
60 cout << "Data and DRS file incompatible (Px*ROI disagree)" << endl;
61 return 1;
62 }
63
64 const UInt_t data_roi = datafile.GetUInt("NROI"); // Value from header
65 const UInt_t data_px = datafile.GetUInt("NPIX");
66
67//-------------------------------------------
68//Link the data to variables
69//-------------------------------------------
70 unsigned int data_num;
71 datafile.SetRefAddress("EventNum", data_num);
72 vector<int16_t> data(data_n);
73 datafile.SetVecAddress("Data", data);
74 vector<int16_t> dataOffset(data_px);
75 datafile.SetVecAddress("StartCellData", dataOffset);
76
77 vector<float> drs_basemean(drs_n);
78 drsfile.SetVecAddress("BaselineMean", drs_basemean);
79 vector<float> drs_gainmean(drs_n);
80 drsfile.SetVecAddress("GainMean", drs_gainmean);
81 vector<float> drs_triggeroffsetmean(drs_n);
82 drsfile.SetVecAddress("TriggerOffsetMean", drs_triggeroffsetmean);
83 drsfile.GetRow(0); //Read the calibration data
84
85//-------------------------------------------
86//Create the canvas, plot and title
87//-------------------------------------------
88 char title[500];
89 std::sprintf(title,"Data: %s, DRS: %s, Px %i",name,drsname,pixelnr);
90 TCanvas *canv = new TCanvas( "canv", "Mean values of the first event", 100, 10, 700, 500 );
91 TProfile *pix = new TProfile("pix", title, 1024, -0.5, 1023.5);
92 pix->GetXaxis()->SetTitle("DRS bin (@2 GHz)");
93 pix->GetYaxis()->SetTitle("Amplitude (mV)");
94
95 char ptitle[500];
96 std::sprintf(ptitle,"Data: %s, DRS: %s, Px %i",name,drsname);
97 TCanvas *pcanv = new TCanvas( "pcanv", title, 800, 10, 700, 500 );
98 TH2F *pers = new TH2F("pers",ptitle,1024,-0.5,1023.5,200,-49.5,150.5);
99 pers->GetXaxis()->SetTitle("DRS bin (@2 GHz)");
100 pers->GetYaxis()->SetTitle("Amplitude (mV)");
101
102//-------------------------------------------
103//Iterate over the events
104//-------------------------------------------
105if(maxevents==0) maxevents = datafile.GetNumRows();
106// for (size_t i=0; i<datafile.GetNumRows(); i++)
107 for (size_t i=0; i<maxevents; i++)
108 {
109// size_t i=eventnr; //Fix the Event
110 datafile.GetRow(i);
111 cout << "Event number: " << data_num << endl;
112
113//-------------------------------------------
114//Iterate over the pixels
115//-------------------------------------------
116// for (int j=0; j<data_px; j++)
117// {
118 size_t j=pixelnr; //Fix the Pixel to a SoftID
119
120//-------------------------------------------
121//Iterate over the slices
122//-------------------------------------------
123 for (UInt_t k=0; k<data_roi; k++)
124 {
125 UInt_t drs_calib_offset = (k+dataOffset[j])%data_roi;
126 //Example values for the various datasets (20110803_015 / 018?9)
127 //data = -1856 +- 30.76
128 //drs_basemean = -906.6 +- 14.73
129 //drs_gainmean = 1696 +- 4
130 //drs_triggeroffsetmean = -0.01 +- 1.532
131
132 //Gain calibration to mV: 50000(DAC) / drs_gainmean(ADC) * 2500(mV) / 65536(DAC)
133 //Target value: 1.907 V für DRS-Calib files beim DAC-Wert 50000
134 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;
135 pix->Fill(k,sample);
136 pers->Fill(k,sample);
137 }
138// }
139 }
140
141//-------------------------------------------
142//Draw the data
143//-------------------------------------------
144 canv->cd();
145 pix->Draw();
146 canv->Modified();
147 canv->Update();
148 pcanv->cd();
149 pers->Draw("col");
150 pcanv->Modified();
151 pcanv->Update();
152 return 0;
153}
Note: See TracBrowser for help on using the repository browser.