source: fact/tools/rootmacros/fpeak_cfd.C@ 12502

Last change on this file since 12502 was 12502, checked in by neise, 13 years ago
runnable with just the first 3 parameters
  • Property svn:executable set to *
File size: 10.2 KB
Line 
1#include <TROOT.h>
2#include <TCanvas.h>
3#include <TProfile.h>
4#include <TTimer.h>
5#include <TH1F.h>
6#include <TH2F.h>
7#include <Getline.h>
8#include <TLine.h>
9#include <TBox.h>
10#include <TMath.h>
11#include <TFile.h>
12#include <TStyle.h>
13
14#include <stdio.h>
15#include <stdint.h>
16#include <cstdio>
17#include <deque>
18
19#define NPIX 1440
20#define NCELL 1024
21#define FAD_MAX_SAMPLES 1024
22
23#define HAVE_ZLIB
24#include "fits.h"
25#include "FOpenCalibFile.c"
26
27#include "discriminator.h"
28#include "discriminator.C"
29#include "zerosearch.h"
30#include "zerosearch.C"
31#include "factfir.C"
32
33#include "FOpenDataFile.h"
34#include "FOpenDataFile.c"
35
36#include "DrsCalibration.C"
37#include "DrsCalibration.h"
38
39#include "SpikeRemoval.h"
40#include "SpikeRemoval.C"
41
42bool breakout=false;
43
44int NEvents;
45vector<int16_t> AllPixelDataVector;
46vector<int16_t> StartCellVector;
47unsigned int CurrentEventID;
48size_t PXLxROI;
49UInt_t RegionOfInterest;
50UInt_t NumberOfPixels;
51
52size_t drs_n;
53vector<float> drs_basemean;
54vector<float> drs_gainmean;
55vector<float> drs_triggeroffsetmean;
56
57vector<float> Ameas(FAD_MAX_SAMPLES); // copy of the data (measured amplitude
58vector<float> Vcorr(FAD_MAX_SAMPLES); // corrected Values
59vector<float> Vslide(FAD_MAX_SAMPLES); // sliding average result
60vector<float> Vcfd(FAD_MAX_SAMPLES); // CDF result
61vector<float> Vcfd2(FAD_MAX_SAMPLES); // CDF result + 2nd sliding average
62
63// histograms
64const int NumberOfDebugHistoTypes = 7;
65const unsigned int
66 Ameas_ = 0,
67 N1mean_ = 1,
68 Vcorr_ = 2,
69 Vtest_ = 3,
70 Vslide_ = 4,
71 Vcfd_ = 5,
72 Vcfd2_ = 6;
73
74TH1F* h;
75TH1F *debugHistos;
76TH2F* hStartCell; // id of the DRS physical pipeline cell where readout starts
77 // x = pixel id, y = DRS cell id
78
79TH2F * hAmplSpek_cfd;
80TObjArray hList;
81
82void BookHistos( TString &title );
83void SaveHistograms(const char * );
84
85///////////////////////////////////////////////////////////////////////////////
86///////////////////////////////////////////////////////////////////////////////
87///////////////////////////////////////////////////////////////////////////////
88// read FACT raw data
89// * remove spikes
90// * calculate baseline
91// * find peaks (CFD and normal discriminator)
92// * compute pulse height and pulse integral spektrum of the peaks
93int fpeak_cfd(
94 char *datafilename = "data/20111016_013.fits.gz",
95 const char *drsfilename = "../../20111016_011.drs.fits.gz",
96 const char *OutRootFileName = "../analysis/fpeak_cdf.Coutput.root",
97 int firstevent = 0,
98 int nevents = -1,
99 int firstpixel = 0,
100 int npixel = -1,
101 bool spikeDebug = false,
102 int avg1 = 8,
103 int avg2 = 8,
104 int verbosityLevel = 1, // different verbosity levels can be implemented here
105 bool ProduceGraphic = true
106 )
107{
108gStyle->SetPalette(1,0);
109gROOT->SetStyle("Plain");
110
111TString histogramtitle;
112histogramtitle += datafilename;
113histogramtitle += " ";
114histogramtitle += drsfilename ;
115histogramtitle += " ";
116
117 histogramtitle += avg1;
118 histogramtitle += " ";
119 histogramtitle += avg2;
120// histogramtitle += " 3.0 100";
121 cout << histogramtitle.Data() << endl;
122
123
124 TCanvas *cFiltered = NULL;
125 if (spikeDebug){
126 cFiltered = new TCanvas("cFiltered","filtered DRS Waveforms",410,410,400,400);
127 cFiltered->Divide(1, 4);
128 }
129
130
131
132// CFD filter settings
133 int k_cfd = 10;
134 vector<double> a_cfd(k_cfd, 0);
135 double b_cfd = 1.;
136 a_cfd[0]=-0.75;
137 a_cfd[k_cfd-1]=1.;
138
139 fits * datafile;
140 // Opens the raw data file and 'binds' the variables given as
141 // Parameters to the data file. So they are filled with
142 // raw data as soon as datafile->GetRow(int) is called.
143 NEvents = OpenDataFile( datafilename, &datafile,
144 AllPixelDataVector, StartCellVector, CurrentEventID,
145 RegionOfInterest, NumberOfPixels, PXLxROI, verbosityLevel);
146 if (NEvents == 0){
147 cout << "return code of OpenDataFile:" << datafilename<< endl;
148 cout << "is zero -> aborting." << endl;
149 return 1;
150 }
151
152 if (verbosityLevel > 0)
153 cout <<"number of events in file: "<< NEvents << "\t";
154 if ( nevents == -1 || nevents > NEvents ) nevents = NEvents; // -1 means all!
155 if (verbosityLevel > 0)
156 cout <<"of, which "<< nevents << "will be processed"<< endl;
157
158 if (verbosityLevel > 0)
159 cout <<"Total # of Pixel: "<< NumberOfPixels << "\t";
160 if ( npixel == -1 || npixel > (int)NumberOfPixels ) npixel = (int)NumberOfPixels; // -1 means all!
161 if (verbosityLevel > 0)
162 cout <<"of, which "<< npixel << "will be processed"<< endl;
163
164
165
166 //Get the DRS calibration
167 FOpenCalibFile( drsfilename,
168 drs_basemean,
169 drs_gainmean,
170 drs_triggeroffsetmean,
171 drs_n);
172
173 BookHistos( histogramtitle );
174//-----------------------------------------------------------------------------
175// Loops over Every Event and Pixel
176//-----------------------------------------------------------------------------
177 for ( int ev = firstevent; ev < firstevent + nevents; ev++) {
178 // Get an Event --> consists of 1440 Pixel ...erm....data
179 datafile->GetRow( ev );
180
181 for ( int pix = firstpixel; pix < firstpixel+npixel; pix++ ){
182 if (verbosityLevel > 0){
183 if (pix == firstpixel){
184 cout << "Processing Event: " << CurrentEventID << "/" << nevents << endl;
185 }
186 }
187
188 applyDrsCalibration( Ameas,pix,12,12,
189 drs_basemean, drs_gainmean, drs_triggeroffsetmean,
190 RegionOfInterest, AllPixelDataVector, StartCellVector);
191
192 // finds spikes in the raw data, and interpolates the value
193 // spikes are: 1 or 2 slice wide, positive non physical artifacts
194 removeSpikes (Ameas, Vcorr);
195
196 // filter Vcorr with sliding average using FIR filter function
197 sliding_avg(Vcorr, Vslide, avg1);
198 // filter Vslide with CFD using FIR filter function
199 factfir(b_cfd , a_cfd, k_cfd, Vslide, Vcfd);
200 // filter Vcfd with sliding average using FIR filter function
201 sliding_avg(Vcfd, Vcfd2, avg2);
202
203
204 // peaks in Ameas[] are found by searching for zero crossings
205 // in Vcfd2
206 // first Argument 1 means ... *rising* edge
207 // second Argument 1 means ... search with stepsize 1 ... 8 is okay as well
208 vector<Region> * zXings = zerosearch( Vcfd2 , 1 , 8);
209 // zXings means "zero cross ings"
210 EnlargeRegion( *zXings, 10, 10);
211 findAbsMaxInRegions( *zXings, Vslide);
212 // removeMaximaBelow( *zXings, 3.0, 0);
213 // removeRegionWithMaxOnEdge( *zXings, 2);
214 // removeRegionOnFallingEdge( *zXings, 100);
215
216
217 // fill maxima in Histogram
218 if (zXings->size() != 0 ){
219 for (unsigned int i=0; i<zXings->size(); i++){
220 if (verbosityLevel > 1){
221 cout << zXings->at(i).maxPos << ":\t"
222 << zXings->at(i).maxVal <<endl;
223 }
224 hAmplSpek_cfd->Fill(pix, zXings->at(i).maxVal);
225 }
226 }
227
228 if ( spikeDebug ){
229
230 // TODO do this correct. The vectors should be the rigt ones... this is just luck
231 debugHistos[Ameas_].GetXaxis()->Set(Ameas.size() , -0.5 , Ameas.size()-0.5);
232 debugHistos[Vcorr_].GetXaxis()->Set(Ameas.size() , -0.5 , Ameas.size()-0.5);
233 debugHistos[Vslide_].GetXaxis()->Set(Ameas.size() , -0.5 , Ameas.size()-0.5);
234 debugHistos[Vcfd_].GetXaxis()->Set(Ameas.size() , -0.5 , Ameas.size()-0.5);
235 debugHistos[Vcfd2_].GetXaxis()->Set(Ameas.size() , -0.5 , Ameas.size()-0.5);
236
237 for ( unsigned int sl = 0; sl < RegionOfInterest; sl++){
238 debugHistos[Ameas_].SetBinContent(sl, Ameas[sl]);
239 debugHistos[Vcorr_].SetBinContent(sl, Vcorr[sl]);
240 debugHistos[Vslide_].SetBinContent( sl, Vslide[sl] );
241 debugHistos[Vcfd_].SetBinContent( sl, Vcfd[sl] );
242 debugHistos[Vcfd2_].SetBinContent( sl, Vcfd2[sl] );
243 }
244
245 cFiltered->cd(1);
246 gPad->SetGrid();
247 debugHistos[Ameas_].Draw();
248
249 cFiltered->cd(2);
250 gPad->SetGrid();
251 debugHistos[Vcorr_].Draw();
252
253 cFiltered->cd(3);
254 gPad->SetGrid();
255 debugHistos[Vslide_].Draw();
256
257 TBox *OneBox;
258 vector<TBox*> MyBoxes;
259 for (unsigned int i=0; i<zXings->size(); i++){
260 OneBox = new TBox(
261 zXings->at(i).maxPos -10 ,
262 zXings->at(i).maxVal -0.5,
263 zXings->at(i).maxPos +10 ,
264 zXings->at(i).maxVal +0.5);
265 OneBox->SetLineColor(kBlue);
266 OneBox->SetLineWidth(1);
267 OneBox->SetFillStyle(0);
268 OneBox->SetFillColor(kRed);
269 MyBoxes.push_back(OneBox);
270 OneBox->Draw();
271 }
272
273 cFiltered->cd(4);
274 gPad->SetGrid();
275 debugHistos[Vcfd2_].Draw();
276 TLine *zeroline = new TLine(0, 0, 1024, 0);
277 zeroline->SetLineColor(kBlue);
278 zeroline->Draw();
279
280 cFiltered->Update();
281 //Process gui events asynchronously during input
282 TTimer timer("gSystem->ProcessEvents();", 50, kFALSE);
283 timer.TurnOn();
284 TString input = Getline("Type 'q' to exit, <return> to go on: ");
285 timer.TurnOff();
286 if (input=="q\n") {
287 breakout=true;
288 }
289
290 //TODO!!!!!!!!!
291 // do some Garbage collection ...
292 }// end of if(spikeDebug)
293
294 delete zXings;
295 if (breakout)
296 break;
297
298 } // end of loop over pixels
299
300 if (breakout)
301 break;
302 } // end of loop over pixels
303if (ProduceGraphic){
304 TCanvas * cSpektrum;
305 cSpektrum = new TCanvas("cSpektrum","Amplitude Spektrum",10,10,400,400);
306 cSpektrum->Divide(1,1);
307 cSpektrum->cd(1);
308 hAmplSpek_cfd->Draw("COLZ");
309 cSpektrum->Modified();
310 cSpektrum->Update();
311}
312 SaveHistograms( OutRootFileName );
313 return( 0 );
314}
315
316// booking and parameter settings for all histos
317void BookHistos( TString& title){
318 // histograms for baseline extraction
319
320
321
322 hAmplSpek_cfd = new TH2F("hAmplSpek_cfd", title,1440,-0.5,1439.5, 256, -27.5, 100.5);
323 hAmplSpek_cfd->GetXaxis()->SetTitle( "pixel" );
324 hAmplSpek_cfd->GetYaxis()->SetTitle( "amplitude in mV" );
325 hList.Add( hAmplSpek_cfd );
326
327
328 debugHistos = new TH1F[ NumberOfDebugHistoTypes ];
329 for ( int type = 0; type < NumberOfDebugHistoTypes; type++){
330 debugHistos[ type ].SetBins(1024, 0, 1024);
331 debugHistos[ type ].SetLineColor(1);
332 debugHistos[ type ].SetLineWidth(2);
333
334 // set X axis paras
335 debugHistos[ type ].GetXaxis()->SetLabelSize(0.1);
336 debugHistos[ type ].GetXaxis()->SetTitleSize(0.1);
337 debugHistos[ type ].GetXaxis()->SetTitleOffset(1.2);
338 debugHistos[ type ].GetXaxis()->SetTitle(Form("Time slice (%.1f ns/slice)", 1./2.));
339
340 // set Y axis paras
341 debugHistos[ type ].GetYaxis()->SetLabelSize(0.1);
342 debugHistos[ type ].GetYaxis()->SetTitleSize(0.1);
343 debugHistos[ type ].GetYaxis()->SetTitleOffset(0.3);
344 debugHistos[ type ].GetYaxis()->SetTitle("Amplitude (a.u.)");
345 }
346}
347
348void SaveHistograms( const char * loc_fname){
349 // create the histogram file (replace if already existing)
350 TFile tf( loc_fname, "RECREATE");
351
352 hList.Write(); // write the major histograms into the top level directory
353
354 tf.Close(); // close the file
355}
Note: See TracBrowser for help on using the repository browser.