MiniballSort
Loading...
Searching...
No Matches
CDCalibrator.cc
Go to the documentation of this file.
1#include "CDCalibrator.hh"
2
3MiniballCDCalibrator::MiniballCDCalibrator( std::shared_ptr<MiniballSettings> myset ){
4
5 // First get the settings
6 set = myset;
7
8 // No calibration file by default
9 overwrite_cal = false;
10
11 // No input file at the start by default
12 flag_input_file = false;
13
14 // Progress bar starts as false
15 _prog_ = false;
16
17 // Start at MBS event 0
18 preveventid = 0;
19
20 // ------------------------------- //
21 // Initialise variables and flags //
22 // ------------------------------- //
23 build_window = set->GetEventWindow();
24
25 // Intialise the hist list
26 histlist = new TList();
27
28}
29
31
32 // Call for every new file
33 // Reset counters etc.
34
35 time_prev = 0;
36 time_min = 0;
37 time_max = 0;
38 time_first = 0;
39
40}
41
42void MiniballCDCalibrator::SetInputFile( std::vector<std::string> input_file_names ) {
43
45 input_tree = new TChain( "mb_sort" );
46 mbsinfo_tree = new TChain( "mbsinfo" );
47 for( unsigned int i = 0; i < input_file_names.size(); i++ ) {
48
49 input_tree->Add( input_file_names[i].data() );
50 mbsinfo_tree->Add( input_file_names[i].data() );
51
52 }
53
54 flag_input_file = true;
55
56 input_tree->SetBranchAddress( "data", &in_data );
57 mbsinfo_tree->SetBranchAddress( "mbsinfo", &mbs_info );
58 mbsinfo_tree->BuildIndex("GetEventID()");
59
60 return;
61
62}
63
64void MiniballCDCalibrator::SetInputFile( std::string input_file_name ) {
65
66 // Open next Root input file.
67 input_file = new TFile( input_file_name.data(), "read" );
68 if( input_file->IsZombie() ) {
69
70 std::cout << "Cannot open " << input_file_name << std::endl;
71 return;
72
73 }
74
75 flag_input_file = true;
76
77 // Set the input tree
78 SetInputTree( (TTree*)input_file->Get("mb_sort") );
79 SetMBSInfoTree( (TTree*)input_file->Get("mbsinfo") );
80 StartFile();
81
82 return;
83
84}
85
86void MiniballCDCalibrator::SetInputTree( TTree *user_tree ){
87
88 // Find the tree and set branch addresses
89 input_tree = (TChain*)user_tree;
90 in_data = nullptr;
91 input_tree->SetBranchAddress( "data", &in_data );
92
93 return;
94
95}
96
97void MiniballCDCalibrator::SetMBSInfoTree( TTree *user_tree ){
98
99 // Find the tree and set branch addresses
100 mbsinfo_tree = (TChain*)user_tree;
101 mbs_info = nullptr;
102 mbsinfo_tree->SetBranchAddress( "mbsinfo", &mbs_info );
103
104 return;
105
106}
107
108void MiniballCDCalibrator::SetOutput( std::string output_file_name, bool cWrite ) {
109
110 // ------------------------------------------------------------------------ //
111 // Create output file and create events tree
112 // ------------------------------------------------------------------------ //
113 output_file = new TFile( output_file_name.data(), "recreate" );
114
115 // Hisograms in separate function
116 MakeHists();
117
118 // Output the calibration coefficients
119 std::string cal_file_name = output_file_name.substr( 0, output_file_name.find_last_of(".") );
120 cal_file_name += ".cal";
121 output_cal.open( cal_file_name.data(), std::ios::trunc );
122
123 // Write once at the start if in spy
124 if( cWrite ) output_file->Write();
125
126}
127
129
131
132 flag_close_event = false;
133 event_open = false;
134
135 hit_ctr = 0;
136
137 std::vector<float>().swap(cd_en_list);
138 std::vector<unsigned int>().swap(cd_Q_list);
139 std::vector<unsigned long long>().swap(cd_ts_list);
140 std::vector<unsigned char>().swap(cd_det_list);
141 std::vector<unsigned char>().swap(cd_sec_list);
142 std::vector<unsigned char>().swap(cd_side_list);
143 std::vector<unsigned char>().swap(cd_strip_list);
144
145 return;
146
147}
148
149
151
152 std::string hname, htitle;
153
154 // ------------- //
155 // CD histograms //
156 // ------------- //
157 cd_pen_nen.resize( set->GetNumberOfCDDetectors() );
158 cd_nen_pen.resize( set->GetNumberOfCDDetectors() );
159 cd_pen_nQ.resize( set->GetNumberOfCDDetectors() );
160 cd_nQ_pQ.resize( set->GetNumberOfCDDetectors() );
161
162 // Get sizes and scales
163 double maxQ = 1073741824;
164 unsigned int Qbins = 8192;
165
166 if( set->GetNumberOfCaenAdcModules() > 0 ) {
167 maxQ = 4096;
168 Qbins = 4096;
169 }
170
171 else if( set->GetNumberOfFebexSfps() > 1 &&
172 set->GetNumberOfFebexBoards() > 0 &&
173 set->GetNumberOfFebexChannels() > 0 ) {
174
175 if( cal->FebexType( 1, 0, 0 ) == "Qshort" ) {
176 maxQ = 65536;
177 }
178 }
179
180 for( unsigned int i = 0; i < set->GetNumberOfCDDetectors(); ++i ) {
181
182 cd_pen_nen[i].resize( set->GetNumberOfCDSectors() );
183 cd_nen_pen[i].resize( set->GetNumberOfCDSectors() );
184 cd_pen_nQ[i].resize( set->GetNumberOfCDSectors() );
185 cd_nQ_pQ[i].resize( set->GetNumberOfCDSectors() );
186
187 for( unsigned int j = 0; j < set->GetNumberOfCDSectors(); ++j ) {
188
189 cd_nen_pen[i][j].resize( set->GetNumberOfCDPStrips() );
190 cd_nQ_pQ[i][j].resize( set->GetNumberOfCDPStrips() );
191
192 for( unsigned int k = 0; k < set->GetNumberOfCDPStrips(); ++k ) {
193
194 hname = "cd_" + std::to_string(i) + "_" + std::to_string(j);
195 hname += "_nen_" + std::to_string(ptag) + "_pen_" + std::to_string(k);
196 htitle = "CD n-side energy vs p-side energy for detector " + std::to_string(i);
197 htitle += ", sector " + std::to_string(j) + ", pid " + std::to_string(k);
198 htitle += ", nid " + std::to_string(ntag);
199 htitle += ";n-side energy (keV);p-side energy (keV);Counts";
200 cd_nen_pen[i][j][k] = new TH2F( hname.data(), htitle.data(),
201 4000, 0, 2000e3, 4000, 0, 2000e3 );
202 histlist->Add(cd_nen_pen[i][j][k]);
203
204 hname = "cd_" + std::to_string(i) + "_" + std::to_string(j);
205 hname += "_nQ_" + std::to_string(ptag) + "_pQ_" + std::to_string(k);
206 htitle = "CD n-side energy vs p-side raw charge for detector " + std::to_string(i);
207 htitle += ", sector " + std::to_string(j) + ", pid " + std::to_string(k);
208 htitle += ", nid " + std::to_string(ntag);
209 htitle += ";n-side raw charge (ADC units);p-side raw charge (ADC units);Counts";
210 cd_nQ_pQ[i][j][k] = new TH2F( hname.data(), htitle.data(),
211 Qbins, 0, maxQ, Qbins, 0, maxQ );
212 histlist->Add(cd_nQ_pQ[i][j][k]);
213
214 } // k
215
216 cd_pen_nen[i][j].resize( set->GetNumberOfCDNStrips() );
217 cd_pen_nQ[i][j].resize( set->GetNumberOfCDNStrips() );
218
219 for( unsigned int k = 0; k < set->GetNumberOfCDNStrips(); ++k ) {
220
221 hname = "cd_" + std::to_string(i) + "_" + std::to_string(j);
222 hname += "_pen_" + std::to_string(ptag) + "_nen_" + std::to_string(k);
223 htitle = "CD p-side energy vs n-side energy for detector " + std::to_string(i);
224 htitle += ", sector " + std::to_string(j) + ", pid " + std::to_string(ptag);
225 htitle += ", nid " + std::to_string(k);
226 htitle += ";p-side energy (keV);n-side energy (keV);Counts";
227 cd_pen_nen[i][j][k] = new TH2F( hname.data(), htitle.data(),
228 4000, 0, 2000e3, 4000, 0, 2000e3 );
229 histlist->Add(cd_pen_nen[i][j][k]);
230
231 hname = "cd_" + std::to_string(i) + "_" + std::to_string(j);
232 hname += "_pen_" + std::to_string(ptag) + "_nQ_" + std::to_string(k);
233 htitle = "CD p-side energy vs n-side raw charge for detector " + std::to_string(i);
234 htitle += ", sector " + std::to_string(j) + ", pid " + std::to_string(ptag);
235 htitle += ", nid " + std::to_string(k);
236 htitle += ";p-side energy (keV);n-side raw charge (ADC units);Counts";
237 cd_pen_nQ[i][j][k] = new TH2F( hname.data(), htitle.data(),
238 4000, 0, 2000e3, Qbins, 0, maxQ );
239 histlist->Add(cd_pen_nQ[i][j][k]);
240
241 } // k
242
243 } // j
244
245 } // i
246
247
248 // flag to denote that hists are ready (used for spy)
249 hists_ready = true;
250
251 return;
252
253}
254
255// Reset histograms in the DataSpy
257
258 // Loop over the hist list
259 TIter next( histlist->MakeIterator() );
260 while( TObject *obj = next() ) {
261
262 if( obj->InheritsFrom( "TH2" ) )
263 ( (TH2*)obj )->Reset("ICESM");
264 else if( obj->InheritsFrom( "TH1" ) )
265 ( (TH1*)obj )->Reset("ICESM");
266
267 }
268
269 return;
270
271}
272
273bool MiniballCDCalibrator::FindCDChannels( int det, int sec, int side, int strip, int &adc, int &ch ) {
274
275 // Loop over ADCs
276 for( unsigned int m = 0; m < set->GetNumberOfCaenAdcModules(); ++m ) {
277
278 // Loop over channels
279 for( unsigned int c = 0; c < set->GetNumberOfCaenAdcChannels(); ++c ) {
280
281 // Check that it's a CD
282 if( !set->IsCD(m,c) ) continue;
283
284 // Check we have the correct CD detector
285 if( set->GetCDDetector(m,c) != det ) continue;
286
287 // Check we have the correct sector
288 if( set->GetCDSector(m,c) != sec ) continue;
289
290 // Check we have an P side (==0)
291 if( set->GetCDSide(m,c) != side ) continue;
292
293 // Check we have the correct strip
294 if( set->GetCDStrip(m,c) != strip ) continue;
295
296 // Then we got the right channel
297 adc = m;
298 ch = c;
299 return true;
300
301 } // c
302
303 } // m
304
305 std::cerr << "CD strip not found, det=" << det << ", sec=" << sec;
306 std::cerr << ", side=" << side << ", strip=" << strip << std::endl;
307 return false;
308
309}
310
311bool MiniballCDCalibrator::FindCDChannels( int det, int sec, int side, int strip, int &sfp, int &board, int &ch ) {
312
313 // Loop over SFPs
314 for( unsigned int s = 0; s < set->GetNumberOfFebexSfps(); ++s ) {
315
316 // Loop over boards
317 for( unsigned int m = 0; m < set->GetNumberOfFebexBoards(); ++m ) {
318
319 // Loop over channels
320 for( unsigned int c = 0; c < set->GetNumberOfFebexChannels(); ++c ) {
321
322 // Check that it's a CD
323 if( !set->IsCD(s,m,c) ) continue;
324
325 // Check we have the correct CD detector
326 if( set->GetCDDetector(s,m,c) != det ) continue;
327
328 // Check we have the correct sector
329 if( set->GetCDSector(s,m,c) != sec ) continue;
330
331 // Check we have the right side
332 if( set->GetCDSide(s,m,c) != side ) continue;
333
334 // Check we have the correct strip
335 if( set->GetCDStrip(s,m,c) != strip ) continue;
336
337 // Then we got the right channel
338 sfp = s;
339 board = m;
340 ch = c;
341 return true;
342
343 } // c
344
345 } // m
346
347 } // s
348
349 std::cerr << "CD strip not found, det=" << det << ", sec=" << sec;
350 std::cerr << ", side=" << side << ", strip=" << strip << std::endl;
351 return false;
352
353}
354
356
357 // Check if we have old or new DAQ
358 bool oldDAQ = false;
359 if( set->GetNumberOfCaenAdcModules() > 0 )
360 oldDAQ = true;
361
362 // Create a TF1 for the linear fit
363 auto pfit = std::make_unique<TF1>( "pfit", "[0]+[1]*x", 0, 1e9 );
364
365 // Some canvases to check fits
366 gErrorIgnoreLevel = kError;
367 std::vector<std::vector<std::unique_ptr<TCanvas>>> canv;
368 canv.resize( set->GetNumberOfCDDetectors() );
369
370 // Loop over detectors
371 for( unsigned int i = 0; i < set->GetNumberOfCDDetectors(); ++i ) {
372
373 canv[i].resize( set->GetNumberOfCDSectors() );
374
375 // Loop over the sectors
376 for( unsigned int j = 0; j < set->GetNumberOfCDSectors(); ++j ) {
377
378 std::string cname = "cdcal_p_" + std::to_string(i) + "_" + std::to_string(j);
379 canv[i][j] = std::make_unique<TCanvas>( cname.data(), cname.data(), 800, 1000 );
380
381 // Loop over all the strips
382 for( unsigned int k = 0; k < set->GetNumberOfCDPStrips(); ++k ) {
383
384 // Get the right histogram to do the fit
385 auto res = cd_nQ_pQ[i][j][k]->Fit( pfit.get(), "QWL" );
386 if( res != 0 ) continue;
387 double fit_gain = ngain / pfit->GetParameter(1);
388 double fit_offset = noffset - pfit->GetParameter(0) * fit_gain;
389
390 // If we have the n-side tag, set the gain and offset
391 if( k == ptag ) {
392 std::cout << "!! This is the p-side tag channel, cross-check check the parameters below !!" << std::endl;
393 pgain = fit_gain;
394 poffset = fit_offset;
395 }
396
397 // Get the output names for the calibration file
398 std::string cal_base;
399 std::string modchstr;
400 int fsfp, fmod, fch;
401 if( oldDAQ ) {
402
403 // Search for the correct ADC and channel combination
404 cal_base = "adc_";
405 if( !FindCDChannels( i, j, 0, k, fmod, fch ) )
406 continue;
407 modchstr = std::to_string(fmod) + "_" + std::to_string(fch);
408
409 } // old DAQ
410
411 else {
412
413 // Search for the correct ADC and channel combination
414 cal_base = "febex_";
415 if( !FindCDChannels( i, j, 0, k, fsfp, fmod, fch ) )
416 continue;
417 modchstr = std::to_string(fsfp) + "_" + std::to_string(fmod);
418 modchstr += "_" + std::to_string(fch);
419
420 } // new DAQ
421
422 // Add gain and offset
423 std::string gainstr = cal_base + modchstr + ".Gain: " + std::to_string( fit_gain );
424 std::string offsetstr = cal_base + modchstr + ".Offset: " + std::to_string( fit_offset );
425
426 // Write them to the file
427 std::cout << gainstr << std::endl;
428 std::cout << offsetstr << std::endl;
429 output_cal << gainstr << std::endl;
430 output_cal << offsetstr << std::endl;
431
432 // Print to a file
433 std::string pdfname = cname + ".pdf";
434 if( k == 0 && set->GetNumberOfCDPStrips() != 1 )
435 pdfname += "(";
436 else if( k > 0 && k == set->GetNumberOfCDPStrips() - 1 )
437 pdfname += ")";
438 canv[i][j]->Print( pdfname.data(), "pdf" );
439
440 } // k
441
442 } // j
443
444 } // i
445
446 // Reset warning level
447 gErrorIgnoreLevel = kInfo;
448
449 return;
450
451}
452
454
455 // Check if we have old or new DAQ
456 bool oldDAQ = false;
457 if( set->GetNumberOfCaenAdcModules() > 0 )
458 oldDAQ = true;
459
460 // Create a TF1 for the linear fit
461 auto nfit = std::make_unique<TF1>( "nfit", "[0]+[1]*x", 0, 1e9 );
462
463 // Some canvases to check fits
464 gErrorIgnoreLevel = kError;
465 std::vector<std::vector<std::unique_ptr<TCanvas>>> canv;
466 canv.resize( set->GetNumberOfCDDetectors() );
467
468 // Loop over detectors
469 for( unsigned int i = 0; i < set->GetNumberOfCDDetectors(); ++i ) {
470
471 canv[i].resize( set->GetNumberOfCDSectors() );
472
473 // Loop over the sectors
474 for( unsigned int j = 0; j < set->GetNumberOfCDSectors(); ++j ) {
475
476 std::string cname = "cdcal_n_" + std::to_string(i) + "_" + std::to_string(j);
477 canv[i][j] = std::make_unique<TCanvas>( cname.data(), cname.data(), 800, 1000 );
478
479 // Loop over all the strips
480 for( unsigned int k = 0; k < set->GetNumberOfCDNStrips(); ++k ) {
481
482 // Get the right histogram to do the fit
483 auto res = cd_pen_nQ[i][j][k]->Fit( nfit.get(), "QWL" );
484 if( res != 0 ) continue;
485 double fit_gain = 1.0 / nfit->GetParameter(1);
486 double fit_offset = -1.0 * nfit->GetParameter(0) * fit_gain;
487 //double fit_gain = 1.0;
488 //double fit_offset = 0.0;
489
490 // If we have the n-side tag, set the gain and offset
491 if( k == ntag ) {
492 ngain = fit_gain;
493 noffset = fit_offset;
494 }
495
496 // Get the output names for the calibration file
497 std::string cal_base;
498 std::string modchstr;
499 int fsfp, fmod, fch;
500 if( oldDAQ ) {
501
502 // Search for the correct ADC and channel combination
503 cal_base = "adc_";
504 if( !FindCDChannels( i, j, 1, k, fmod, fch ) )
505 continue;
506 modchstr = std::to_string(fmod) + "_" + std::to_string(fch);
507
508 } // old DAQ
509
510 else {
511
512 // Search for the correct ADC and channel combination
513 cal_base = "febex_";
514 if( !FindCDChannels( i, j, 1, k, fsfp, fmod, fch ) )
515 continue;
516 modchstr = std::to_string(fsfp) + "_" + std::to_string(fmod);
517 modchstr += "_" + std::to_string(fch);
518
519 } // new DAQ
520
521 // Add gain and offset
522 std::string gainstr = cal_base + modchstr + ".Gain: " + std::to_string( fit_gain );
523 std::string offsetstr = cal_base + modchstr + ".Offset: " + std::to_string( fit_offset );
524
525 // Write them to the file
526 std::cout << gainstr << std::endl;
527 std::cout << offsetstr << std::endl;
528 output_cal << gainstr << std::endl;
529 output_cal << offsetstr << std::endl;
530
531 // Print to a file
532 std::string pdfname = cname + ".pdf";
533 if( k == 0 && set->GetNumberOfCDNStrips() != 1 )
534 pdfname += "(";
535 else if( k > 0 && k == set->GetNumberOfCDNStrips() - 1 )
536 pdfname += ")";
537 canv[i][j]->Print( pdfname.data(), "pdf" );
538
539 } // k
540
541 } // j
542
543 } // i
544
545 // Reset warning level
546 gErrorIgnoreLevel = kInfo;
547
548 return;
549
550}
551
553
554 // Variables for the finder algorithm
555 std::vector<unsigned char> pindex;
556 std::vector<unsigned char> nindex;
557
558 // Loop over each detector and sector
559 for( unsigned int i = 0; i < set->GetNumberOfCDDetectors(); ++i ){
560
561 for( unsigned int j = 0; j < set->GetNumberOfCDSectors(); ++j ){
562
563 // Reset variables for a new detector element
564 pindex.clear();
565 nindex.clear();
566 std::vector<unsigned char>().swap(pindex);
567 std::vector<unsigned char>().swap(nindex);
568
569 // Calculate p/n side multiplicities and get indicies
570 for( unsigned int k = 0; k < cd_en_list.size(); ++k ){
571
572 // Test that we have the correct detector and quadrant
573 if( i != cd_det_list.at(k) || j != cd_sec_list.at(k) )
574 continue;
575
576 // Check max energy and push back the multiplicity
577 if( cd_side_list.at(k) == 0 )
578 pindex.push_back(k);
579
580 else if( cd_side_list.at(k) == 1 )
581 nindex.push_back(k);
582
583 } // k: all CD events
584
585 // Keep only multiplicity 1v1 events
586 if( pindex.size() != 1 || nindex.size() != 1 )
587 continue;
588
589 // Fill the hit in the right pixel
590 int pid = cd_strip_list[pindex[0]];
591 int nid = cd_strip_list[nindex[0]];
592 double pen = cd_en_list[pindex[0]];
593 double nen = cd_en_list[nindex[0]];
594 unsigned int pQ = cd_Q_list[pindex[0]];
595 unsigned int nQ = cd_Q_list[nindex[0]];
596
597 // skip events with very diiferent energies
598 if( nQ / pQ > 1.5 || pQ / nQ > 1.5 ) continue;
599
600 // For p-side tags
601 if( pid == ptag ) {
602
603 cd_pen_nen[i][j][nid]->Fill( pen, nen );
604 cd_pen_nQ[i][j][nid]->Fill( pen, nQ );
605
606 }
607
608 // For n-side tags
609 if( nid == ntag ) {
610
611 cd_nen_pen[i][j][pid]->Fill( nen, pen );
612 cd_nQ_pQ[i][j][pid]->Fill( nQ, pQ );
613
614 }
615
616 } // j
617
618 } // i
619
620
621}
622
624
626
627 if( input_tree->LoadTree(0) < 0 ){
628
629 std::cout << " CD Calibrator: nothing to do" << std::endl;
630 return 0;
631
632 }
633
634 // Get ready and go
635 Initialise();
636 n_entries = input_tree->GetEntries();
637 n_mbs_entries = mbsinfo_tree->GetEntries();
638
639 std::cout << " CD Calibrator: number of entries in input tree = ";
640 std::cout << n_entries << std::endl;
641
642 std::cout << "\tnumber of MBS Events/triggers in input tree = ";
643 std::cout << n_mbs_entries << std::endl;
644
645 // ------------------------------------------------------------------------ //
646 // Main loop over TTree to find events
647 // ------------------------------------------------------------------------ //
648 for( unsigned long i = 0; i < n_entries; ++i ) {
649
650 // First event, yes please!
651 if( i == 0 ){
652
653 input_tree->GetEntry(i);
656
657 // Try to get the MBS info event with the index
658 if( mbsinfo_tree->GetEntryWithIndex( myeventid ) < 0 &&
659 n_mbs_entries > 0 ) {
660
661 // Look for the matches MBS Info event if we didn't match automatically
662 for( unsigned long j = 0; j < n_mbs_entries; ++j ){
663
664 mbsinfo_tree->GetEntry(j);
665 if( mbs_info->GetEventID() == myeventid ) {
667 break;
668 }
669
670 // Panic if we failed!
671 if( j+1 == n_mbs_entries ) {
672 std::cerr << "Didn't find matching MBS Event IDs at start of the file: ";
673 std::cerr << myeventid << std::endl;
674 }
675
676 }
677
678 }
679
680 //std::cout << "MBS Trigger time = " << myeventtime << std::endl;
681
682 }
683
684 // Get the time of the event
685 if( set->GetMbsEventMode() ) {
686
689
690 }
691
692 else mytime = in_data->GetTime();
693
694 // check time stamp monotonically increases in time-ordered mode!
695 if( time_prev > mytime && !set->GetMbsEventMode() ) {
696
697 std::cout << "Out of order event in ";
698 std::cout << input_tree->GetName() << std::endl;
699
700 }
701
702 // check event id is increasing in MBS event ordered mode
703 if( preveventid > myeventid && set->GetMbsEventMode() ) {
704
705 std::cout << "Out of order MBS event " << myeventid;
706 std::cout << " < " << preveventid << std::endl;
707
708 }
709
710 // record time of this event
712
713 // assume this is above threshold initially
714 mythres = true;
715
716 // ------------------------------------------ //
717 // Find FEBEX data
718 // ------------------------------------------ //
719 if( in_data->IsFebex() ) {
720
721 // Get the data
723 mysfp = febex_data->GetSfp();
724 myboard = febex_data->GetBoard();
725 mych = febex_data->GetChannel();
726 mypileup = febex_data->IsPileup();
727 myclipped = febex_data->IsClipped();
728
729 // Update calibration always for CD calibrator
730 unsigned int adc_tmp_value;
731 if( cal->FebexType( mysfp, myboard, mych ) == "Qshort" )
732 adc_tmp_value = febex_data->GetQshort();
733 else if( cal->FebexType( mysfp, myboard, mych ) == "Qint" )
734 adc_tmp_value = febex_data->GetQint();
735 else adc_tmp_value = febex_data->GetQshort();
736
737 myenergy = cal->FebexEnergy( mysfp, myboard, mych, adc_tmp_value );
738
739 if( adc_tmp_value > cal->FebexThreshold( mysfp, myboard, mych ) )
740 mythres = true;
741 else mythres = false;
742
743 // Is it a particle from the CD?
744 if( set->IsCD( mysfp, myboard, mych ) && mythres ) {
745
746 // Increment counts and open the event
747 hit_ctr++;
748
749 // Clipped rejection and pileup rejection
750 if( ( !myclipped || !set->GetClippedRejection() ) &&
751 ( !mypileup || !set->GetPileupRejection() ) ) {
752
753 event_open = true;
754 cd_en_list.push_back( myenergy );
755 cd_Q_list.push_back( adc_tmp_value );
756 cd_ts_list.push_back( mytime );
757 cd_det_list.push_back( set->GetCDDetector( mysfp, myboard, mych ) );
758 cd_sec_list.push_back( set->GetCDSector( mysfp, myboard, mych ) );
759 cd_side_list.push_back( set->GetCDSide( mysfp, myboard, mych ) );
760 cd_strip_list.push_back( set->GetCDStrip( mysfp, myboard, mych ) );
761
762 }
763
764 }
765
766 }
767
768 // ------------------------------------------ //
769 // Find ADC data
770 // ------------------------------------------ //
771 if( in_data->IsAdc() ) {
772
773 // Get the data
775 myadc = adc_data->GetModule();
776 mych = adc_data->GetChannel();
777 myclipped = adc_data->IsClipped();
778
779 // Update calibration always for CD calibrator
780 unsigned int adc_tmp_value = adc_data->GetQshort();
781 myenergy = cal->AdcEnergy( myadc, mych, adc_tmp_value );
782
783 if( adc_tmp_value > cal->AdcThreshold( myadc, mych ) )
784 mythres = true;
785 else mythres = false;
786
787 // Is it a particle from the CD?
788 if( set->IsCD( myadc, mych ) && mythres ) {
789
790 // Increment counts and open the event
791 hit_ctr++;
792
793 if( !myclipped || !set->GetClippedRejection() ) {
794
795 event_open = true;
796 cd_en_list.push_back( myenergy );
797 cd_Q_list.push_back( adc_tmp_value );
798 cd_ts_list.push_back( mytime );
799 cd_det_list.push_back( set->GetCDDetector( myadc, mych ) );
800 cd_sec_list.push_back( set->GetCDSector( myadc, mych ) );
801 cd_side_list.push_back( set->GetCDSide( myadc, mych ) );
802 cd_strip_list.push_back( set->GetCDStrip( myadc, mych ) );
803
804 }
805
806 }
807
808 }
809
810 // Sort out the timing for the event window
811 // but only if it isn't an info event, i.e only for real data
812 if( !in_data->IsInfo() ) {
813
814 // if this is first datum included in Event
815 if( hit_ctr == 1 && mythres ) {
816
820
821 }
822
823 // Update min and max
824 if( mytime > time_max ) time_max = mytime;
825 else if( mytime < time_min ) time_min = mytime;
826
827 } // not info data
828
829 //------------------------------
830 // check if last datum from this event and do some cleanup
831 //------------------------------
832
833 if( input_tree->GetEntry(i+1) ) {
834
835 // Get the next MBS event ID
838
839 // If the next MBS event ID is the same, carry on
840 // If not, we have to go look for the next trigger time
841 if( myeventid != preveventid ) {
842
843 // Close the event
844 flag_close_event = true;
845
846 // And find the next MBS event ID
847 if( mbsinfo_tree->GetEntryWithIndex( myeventid ) < 0 &&
848 n_mbs_entries > 0 ) {
849
850 std::cerr << "MBS Event " << myeventid << " not found by index, looking up manually" << std::endl;
851
852 // Look for the matches MBS Info event if we didn't match automatically
853 for( unsigned long j = 0; j < n_mbs_entries; ++j ){
854
855 mbsinfo_tree->GetEntry(j);
856 if( mbs_info->GetEventID() == myeventid ) {
858 break;
859 }
860
861 // Panic if we failed!
862 if( j+1 == n_mbs_entries ) {
863 std::cerr << "Didn't find matching MBS Event IDs at start of the file: ";
864 std::cerr << myeventid << std::endl;
865 }
866 }
867
868 }
869
870 else myeventtime = mbs_info->GetTime();
871
872 }
873
874 // BELOW IS THE TIME-ORDERED METHOD!
875
876 // Get the time of the next event
877 if( set->GetMbsEventMode() ) {
878
881
882 }
883
884 else mytime = in_data->GetTime();
885
886 // Calculate time diff
888
889 // window = time_stamp_first + time_window
890 if( time_diff > build_window )
891 flag_close_event = true; // set flag to close this event
892
893 // we've gone on to the next file in the chain
894 else if( time_diff < 0 )
895 flag_close_event = true; // set flag to close this event
896
897 } // if next entry beyond time window: close event!
898
899
900 //----------------------------
901 // if close this event or last entry
902 //----------------------------
903 if( flag_close_event || (i+1) == n_entries ) {
904
905 //--------------------------------------------------
906 // clear values of arrays to store intermediate info
907 //--------------------------------------------------
909 Initialise();
910
911 } // if close event && hit_ctr > 0
912
913 // Progress bar
914 bool update_progress = false;
915 if( n_entries < 200 )
916 update_progress = true;
917 else if( i % (n_entries/100) == 0 || i+1 == n_entries )
918 update_progress = true;
919
920 if( update_progress ) {
921
922 // Percent complete
923 float percent = (float)(i+1)*100.0/(float)n_entries;
924
925 // Progress bar in GUI
926 if( _prog_ ) {
927
928 prog->SetPosition( percent );
929 gSystem->ProcessEvents();
930
931 }
932
933 // Progress bar in terminal
934 std::cout << " " << std::setw(6) << std::setprecision(4);
935 std::cout << percent << "% \r";
936 std::cout.flush();
937
938 }
939
940 } // End of main loop over TTree to process raw FEBEX data entries (for n_entries)
941
942
943 //--------------------------
944 // Do the fitting to get calibration coefficients
945 //--------------------------
946
947 std::cout << "\n\nUsing p-side strip " << (int)ptag << " as reference for calibrating n-sides" << std::endl;
949 std::cout << "\n\nUsing n-side strip " << (int)ntag << " as reference for calibrating p-sides" << std::endl;
951
952 //--------------------------
953 // Clean up
954 //--------------------------
955
956 std::cout << "\n MiniballCDCalibrator finished..." << std::endl;
957
958 return n_entries;
959
960}
long long int GetTime()
unsigned long long int GetEventID()
unsigned long long myeventid
MBS event id.
unsigned long long mytime
absolute timestamp
std::shared_ptr< TGProgressBar > prog
unsigned long long time_max
void SetOutput(std::string output_file_name, bool cWrite=false)
unsigned char myboard
febex board number
bool mythres
above threshold?
unsigned long hit_ctr
unsigned char mysfp
sfp number
float myenergy
calibrated energy
unsigned long long preveventid
previous MBS event id
std::vector< unsigned int > cd_Q_list
list of CD uncalibrated energies for ParticleFinder
long myhittime
hit time with respect to event time
unsigned long long time_first
std::vector< std::vector< std::vector< TH2F * > > > cd_pen_nQ
std::ofstream output_cal
MBSInfoPackets * mbs_info
unsigned char myadc
ADC module number.
std::shared_ptr< FebexData > febex_data
bool FindCDChannels(int det, int sec, int side, int strip, int &adc, int &ch)
std::vector< unsigned char > cd_side_list
list of CD side IDs; 0 = p, 1 = n
unsigned long long time_prev
std::vector< std::vector< std::vector< TH2F * > > > cd_pen_nen
bool mypileup
pileup flag?
MiniballCDCalibrator(std::shared_ptr< MiniballSettings > myset)
unsigned long n_entries
unsigned long FillHists()
unsigned char mych
channel number
void Initialise()
called for every event
std::vector< std::vector< std::vector< TH2F * > > > cd_nQ_pQ
std::shared_ptr< MiniballCalibration > cal
std::vector< unsigned char > cd_det_list
list of CD detector IDs
void SetInputFile(std::vector< std::string > input_file_names)
unsigned long n_mbs_entries
std::vector< float > cd_en_list
list of CD energies for ParticleFinder
std::vector< unsigned char > cd_strip_list
list of CD strip IDs
std::vector< unsigned char > cd_sec_list
list of CD sector IDs
void SetInputTree(TTree *user_tree)
bool myclipped
clipped flag?
void StartFile()
called for every file
unsigned long long myeventtime
MBS event time.
std::vector< unsigned long long > cd_ts_list
list of CD timestamps for ParticleFinder
std::shared_ptr< MiniballSettings > set
MiniballDataPackets * in_data
void SetMBSInfoTree(TTree *user_tree)
std::shared_ptr< AdcData > adc_data
std::vector< std::vector< std::vector< TH2F * > > > cd_nen_pen
unsigned long long time_min
TFile * output_file
Outputs.
TFile * input_file
Input tree.
bool flag_close_event
length of build window in ns
unsigned long long int GetEventID() const
bool IsInfo() const
std::shared_ptr< FebexData > GetFebexData() const
bool IsFebex() const
std::shared_ptr< AdcData > GetAdcData() const
long long int GetTime() const
std::shared_ptr< MiniballSettings > myset
Definition mb_sort.cc:123