ATLAS Offline Software
Loading...
Searching...
No Matches
TRTCalib_bhadd.cxx File Reference
#include <fstream>
#include <sstream>
#include <iomanip>
#include <cstdlib>
#include <cstdio>
#include <cmath>
#include <TMath.h>
#include <iostream>
#include <sys/time.h>
#include <sys/resource.h>
#include <sys/sysinfo.h>
#include <unistd.h>
#include <map>
#include <vector>
#include <TFile.h>
#include <TNtuple.h>
#include <TH2F.h>
#include <TH1D.h>
#include <TVectorD.h>

Go to the source code of this file.

Classes

class  CalHist
class  CompBHist
class  trackdata
struct  calibpars

Functions

int Simple1dHist (float min, float max, int nbins, float value)
int Simple2dHist (float minx, float maxx, int nbinsx, float miny, float maxy, int nbinsy, float valuex, float valuey)
void process_mem_usage (double &vm_usage, double &resident_set)
map< string, calibparsreadoldpars ()
int dump_tracktuple (map< string, trackdata > *trackmap)
int dump_hists (map< int, CalHist * > *histmap, int ntbins, int nrbins, int, char *fname, int fileno)
int main (int argc, char *argv[])

Function Documentation

◆ dump_hists()

int dump_hists ( map< int, CalHist * > * histmap,
int ntbins,
int nrbins,
int ,
char * fname,
int fileno )

Definition at line 294 of file TRTCalib_bhadd.cxx.

295{
296
297 std::vector<int> temparray(ntbins * nrbins + 200, 0);
298 int maxnpop = 0;
299 int totnpop = 0;
300 ofstream ofile(Form("%s.part%i", fname, fileno), ios::out | ios::binary | ios::app);
301 for (map<int, CalHist *>::iterator it = histmap->begin(); it != histmap->end(); ++it)
302 {
303 it->second->GetArray(temparray, ntbins * nrbins + 200);
304 CompBHist *cbhist = new CompBHist(it->first, temparray, ntbins, nrbins);
305 cbhist->Write(&ofile);
306 if (cbhist->npop > maxnpop)
307 maxnpop = cbhist->npop;
308 totnpop += cbhist->npop;
309 delete cbhist;
310 }
311
312 cout << "DUMPING " << histmap->size() << " HISTOGRAMS, MAX NONZERO BINS: " << maxnpop << ", AVERAGE NONZERO BINS: " << (float)totnpop / histmap->size() << endl;
313
314 ofile.close();
315
316 return 0;
317}
void Write(ofstream *)
STL iterator class.

◆ dump_tracktuple()

int dump_tracktuple ( map< string, trackdata > * trackmap)

Definition at line 279 of file TRTCalib_bhadd.cxx.

280{
281
282 TFile *ttfile = new TFile("tracktuple.root", "UPDATE");
283 TNtuple *tracktup = new TNtuple("tracktuple", "track data", "run:evt:track:epnew:epold:nhits:t0:t:ttrack:theta:phi:d0:pt:trackres:trackresMean");
284 for (std::map<std::string, trackdata>::iterator iep = trackmap->begin(); iep != trackmap->end(); ++iep)
285 {
286 tracktup->Fill(iep->second.run, iep->second.event, iep->second.track, iep->second.epnew / iep->second.nhits, iep->second.epold, iep->second.nhits, iep->second.t0 / iep->second.nhits, iep->second.t / iep->second.nhits, iep->second.ttrack / iep->second.nhits, iep->second.theta, iep->second.phi, iep->second.d0, iep->second.pt, (iep->second.trackres / iep->second.nhits), (iep->second.trackresMean / iep->second.nhits));
287 }
288 ttfile->Write();
289 ttfile->Close();
290 // tracktup->Delete();
291 return 0;
292}

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 321 of file TRTCalib_bhadd.cxx.

322{
323
324 bool existAr = false;
325
326 float ptcutforerrors = 2;
327 // int who = RUSAGE_SELF;
328 // struct rusage usage;
329 // struct sysinfo sinfo;
330 // int ret;
331 double vm, rss;
332
333 float run, evt, trk, det, lay, mod, stl, stw, brd, chp, fsid, r, dr, t, rtrack, ttrack, drtrack, t0, ephase, theta, phi, pt, d0, ToT, HT, qoverp, isArgonStraw;
334 int sid;
335 pt = 0;
336 d0 = 0;
337 phi = 0;
338 theta = 0;
339 map<int, CalHist *> histmap;
340 map<string, trackdata> trackmap;
341
342 int ntbins = 55, nrbins = 100, ntresbins = 100, nresbins = 100, nhistbins, maxvalue = 0;
343 // int ntbins=64,nrbins=64,ntresbins=100,nresbins=100,nhistbins,maxvalue=0;
344 double minr = 0, maxr = 2, mint = -5, maxt = 50, mintres = -10, maxtres = 10, minres = -0.6, maxres = 0.6;
345 // double minr=0,maxr=2,mint=-10,maxt=80,mintres=-25,maxtres=25,minres=-1.0,maxres=1.0;
346
347 // int ntbins=64,nrbins=64,ntresbins=100,nresbins=100,nhistbins,maxvalue=0;
348 // double minr=0,maxr=2,mint=-10,maxt=80,mintres=-25,maxtres=25,minres=-1.0,maxres=1.0;
349
350 // sscanf (argv[1],"%i",&ntbins);
351 // sscanf (argv[2],"%i",&nrbins);
352 nhistbins = ntbins * nrbins + 200;
353
354 int npop;
355 // char header[50];
356 std::vector<int> histdata(nhistbins, 0);
357
358 int nhits = 0, nhists = 0;
359 int ntres = 0, nres = 0, nrt = 0;
360 int nhistsadd = 0, nbinsadd = 0, nhiststmp = 0;
361
362 int nfiles = argc - 2;
363 cout << "PROCESSING " << nfiles << " FILE(S)" << endl
364 << endl;
365
366 ofstream ofilestat((string(argv[1]) + ".stat").data(), ios::out);
367
368 TH1F *residual_ref = new TH1F("residual_-2_0_0_1", "residual_-2_0_0_1", nresbins, minres, maxres);
369 TH1F *timeresidual_ref = new TH1F("timeresidual_1_0_0_0_1", "timeresidual_1_0_0_0_1", ntresbins, mintres, maxtres);
370 TH2F *rt_ref = new TH2F("rt_-1_2", "rt_-1_2", ntbins, mint, maxt, nrbins, minr, maxr);
371 TH2F *rt_ref2 = new TH2F("rt_-1", "rt_-1", ntbins, mint, maxt, nrbins, minr, maxr);
372
373 TH2F *reshist_trt = new TH2F("residual_trt", "residual_trt", ntbins, mint, maxt, 100, -0.4, 0.4);
374 TH2F *reshist_bar = new TH2F("residual_bar", "residual_bar", ntbins, mint, maxt, 100, -0.4, 0.4);
375 TH2F *reshist1 = new TH2F("residual_1", "residual_1", ntbins, mint, maxt, 100, -0.4, 0.4);
376 TH2F *reshist2 = new TH2F("residual_-1", "residual_-1", ntbins, mint, maxt, 100, -0.4, 0.4);
377 TH2F *reshist3 = new TH2F("residual_2", "residual_2", ntbins, mint, maxt, 100, -0.4, 0.4);
378 TH2F *reshist4 = new TH2F("residual_-2", "residual_-2", ntbins, mint, maxt, 100, -0.4, 0.4);
379
380 TH2F *treshist_trt = new TH2F("timeresidual_trt", "timeresidual_trt", 50, mintres, maxtres, nrbins, minr, maxr);
381 TH2F *treshist_bar = new TH2F("timeresidual_bar", "timeresidual_bar", 50, mintres, maxtres, nrbins, minr, maxr);
382 TH2F *treshist1 = new TH2F("timeresidual_1", "timeresidual_1", 50, mintres, maxtres, nrbins, minr, maxr);
383 TH2F *treshist2 = new TH2F("timeresidual_-1", "timeresidual_-1", 50, mintres, maxtres, nrbins, minr, maxr);
384 TH2F *treshist3 = new TH2F("timeresidual_2", "timeresidual_2", 50, mintres, maxtres, nrbins, minr, maxr);
385 TH2F *treshist4 = new TH2F("timeresidual_-2", "timeresidual_-2", 50, mintres, maxtres, nrbins, minr, maxr);
386
387 map<string, calibpars> pt0map = readoldpars();
388
389 // Histograms to do error study: Barrel/Ecs:
390
391 // Pt dependent:
392
393 // Pull All hits:
394 TH2F *pull_trt = new TH2F("pull_trt", "pull_TRT_allhits", 16, 0, 16, 200, -4, 4);
395 TH2F *pull_ba = new TH2F("pull_ba", "pull_BarrelA_allhits", 16, 0, 16, 200, -4, 4);
396 TH2F *pull_bc = new TH2F("pull_bc", "pull_BarrelC_allhits", 16, 0, 16, 200, -4, 4);
397 TH2F *pull_ea = new TH2F("pull_ea", "pull_EndcapA_allhits", 16, 0, 16, 200, -4, 4);
398 TH2F *pull_ec = new TH2F("pull_ec", "pull_EndcapC_allhits", 16, 0, 16, 200, -4, 4);
399
400 // Residual All hits:
401 TH2F *residual_trt = new TH2F("residual_allhits", "residual_TRT_allhits", 16, 0, 16, 200, -2, 2);
402 TH2F *residual_ba = new TH2F("residual_ba", "residual_BarrelA_allhits", 16, 0, 16, 200, -2, 2);
403 TH2F *residual_bc = new TH2F("residual_bc", "residual_BarrelC_allhits", 16, 0, 16, 200, -2, 2);
404 TH2F *residual_ea = new TH2F("residual_ea", "residual_EndcapA_allhits", 16, 0, 16, 200, -2, 2);
405 TH2F *residual_ec = new TH2F("residual_ec", "residual_EndcapC_allhits", 16, 0, 16, 200, -2, 2);
406
407 // Time Residual All hits:
408 TH2F *tresidual_trt = new TH2F("tresidual_trt", "tresidual_TRT_allhits", 16, 0, 16, 200, -20, 20);
409 TH2F *tresidual_ba = new TH2F("tresidual_ba", "tresidual_BarrelA_allhits", 16, 0, 16, 200, -20, 20);
410 TH2F *tresidual_bc = new TH2F("tresidual_bc", "tresidual_BarrelC_allhits", 16, 0, 16, 200, -20, 20);
411 TH2F *tresidual_ea = new TH2F("tresidual_ea", "tresidual_EndcapA_allhits", 16, 0, 16, 200, -20, 20);
412 TH2F *tresidual_ec = new TH2F("tresidual_ec", "tresidual_EndcapC_allhits", 16, 0, 16, 200, -20, 20);
413
414 // Pull Precision hits:
415 TH2F *pull_trtP = new TH2F("pull_trtP", "pull_TRT_Phits", 16, 0, 16, 200, -4, 4);
416 TH2F *pull_baP = new TH2F("pull_baP", "pull_BarrelA_Phits", 16, 0, 16, 200, -4, 4);
417 TH2F *pull_bcP = new TH2F("pull_bcP", "pull_BarrelC_Phits", 16, 0, 16, 200, -4, 4);
418 TH2F *pull_eaP = new TH2F("pull_eaP", "pull_EndcapA_Phits", 16, 0, 16, 200, -4, 4);
419 TH2F *pull_ecP = new TH2F("pull_ecP", "pull_EndcapC_Phits", 16, 0, 16, 200, -4, 4);
420 // Residaul Precision hits:
421 TH2F *residual_trtP = new TH2F("residual_trtP", "residual_TRT_Phits", 16, 0, 16, 200, -2, 2);
422 TH2F *residual_baP = new TH2F("residual_baP", "residual_BarrelA_Phits", 16, 0, 16, 200, -2, 2);
423 TH2F *residual_bcP = new TH2F("residual_bcP", "residual_BarrelC_Phits", 16, 0, 16, 200, -2, 2);
424 TH2F *residual_eaP = new TH2F("residual_eaP", "residual_EndcapA_Phits", 16, 0, 16, 200, -2, 2);
425 TH2F *residual_ecP = new TH2F("residual_ecP", "residual_EndcapC_Phits", 16, 0, 16, 200, -2, 2);
426 // Residaul Precision hits:
427 TH2F *tresidual_trtP = new TH2F("tresidual_trtP", "residual_TRT_Phits", 16, 0, 16, 200, -20, 20);
428 TH2F *tresidual_baP = new TH2F("tresidual_baP", "residual_BarrelA_Phits", 16, 0, 16, 200, -20, 20);
429 TH2F *tresidual_bcP = new TH2F("tresidual_bcP", "residual_BarrelC_Phits", 16, 0, 16, 200, -20, 20);
430 TH2F *tresidual_eaP = new TH2F("tresidual_eaP", "residual_EndcapA_Phits", 16, 0, 16, 200, -20, 20);
431 TH2F *tresidual_ecP = new TH2F("tresidual_ecP", "residual_EndcapC_Phits", 16, 0, 16, 200, -20, 20);
432 // Track Errors :
433 TH2F *trackerrors_trtP = new TH2F("trackerrors_trt", "trackerrors_vs_pt_trt", 16, 0, 16, 200, 0, 0.2);
434 TH2F *trackerrors_baP = new TH2F("trackerrors_ba", "trackerrors_vs_pt_BarrelA", 16, 0, 16, 200, 0, .2);
435 TH2F *trackerrors_bcP = new TH2F("trackerrors_bc", "trackerrors_vs_pt_BarrelC", 16, 0, 16, 200, 0, .2);
436 TH2F *trackerrors_eaP = new TH2F("trackerrors_ea", "trackerrors_vs_pt_EndcapA", 16, 0, 16, 200, 0, .2);
437 TH2F *trackerrors_ecP = new TH2F("trackerrors_ec", "trackerrors_vs_pt_EndcapC", 16, 0, 16, 200, 0, .2);
438
439 // for (map<string,calibpars>::iterator iot=pt0map.begin(); iot!=pt0map.end(); iot++) {
440
441 // Time Bin Dependent:
442
443 int ntbins2 = 30;
444 float tmin = -10;
445 float tmax = 75;
446
447 // TH2 to Store the Pull vs width:
448 TH2F *pull_vs_tb_trt = new TH2F("pull_vs_tb_trt", "pull_vs_timebin_trt", ntbins2, tmin, tmax, 200, -4, 4);
449 TH2F *pull_vs_tb_ba = new TH2F("pull_vs_tb_ba", "pull_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, -4, 4);
450 TH2F *pull_vs_tb_bc = new TH2F("pull_vs_tb_bc", "pull_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, -4, 4);
451 TH2F *pull_vs_tb_ea = new TH2F("pull_vs_tb_ea", "pull_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, -4, 4);
452 TH2F *pull_vs_tb_ec = new TH2F("pull_vs_tb_ec", "pull_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, -4, 4);
453
454 // TH2 to Store the ErrorsUsed vs width:
455 TH2F *errors_vs_tb_trt = new TH2F("errors_vs_tb_trt", "errors_vs_timebin_trt", ntbins2, tmin, tmax, 200, 0, 0.4);
456 TH2F *errors_vs_tb_ba = new TH2F("errors_vs_tb_ba", "errors_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, 0, 0.4);
457 TH2F *errors_vs_tb_bc = new TH2F("errors_vs_tb_bc", "errors_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, 0, 0.4);
458 TH2F *errors_vs_tb_ea = new TH2F("errors_vs_tb_ea", "errors_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, 0, 0.4);
459 TH2F *errors_vs_tb_ec = new TH2F("errors_vs_tb_ec", "errors_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, 0, 0.4);
460
461 // TH2 to Store the ErrorsUsed vs width:
462 TH2F *residual_vs_tb_trt = new TH2F("residual_vs_tb_trt", "residual_vs_timebin_trt", ntbins2, tmin, tmax, 200, -2, .2);
463 TH2F *residual_vs_tb_ba = new TH2F("residual_vs_tb_ba", "residual_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, -2, .2);
464 TH2F *residual_vs_tb_bc = new TH2F("residual_vs_tb_bc", "residual_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, -2, .2);
465 TH2F *residual_vs_tb_ea = new TH2F("residual_vs_tb_ea", "residual_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, -2, .2);
466 TH2F *residual_vs_tb_ec = new TH2F("residual_vs_tb_ec", "residual_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, -2, .2);
467
468 // TH2 to Store the Track Errors vs width:
469 TH2F *trackerrors_vs_tb_trt = new TH2F("trackerrors_vs_tb_trt", "trackerrors_vs_timebin_trt", ntbins2, tmin, tmax, 200, 0, 0.2);
470 TH2F *trackerrors_vs_tb_ba = new TH2F("trackerrors_vs_tb_ba", "trackerrors_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, 0, .2);
471 TH2F *trackerrors_vs_tb_bc = new TH2F("trackerrors_vs_tb_bc", "trackerrors_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, 0, .2);
472 TH2F *trackerrors_vs_tb_ea = new TH2F("trackerrors_vs_tb_ea", "trackerrors_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, 0, .2);
473 TH2F *trackerrors_vs_tb_ec = new TH2F("trackerrors_vs_tb_ec", "trackerrors_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, 0, .2);
474
475 // Histograms to store the ToT Corrections, TIME RESIDUAL:
476 TH2F *tres_vs_ToT_trt = new TH2F("tres_vs_ToT_trt", "tres_vs_ToT_trt", 24, -0.5, 23.5, 200, -20, 20);
477 TH2F *tres_vs_ToT_ba = new TH2F("tres_vs_ToT_ba", "tres_vs_ToT_ba", 24, -0.5, 23.5, 200, -20, 20);
478 TH2F *tres_vs_ToT_bc = new TH2F("tres_vs_ToT_bc", "tres_vs_ToT_bc", 24, -0.5, 23.5, 200, -20, 20);
479 TH2F *tres_vs_ToT_ea = new TH2F("tres_vs_ToT_ea", "tres_vs_ToT_ea", 24, -0.5, 23.5, 200, -20, 20);
480 TH2F *tres_vs_ToT_ec = new TH2F("tres_vs_ToT_ec", "tres_vs_ToT_ec", 24, -0.5, 23.5, 200, -20, 20);
481
482 // Residual:
483 TH2F *res_vs_ToT_trt = new TH2F("res_vs_ToT_trt", "res_vs_ToT_trt", 24, -0.5, 23.5, 200, -2, 2);
484 TH2F *res_vs_ToT_ba = new TH2F("res_vs_ToT_ba", "res_vs_ToT_ba", 24, -0.5, 23.5, 200, -2, 2);
485 TH2F *res_vs_ToT_bc = new TH2F("res_vs_ToT_bc", "res_vs_ToT_bc", 24, -0.5, 23.5, 200, -2, 2);
486 TH2F *res_vs_ToT_ea = new TH2F("res_vs_ToT_ea", "res_vs_ToT_ea", 24, -0.5, 23.5, 200, -2, 2);
487 TH2F *res_vs_ToT_ec = new TH2F("res_vs_ToT_ec", "res_vs_ToT_ec", 24, -0.5, 23.5, 200, -2, 2);
488
489 // Histograms to store the HT Corrections, TIME RESIDUAL:
490 TH2F *tres_vs_HT_trt = new TH2F("tres_vs_HT_trt", "tres_vs_HT_trt", 2, -0.5, 1.5, 200, -20, 20);
491 TH2F *tres_vs_HT_ba = new TH2F("tres_vs_HT_ba", "tres_vs_HT_ba", 2, -0.5, 1.5, 200, -20, 20);
492 TH2F *tres_vs_HT_bc = new TH2F("tres_vs_HT_bc", "tres_vs_HT_bc", 2, -0.5, 1.5, 200, -20, 20);
493 TH2F *tres_vs_HT_ea = new TH2F("tres_vs_HT_ea", "tres_vs_HT_ea", 2, -0.5, 1.5, 200, -20, 20);
494 TH2F *tres_vs_HT_ec = new TH2F("tres_vs_HT_ec", "tres_vs_HT_ec", 2, -0.5, 1.5, 200, -20, 20);
495
496 // Residual:
497 TH2F *res_vs_HT_trt = new TH2F("res_vs_HT_trt", "res_vs_HT_trt", 2, -0.5, 1.5, 200, -2, 2);
498 TH2F *res_vs_HT_ba = new TH2F("res_vs_HT_ba", "res_vs_HT_ba", 2, -0.5, 1.5, 200, -2, 2);
499 TH2F *res_vs_HT_bc = new TH2F("res_vs_HT_bc", "res_vs_HT_bc", 2, -0.5, 1.5, 200, -2, 2);
500 TH2F *res_vs_HT_ea = new TH2F("res_vs_HT_ea", "res_vs_HT_ea", 2, -0.5, 1.5, 200, -2, 2);
501 TH2F *res_vs_HT_ec = new TH2F("res_vs_HT_ec", "res_vs_HT_ec", 2, -0.5, 1.5, 200, -2, 2);
502
503 // Histograms to store the sin(theta)/p Corrections, TIME RESIDUAL:
504 TH2F *tres_vs_SinOverP_trt = new TH2F("tres_vs_SinOverP_trt", "tres_vs_SinOverP_trt", 24, 0., 1., 200, -20., 20.);
505 TH2F *tres_vs_SinOverP_ba = new TH2F("tres_vs_SinOverP_ba", "tres_vs_SinOverP_ba", 24, 0., 1., 200, -20., 20.);
506 TH2F *tres_vs_SinOverP_bc = new TH2F("tres_vs_SinOverP_bc", "tres_vs_SinOverP_bc", 24, 0., 1., 200, -20., 20.);
507 TH2F *tres_vs_SinOverP_ea = new TH2F("tres_vs_SinOverP_ea", "tres_vs_SinOverP_ea", 24, 0., 1., 200, -20., 20.);
508 TH2F *tres_vs_SinOverP_ec = new TH2F("tres_vs_SinOverP_ec", "tres_vs_SinOverP_ec", 24, 0., 1., 200, -20., 20.);
509
510 // Residual:
511 TH2F *res_vs_SinOverP_trt = new TH2F("res_vs_SinOverP_trt", "res_vs_SinOverP_trt", 24, 0., 1., 200, -2., 2.);
512 TH2F *res_vs_SinOverP_ba = new TH2F("res_vs_SinOverP_ba", "res_vs_SinOverP_ba", 24, 0., 1., 200, -2., 2.);
513 TH2F *res_vs_SinOverP_bc = new TH2F("res_vs_SinOverP_bc", "res_vs_SinOverP_bc", 24, 0., 1., 200, -2., 2.);
514 TH2F *res_vs_SinOverP_ea = new TH2F("res_vs_SinOverP_ea", "res_vs_SinOverP_ea", 24, 0., 1., 200, -2., 2.);
515 TH2F *res_vs_SinOverP_ec = new TH2F("res_vs_SinOverP_ec", "res_vs_SinOverP_ec", 24, 0., 1., 200, -2., 2.);
516
517 // Histograms to store the sin(theta)/p Corrections, TIME RESIDUAL:
518 TH2F *tres_vs_CosOverP_trt = new TH2F("tres_vs_CosOverP_trt", "tres_vs_CosOverP_trt", 24, 0., 1., 200, -20., 20.);
519 TH2F *tres_vs_CosOverP_ba = new TH2F("tres_vs_CosOverP_ba", "tres_vs_CosOverP_ba", 24, 0., 1., 200, -20., 20.);
520 TH2F *tres_vs_CosOverP_bc = new TH2F("tres_vs_CosOverP_bc", "tres_vs_CosOverP_bc", 24, 0., 1., 200, -20., 20.);
521 TH2F *tres_vs_CosOverP_ea = new TH2F("tres_vs_CosOverP_ea", "tres_vs_CosOverP_ea", 24, 0., 1., 200, -20., 20.);
522 TH2F *tres_vs_CosOverP_ec = new TH2F("tres_vs_CosOverP_ec", "tres_vs_CosOverP_ec", 24, 0., 1., 200, -20., 20.);
523
524 // Residual:
525 TH2F *res_vs_CosOverP_trt = new TH2F("res_vs_CosOverP_trt", "res_vs_CosOverP_trt", 24, 0., 1., 200, -2., 2.);
526 TH2F *res_vs_CosOverP_ba = new TH2F("res_vs_CosOverP_ba", "res_vs_CosOverP_ba", 24, 0., 1., 200, -2., 2.);
527 TH2F *res_vs_CosOverP_bc = new TH2F("res_vs_CosOverP_bc", "res_vs_CosOverP_bc", 24, 0., 1., 200, -2., 2.);
528 TH2F *res_vs_CosOverP_ea = new TH2F("res_vs_CosOverP_ea", "res_vs_CosOverP_ea", 24, 0., 1., 200, -2., 2.);
529 TH2F *res_vs_CosOverP_ec = new TH2F("res_vs_CosOverP_ec", "res_vs_CosOverP_ec", 24, 0., 1., 200, -2., 2.);
530
531 // HISTOGRAMS FOR ARGON STRAWS!!!!
532 TH2F *resArhist_trt = new TH2F("residualAr_trt", "residualAr_trt", ntbins, mint, maxt, 100, -0.4, 0.4);
533 TH2F *resArhist_bar = new TH2F("residualAr_bar", "residualAr_bar", ntbins, mint, maxt, 100, -0.4, 0.4);
534 TH2F *resArhist1 = new TH2F("residualAr_1", "residualAr_1", ntbins, mint, maxt, 100, -0.4, 0.4);
535 TH2F *resArhist2 = new TH2F("residualAr_-1", "residualAr_-1", ntbins, mint, maxt, 100, -0.4, 0.4);
536 TH2F *resArhist3 = new TH2F("residualAr_2", "residualAr_2", ntbins, mint, maxt, 100, -0.4, 0.4);
537 TH2F *resArhist4 = new TH2F("residualAr_-2", "residualAr_-2", ntbins, mint, maxt, 100, -0.4, 0.4);
538
539 TH2F *tresArhist_trt = new TH2F("timeresidualAr_trt", "timeresidualAr_trt", 50, mintres, maxtres, nrbins, minr, maxr);
540 TH2F *tresArhist_bar = new TH2F("timeresidualAr_bar", "timeresidualAr_bar", 50, mintres, maxtres, nrbins, minr, maxr);
541 TH2F *tresArhist1 = new TH2F("timeresidualAr_1", "timeresidualAr_1", 50, mintres, maxtres, nrbins, minr, maxr);
542 TH2F *tresArhist2 = new TH2F("timeresidualAr_-1", "timeresidualAr_-1", 50, mintres, maxtres, nrbins, minr, maxr);
543 TH2F *tresArhist3 = new TH2F("timeresidualAr_2", "timeresidualAr_2", 50, mintres, maxtres, nrbins, minr, maxr);
544 TH2F *tresArhist4 = new TH2F("timeresidualAr_-2", "timeresidualAr_-2", 50, mintres, maxtres, nrbins, minr, maxr);
545 // Histograms to do error study: Barrel/Ecs:
546
547 // Pt dependent:
548 // Pull All hits:
549 TH2F *pullAr_trt = new TH2F("pullAr_trt", "pullAr_TRT_allhits", 16, 0, 16, 200, -4, 4);
550 TH2F *pullAr_ba = new TH2F("pullAr_ba", "pullAr_BarrelA_allhits", 16, 0, 16, 200, -4, 4);
551 TH2F *pullAr_bc = new TH2F("pullAr_bc", "pullAr_BarrelC_allhits", 16, 0, 16, 200, -4, 4);
552 TH2F *pullAr_ea = new TH2F("pullAr_ea", "pullAr_EndcapA_allhits", 16, 0, 16, 200, -4, 4);
553 TH2F *pullAr_ec = new TH2F("pullAr_ec", "pullAr_EndcapC_allhits", 16, 0, 16, 200, -4, 4);
554
555 // Residual All hits:
556 TH2F *residualAr_trt = new TH2F("residualAr_allhits", "residualAr_TRT_allhits", 16, 0, 16, 200, -2, 2);
557 TH2F *residualAr_ba = new TH2F("residualAr_ba", "residualAr_BarrelA_allhits", 16, 0, 16, 200, -2, 2);
558 TH2F *residualAr_bc = new TH2F("residualAr_bc", "residualAr_BarrelC_allhits", 16, 0, 16, 200, -2, 2);
559 TH2F *residualAr_ea = new TH2F("residualAr_ea", "residualAr_EndcapA_allhits", 16, 0, 16, 200, -2, 2);
560 TH2F *residualAr_ec = new TH2F("residualAr_ec", "residualAr_EndcapC_allhits", 16, 0, 16, 200, -2, 2);
561
562 // Time Residual All hits:
563 TH2F *tresidualAr_trt = new TH2F("tresidualAr_trt", "tresidualAr_TRT_allhits", 16, 0, 16, 200, -20, 20);
564 TH2F *tresidualAr_ba = new TH2F("tresidualAr_ba", "tresidualAr_BarrelA_allhits", 16, 0, 16, 200, -20, 20);
565 TH2F *tresidualAr_bc = new TH2F("tresidualAr_bc", "tresidualAr_BarrelC_allhits", 16, 0, 16, 200, -20, 20);
566 TH2F *tresidualAr_ea = new TH2F("tresidualAr_ea", "tresidualAr_EndcapA_allhits", 16, 0, 16, 200, -20, 20);
567 TH2F *tresidualAr_ec = new TH2F("tresidualAr_ec", "tresidualAr_EndcapC_allhits", 16, 0, 16, 200, -20, 20);
568
569 // Pull Precision hits:
570 TH2F *pullAr_trtP = new TH2F("pullAr_trtP", "pullAr_TRT_Phits", 16, 0, 16, 200, -4, 4);
571 TH2F *pullAr_baP = new TH2F("pullAr_baP", "pullAr_BarrelA_Phits", 16, 0, 16, 200, -4, 4);
572 TH2F *pullAr_bcP = new TH2F("pullAr_bcP", "pullAr_BarrelC_Phits", 16, 0, 16, 200, -4, 4);
573 TH2F *pullAr_eaP = new TH2F("pullAr_eaP", "pullAr_EndcapA_Phits", 16, 0, 16, 200, -4, 4);
574 TH2F *pullAr_ecP = new TH2F("pullAr_ecP", "pullAr_EndcapC_Phits", 16, 0, 16, 200, -4, 4);
575 // Residaul Precision hits:
576 TH2F *residualAr_trtP = new TH2F("residualAr_trtP", "residualAr_TRT_Phits", 16, 0, 16, 200, -2, 2);
577 TH2F *residualAr_baP = new TH2F("residualAr_baP", "residualAr_BarrelA_Phits", 16, 0, 16, 200, -2, 2);
578 TH2F *residualAr_bcP = new TH2F("residualAr_bcP", "residualAr_BarrelC_Phits", 16, 0, 16, 200, -2, 2);
579 TH2F *residualAr_eaP = new TH2F("residualAr_eaP", "residualAr_EndcapA_Phits", 16, 0, 16, 200, -2, 2);
580 TH2F *residualAr_ecP = new TH2F("residualAr_ecP", "residualAr_EndcapC_Phits", 16, 0, 16, 200, -2, 2);
581 // Residaul Precision hits:
582 TH2F *tresidualAr_trtP = new TH2F("tresidualAr_trtP", "residualAr_TRT_Phits", 16, 0, 16, 200, -20, 20);
583 TH2F *tresidualAr_baP = new TH2F("tresidualAr_baP", "residualAr_BarrelA_Phits", 16, 0, 16, 200, -20, 20);
584 TH2F *tresidualAr_bcP = new TH2F("tresidualAr_bcP", "residualAr_BarrelC_Phits", 16, 0, 16, 200, -20, 20);
585 TH2F *tresidualAr_eaP = new TH2F("tresidualAr_eaP", "residualAr_EndcapA_Phits", 16, 0, 16, 200, -20, 20);
586 TH2F *tresidualAr_ecP = new TH2F("tresidualAr_ecP", "residualAr_EndcapC_Phits", 16, 0, 16, 200, -20, 20);
587 // Track Errors :
588 TH2F *trackerrorsAr_trtP = new TH2F("trackerrorsAr_trt", "trackerrorsAr_vs_pt_trt", 16, 0, 16, 200, 0, 0.2);
589 TH2F *trackerrorsAr_baP = new TH2F("trackerrorsAr_ba", "trackerrorsAr_vs_pt_BarrelA", 16, 0, 16, 200, 0, .2);
590 TH2F *trackerrorsAr_bcP = new TH2F("trackerrorsAr_bc", "trackerrorsAr_vs_pt_BarrelC", 16, 0, 16, 200, 0, .2);
591 TH2F *trackerrorsAr_eaP = new TH2F("trackerrorsAr_ea", "trackerrorsAr_vs_pt_EndcapA", 16, 0, 16, 200, 0, .2);
592 TH2F *trackerrorsAr_ecP = new TH2F("trackerrorsAr_ec", "trackerrorsAr_vs_pt_EndcapC", 16, 0, 16, 200, 0, .2);
593
594 // Time Bin Dependent:
595
596 // TH2 to Store the Pull vs width:
597 TH2F *pullAr_vs_tb_trt = new TH2F("pullAr_vs_tb_trt", "pullAr_vs_timebin_trt", ntbins2, tmin, tmax, 200, -4, 4);
598 TH2F *pullAr_vs_tb_ba = new TH2F("pullAr_vs_tb_ba", "pullAr_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, -4, 4);
599 TH2F *pullAr_vs_tb_bc = new TH2F("pullAr_vs_tb_bc", "pullAr_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, -4, 4);
600 TH2F *pullAr_vs_tb_ea = new TH2F("pullAr_vs_tb_ea", "pullAr_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, -4, 4);
601 TH2F *pullAr_vs_tb_ec = new TH2F("pullAr_vs_tb_ec", "pullAr_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, -4, 4);
602
603 // TH2 to Store the ErrorsUsed vs width:
604 TH2F *errorsAr_vs_tb_trt = new TH2F("errorsAr_vs_tb_trt", "errorsAr_vs_timebin_trt", ntbins2, tmin, tmax, 200, 0, 0.4);
605 TH2F *errorsAr_vs_tb_ba = new TH2F("errorsAr_vs_tb_ba", "errorsAr_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, 0, 0.4);
606 TH2F *errorsAr_vs_tb_bc = new TH2F("errorsAr_vs_tb_bc", "errorsAr_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, 0, 0.4);
607 TH2F *errorsAr_vs_tb_ea = new TH2F("errorsAr_vs_tb_ea", "errorsAr_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, 0, 0.4);
608 TH2F *errorsAr_vs_tb_ec = new TH2F("errorsAr_vs_tb_ec", "errorsAr_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, 0, 0.4);
609
610 // TH2 to Store the ErrorsUsed vs width:
611 TH2F *residualAr_vs_tb_trt = new TH2F("residualAr_vs_tb_trt", "residualAr_vs_timebin_trt", ntbins2, tmin, tmax, 200, -2, .2);
612 TH2F *residualAr_vs_tb_ba = new TH2F("residualAr_vs_tb_ba", "residualAr_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, -2, .2);
613 TH2F *residualAr_vs_tb_bc = new TH2F("residualAr_vs_tb_bc", "residualAr_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, -2, .2);
614 TH2F *residualAr_vs_tb_ea = new TH2F("residualAr_vs_tb_ea", "residualAr_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, -2, .2);
615 TH2F *residualAr_vs_tb_ec = new TH2F("residualAr_vs_tb_ec", "residualAr_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, -2, .2);
616
617 // TH2 to Store the Track Errors vs width:
618 TH2F *trackerrorsAr_vs_tb_trt = new TH2F("trackerrorsAr_vs_tb_trt", "trackerrorsAr_vs_timebin_trt", ntbins2, tmin, tmax, 200, 0, 0.2);
619 TH2F *trackerrorsAr_vs_tb_ba = new TH2F("trackerrorsAr_vs_tb_ba", "trackerrorsAr_vs_timebin_BarrelA", ntbins2, tmin, tmax, 200, 0, .2);
620 TH2F *trackerrorsAr_vs_tb_bc = new TH2F("trackerrorsAr_vs_tb_bc", "trackerrorsAr_vs_timebin_BarrelC", ntbins2, tmin, tmax, 200, 0, .2);
621 TH2F *trackerrorsAr_vs_tb_ea = new TH2F("trackerrorsAr_vs_tb_ea", "trackerrorsAr_vs_timebin_EndcapA", ntbins2, tmin, tmax, 200, 0, .2);
622 TH2F *trackerrorsAr_vs_tb_ec = new TH2F("trackerrorsAr_vs_tb_ec", "trackerrorsAr_vs_timebin_EndcapC", ntbins2, tmin, tmax, 200, 0, .2);
623
624 // Histograms to store the ToT Corrections, TIME RESIDUAL:
625 TH2F *tresAr_vs_ToT_trt = new TH2F("tresAr_vs_ToT_trt", "tresAr_vs_ToT_trt", 24, -0.5, 23.5, 200, -20, 20);
626 TH2F *tresAr_vs_ToT_ba = new TH2F("tresAr_vs_ToT_ba", "tresAr_vs_ToT_ba", 24, -0.5, 23.5, 200, -20, 20);
627 TH2F *tresAr_vs_ToT_bc = new TH2F("tresAr_vs_ToT_bc", "tresAr_vs_ToT_bc", 24, -0.5, 23.5, 200, -20, 20);
628 TH2F *tresAr_vs_ToT_ea = new TH2F("tresAr_vs_ToT_ea", "tresAr_vs_ToT_ea", 24, -0.5, 23.5, 200, -20, 20);
629 TH2F *tresAr_vs_ToT_ec = new TH2F("tresAr_vs_ToT_ec", "tresAr_vs_ToT_ec", 24, -0.5, 23.5, 200, -20, 20);
630
631 // Residual:
632 TH2F *resAr_vs_ToT_trt = new TH2F("resAr_vs_ToT_trt", "resAr_vs_ToT_trt", 24, -0.5, 23.5, 200, -2, 2);
633 TH2F *resAr_vs_ToT_ba = new TH2F("resAr_vs_ToT_ba", "resAr_vs_ToT_ba", 24, -0.5, 23.5, 200, -2, 2);
634 TH2F *resAr_vs_ToT_bc = new TH2F("resAr_vs_ToT_bc", "resAr_vs_ToT_bc", 24, -0.5, 23.5, 200, -2, 2);
635 TH2F *resAr_vs_ToT_ea = new TH2F("resAr_vs_ToT_ea", "resAr_vs_ToT_ea", 24, -0.5, 23.5, 200, -2, 2);
636 TH2F *resAr_vs_ToT_ec = new TH2F("resAr_vs_ToT_ec", "resAr_vs_ToT_ec", 24, -0.5, 23.5, 200, -2, 2);
637
638 // Histograms to store the HT Corrections, TIME RESIDUAL:
639 TH2F *tresAr_vs_HT_trt = new TH2F("tresAr_vs_HT_trt", "tresAr_vs_HT_trt", 2, -0.5, 1.5, 200, -20, 20);
640 TH2F *tresAr_vs_HT_ba = new TH2F("tresAr_vs_HT_ba", "tresAr_vs_HT_ba", 2, -0.5, 1.5, 200, -20, 20);
641 TH2F *tresAr_vs_HT_bc = new TH2F("tresAr_vs_HT_bc", "tresAr_vs_HT_bc", 2, -0.5, 1.5, 200, -20, 20);
642 TH2F *tresAr_vs_HT_ea = new TH2F("tresAr_vs_HT_ea", "tresAr_vs_HT_ea", 2, -0.5, 1.5, 200, -20, 20);
643 TH2F *tresAr_vs_HT_ec = new TH2F("tresAr_vs_HT_ec", "tresAr_vs_HT_ec", 2, -0.5, 1.5, 200, -20, 20);
644
645 // Residual:
646 TH2F *resAr_vs_HT_trt = new TH2F("resAr_vs_HT_trt", "resAr_vs_HT_trt", 2, -0.5, 1.5, 200, -2, 2);
647 TH2F *resAr_vs_HT_ba = new TH2F("resAr_vs_HT_ba", "resAr_vs_HT_ba", 2, -0.5, 1.5, 200, -2, 2);
648 TH2F *resAr_vs_HT_bc = new TH2F("resAr_vs_HT_bc", "resAr_vs_HT_bc", 2, -0.5, 1.5, 200, -2, 2);
649 TH2F *resAr_vs_HT_ea = new TH2F("resAr_vs_HT_ea", "resAr_vs_HT_ea", 2, -0.5, 1.5, 200, -2, 2);
650 TH2F *resAr_vs_HT_ec = new TH2F("resAr_vs_HT_ec", "resAr_vs_HT_ec", 2, -0.5, 1.5, 200, -2, 2);
651
652 // Histograms to store the sin(theta)/p Corrections, TIME RESIDUAL:
653 TH2F *tresAr_vs_SinOverP_trt = new TH2F("tresAr_vs_SinOverP_trt", "tresAr_vs_SinOverP_trt", 24, 0., 1., 200, -20., 20.);
654 TH2F *tresAr_vs_SinOverP_ba = new TH2F("tresAr_vs_SinOverP_ba", "tresAr_vs_SinOverP_ba", 24, 0., 1., 200, -20., 20.);
655 TH2F *tresAr_vs_SinOverP_bc = new TH2F("tresAr_vs_SinOverP_bc", "tresAr_vs_SinOverP_bc", 24, 0., 1., 200, -20., 20.);
656 TH2F *tresAr_vs_SinOverP_ea = new TH2F("tresAr_vs_SinOverP_ea", "tresAr_vs_SinOverP_ea", 24, 0., 1., 200, -20., 20.);
657 TH2F *tresAr_vs_SinOverP_ec = new TH2F("tresAr_vs_SinOverP_ec", "tresAr_vs_SinOverP_ec", 24, 0., 1., 200, -20., 20.);
658
659 // Residual:
660 TH2F *resAr_vs_SinOverP_trt = new TH2F("resAr_vs_SinOverP_trt", "resAr_vs_SinOverP_trt", 24, 0., 1., 200, -2., 2.);
661 TH2F *resAr_vs_SinOverP_ba = new TH2F("resAr_vs_SinOverP_ba", "resAr_vs_SinOverP_ba", 24, 0., 1., 200, -2., 2.);
662 TH2F *resAr_vs_SinOverP_bc = new TH2F("resAr_vs_SinOverP_bc", "resAr_vs_SinOverP_bc", 24, 0., 1., 200, -2., 2.);
663 TH2F *resAr_vs_SinOverP_ea = new TH2F("resAr_vs_SinOverP_ea", "resAr_vs_SinOverP_ea", 24, 0., 1., 200, -2., 2.);
664 TH2F *resAr_vs_SinOverP_ec = new TH2F("resAr_vs_SinOverP_ec", "resAr_vs_SinOverP_ec", 24, 0., 1., 200, -2., 2.);
665
666 // Histograms to store the sin(theta)/p Corrections, TIME RESIDUAL:
667 TH2F *tresAr_vs_CosOverP_trt = new TH2F("tresAr_vs_CosOverP_trt", "tresAr_vs_CosOverP_trt", 24, 0., 1., 200, -20., 20.);
668 TH2F *tresAr_vs_CosOverP_ba = new TH2F("tresAr_vs_CosOverP_ba", "tresAr_vs_CosOverP_ba", 24, 0., 1., 200, -20., 20.);
669 TH2F *tresAr_vs_CosOverP_bc = new TH2F("tresAr_vs_CosOverP_bc", "tresAr_vs_CosOverP_bc", 24, 0., 1., 200, -20., 20.);
670 TH2F *tresAr_vs_CosOverP_ea = new TH2F("tresAr_vs_CosOverP_ea", "tresAr_vs_CosOverP_ea", 24, 0., 1., 200, -20., 20.);
671 TH2F *tresAr_vs_CosOverP_ec = new TH2F("tresAr_vs_CosOverP_ec", "tresAr_vs_CosOverP_ec", 24, 0., 1., 200, -20., 20.);
672
673 // Residual:
674 TH2F *resAr_vs_CosOverP_trt = new TH2F("resAr_vs_CosOverP_trt", "resAr_vs_CosOverP_trt", 24, 0., 1., 200, -2., 2.);
675 TH2F *resAr_vs_CosOverP_ba = new TH2F("resAr_vs_CosOverP_ba", "resAr_vs_CosOverP_ba", 24, 0., 1., 200, -2., 2.);
676 TH2F *resAr_vs_CosOverP_bc = new TH2F("resAr_vs_CosOverP_bc", "resAr_vs_CosOverP_bc", 24, 0., 1., 200, -2., 2.);
677 TH2F *resAr_vs_CosOverP_ea = new TH2F("resAr_vs_CosOverP_ea", "resAr_vs_CosOverP_ea", 24, 0., 1., 200, -2., 2.);
678 TH2F *resAr_vs_CosOverP_ec = new TH2F("resAr_vs_CosOverP_ec", "resAr_vs_CosOverP_ec", 24, 0., 1., 200, -2., 2.);
679
680 // for (map<string,calibpars>::iterator iot=pt0map.begin(); iot!=pt0map.end(); iot++) {
681 // cout << iot->second.t0 << " " << iot->second.dt0 << endl;
682 // }
683
684 //*********************************************************************
685 // FILE LOOP
686 for (int ifiles = 0; ifiles < nfiles; ifiles++)
687 {
688 //*********************************************************************
689
690 bool isntuple = false;
691 char filetype[5];
692 ifstream myFile(argv[ifiles + 2], ios::in | ios::binary);
693 if (myFile.is_open())
694 {
695 myFile.read((char *)filetype, 4);
696 filetype[4] = 0;
697 if (strcmp(filetype, "root") == 0)
698 isntuple = true;
699 // if (strcmp(filetype,"bina")==0) return 3;
700 myFile.close();
701 }
702 else
703 {
704 std::cerr << "WRONG INPUT FILE!\n";
705 return -1;
706 }
707
708 //*********************************************************************
709 // SCANNING ROOT FILE
710 if (isntuple)
711 {
712 //*********************************************************************
713
714 nhistsadd = 0;
715
716 TFile *ntfile = new TFile(argv[ifiles + 2]);
717 TNtuple *hittuple = (TNtuple *)ntfile->Get("ntuple");
718
719 hittuple->SetBranchAddress("run", &run);
720 hittuple->SetBranchAddress("evt", &evt);
721 hittuple->SetBranchAddress("trk", &trk);
722 hittuple->SetBranchAddress("det", &det);
723 hittuple->SetBranchAddress("lay", &lay);
724 hittuple->SetBranchAddress("mod", &mod);
725 hittuple->SetBranchAddress("brd", &brd);
726 hittuple->SetBranchAddress("chp", &chp);
727 hittuple->SetBranchAddress("sid", &fsid);
728 hittuple->SetBranchAddress("stl", &stl);
729 hittuple->SetBranchAddress("stw", &stw);
730 hittuple->SetBranchAddress("r", &r);
731 hittuple->SetBranchAddress("dr", &dr);
732 hittuple->SetBranchAddress("dr", &dr);
733 hittuple->SetBranchAddress("t", &t);
734 hittuple->SetBranchAddress("t0", &t0);
735 // hittuple->SetBranchAddress("rtrack",&rtrack);
736 // hittuple->SetBranchAddress("ttrack",&ttrack);
737 hittuple->SetBranchAddress("rtrackunbias", &rtrack);
738 hittuple->SetBranchAddress("drrtrackunbias", &drtrack);
739 hittuple->SetBranchAddress("ttrackunbias", &ttrack);
740 hittuple->SetBranchAddress("ephase", &ephase);
741
742 if (hittuple->GetListOfBranches()->FindObject("theta"))
743 {
744 hittuple->SetBranchAddress("theta", &theta);
745 }
746 if (hittuple->GetListOfBranches()->FindObject("phi"))
747 {
748 hittuple->SetBranchAddress("phi", &phi);
749 }
750 if (hittuple->GetListOfBranches()->FindObject("d0"))
751 {
752 hittuple->SetBranchAddress("d0", &d0);
753 }
754 if (hittuple->GetListOfBranches()->FindObject("pt"))
755 {
756 hittuple->SetBranchAddress("pt", &pt);
757 }
758 if (hittuple->GetListOfBranches()->FindObject("ToT"))
759 {
760 hittuple->SetBranchAddress("ToT", &ToT);
761 }
762 else
763 ToT = 0;
764 if (hittuple->GetListOfBranches()->FindObject("HT"))
765 {
766 hittuple->SetBranchAddress("HT", &HT);
767 }
768 else
769 HT = 0;
770 if (hittuple->GetListOfBranches()->FindObject("qoverp"))
771 {
772 hittuple->SetBranchAddress("qoverp", &qoverp);
773 }
774 else
775 qoverp = 0;
776
777 if (hittuple->GetListOfBranches()->FindObject("isArgonStraw"))
778 {
779 hittuple->SetBranchAddress("isArgonStraw", &isArgonStraw);
780 }
781 else
782 isArgonStraw = 0;
783
784 int nevents;
785 int npt0 = 0;
786 float diffpt0 = 0;
787 nevents = hittuple->GetEntries();
788
789 cout << "SCANNING ROOT FILE " << argv[ifiles + 2] << " (" << nevents << " HITS)" << endl;
790
791 //*********************************************************************
792 // EVENT LOOP
793 for (int ievt = 0; ievt < nevents; ievt++)
794 {
795 // for (int ievt=25000000;ievt<nevents;ievt++){
796 // for (int ievt=0;ievt<2000100;ievt++){
797 // for (int ievt=0;ievt<1000;ievt++){
798 //*********************************************************************
799
800 hittuple->GetEntry(ievt);
801
802 if (dr < 1.0)
803 { // event selection
804 // if (true){ //event selection
805
806 nhits++;
807
808 string trkkey = string(Form("_%i_%i_%i", (int)run, (int)evt, (int)trk));
809
810 string pt0key = string(Form("_%i_%i_%i_%i_%i", (int)det, (int)lay, (int)mod, (int)stl, (int)stw));
811 if (pt0map.find(pt0key) != pt0map.end())
812 {
813 diffpt0 += fabs(t0 - pt0map[string(Form("_%i_%i_%i_%i_%i", (int)det, (int)lay, (int)mod, (int)stl, (int)stw))].t0);
814 t0 = pt0map[string(Form("_%i_%i_%i_%i_%i", (int)det, (int)lay, (int)mod, (int)stl, (int)stw))].t0;
815 npt0++;
816 // cout << string(Form("_%i_%i_%i_%i_%i",(int)det,(int)lay,(int)mod,(int)stl,(int)stw)) << " " << t0 << " -> " << pt0map[string(Form("_%i_%i_%i_%i_%i",(int)det,(int)lay,(int)mod,(int)stl,(int)stw))].t0 << endl;
817 }
818
819 if (trackmap.find(trkkey) == trackmap.end())
820 {
821 trackmap[trkkey].run = run;
822 trackmap[trkkey].event = evt;
823 trackmap[trkkey].track = trk;
824 trackmap[trkkey].epnew = (t + ephase) - t0 - ttrack;
825 trackmap[trkkey].epold = ephase;
826 trackmap[trkkey].nhits = 1.0;
827 trackmap[trkkey].t = t;
828 trackmap[trkkey].ttrack = ttrack;
829 trackmap[trkkey].t0 = t0;
830 trackmap[trkkey].theta = theta;
831 trackmap[trkkey].phi = phi;
832 trackmap[trkkey].d0 = d0;
833 trackmap[trkkey].pt = pt;
834 trackmap[trkkey].trackres = (r - rtrack) * (r - rtrack);
835 trackmap[trkkey].trackresMean = (r - rtrack);
836 }
837 else
838 {
839 trackmap[trkkey].epnew += (t + ephase) - t0 - ttrack;
840 trackmap[trkkey].t += t;
841 trackmap[trkkey].ttrack += ttrack;
842 trackmap[trkkey].t0 += t0;
843 trackmap[trkkey].nhits++;
844 trackmap[trkkey].trackres += (r - rtrack) * (r - rtrack);
845 trackmap[trkkey].trackresMean += (r - rtrack);
846 }
847 sid = (int)fsid;
848
849 short tresbin = Simple1dHist(mintres, maxtres, ntresbins, t - ttrack - t0 + ephase);
850 short rresbin = Simple1dHist(minres, maxres, nresbins, r - rtrack);
851 short trbin = Simple2dHist(mint, maxt, ntbins, minr, maxr, nrbins, t - t0 + ephase, fabs(rtrack));
852
853 if ((int)det == -2 && (int)lay == 0 && (int)mod == 0 && (int)brd == 1)
854 residual_ref->Fill(r - rtrack);
855 if ((int)det == 1 && (int)lay == 0 && (int)mod == 0 && (int)brd == 0 && (int)chp == 1)
856 timeresidual_ref->Fill(t - ttrack - t0 + ephase);
857 if ((int)det == -1 && (int)lay == 2)
858 rt_ref->Fill(t - t0 + ephase, fabs(rtrack));
859 if ((int)det == -1)
860 rt_ref2->Fill(t - t0 + ephase, fabs(rtrack));
861
862 if (histmap.find(sid) == histmap.end())
863 {
864 nhistsadd++;
865 nhists++;
866 histmap[sid] = new CalHist();
867 }
868 if (tresbin >= 0)
869 {
870 nbinsadd += histmap[sid]->IncreaseBin(tresbin, 1);
871 ntres++;
872 }
873 if (rresbin >= 0)
874 {
875 nbinsadd += histmap[sid]->IncreaseBin(rresbin + 100, 1);
876 nres++;
877 }
878 if (trbin >= 0)
879 {
880 nbinsadd += histmap[sid]->IncreaseBin(trbin + 200, 1);
881 nrt++;
882 }
883
884 if (histmap[sid]->maxvalue > maxvalue)
885 maxvalue = histmap[sid]->maxvalue;
886
887 if ((nhits % 1000000 == 0) | (ievt == nevents - 1))
888 {
889 process_mem_usage(vm, rss);
890 // cppcheck-suppress invalidPrintfArgType_uint
891 printf("%8i HITS READ, %8i HISTOGRAMS (%8lu ADDED), %8i BINS ADDED, MAXVALUE: %3i, VM: %8.0f, RSS: %8.0f\n", nhits, nhists, histmap.size() - nhiststmp, nbinsadd, maxvalue, vm, rss);
892 nhiststmp = histmap.size();
893 // nhistsadd=0;
894 nbinsadd = 0;
895 }
896
897 if (nhits % 20000000 == 0)
898 {
899 dump_hists(&histmap, ntbins, nrbins, ntres, argv[1], 0);
900 for (map<int, CalHist *>::iterator it = histmap.begin(); it != histmap.end(); ++it)
901 {
902 delete it->second;
903 }
904 histmap.clear();
905 dump_tracktuple(&trackmap);
906 trackmap.clear();
907 nhiststmp = 0;
908 maxvalue = 0;
909 }
910
911 if (isArgonStraw == 0)
912 { // ONLY STORE THE XENON CIRCLES:
913 reshist_trt->Fill(t - t0, fabs(r) - fabs(rtrack));
914 if ((int)det == -2)
915 reshist1->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
916 if ((int)det == -1)
917 reshist2->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
918 if ((int)det == 1)
919 reshist3->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
920 if ((int)det == 2)
921 reshist4->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
922 if (abs((int)det) == 1)
923 reshist_bar->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
924
925 treshist_trt->Fill(t - t0 - ttrack, fabs(rtrack));
926 if ((int)det == -2)
927 treshist1->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
928 if ((int)det == -1)
929 treshist2->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
930 if ((int)det == 1)
931 treshist3->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
932 if ((int)det == 2)
933 treshist4->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
934 if (abs((int)det) == 1)
935 treshist_bar->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
936
937 float residual = r - rtrack;
938 // float timebin = int((t-t0 + ephase)/3.125 );
939 float time = t - t0 + ephase;
940 float pull = (r - rtrack) / TMath::Sqrt(dr * dr + drtrack * drtrack);
941 float tresidual = t - t0 - ttrack + ephase;
942
943 ToT = ToT / 3.125;
944
945 float sinoverp = 0.0;
946 float cosoverp = 0.0;
947 if (pt != 0.0)
948 {
949 sinoverp = 1000. * sin(theta) * sin(theta) / pt;
950 cosoverp = 1000. * cos(theta) * cos(theta) / pt;
951 }
952
953 float ppt = pt / 1000.;
954 if (ppt > 16)
955 ppt = 15.5;
956
957 tres_vs_ToT_trt->Fill(ToT, tresidual);
958 res_vs_ToT_trt->Fill(ToT, residual);
959 tres_vs_HT_trt->Fill(HT, tresidual);
960 res_vs_HT_trt->Fill(HT, residual);
961 tres_vs_SinOverP_trt->Fill(sinoverp, tresidual);
962 res_vs_SinOverP_trt->Fill(sinoverp, residual);
963 tres_vs_CosOverP_trt->Fill(cosoverp, tresidual);
964 res_vs_CosOverP_trt->Fill(cosoverp, residual);
965
966 pull_trtP->Fill(ppt, pull);
967 residual_trtP->Fill(ppt, residual);
968 tresidual_trtP->Fill(ppt, tresidual);
969 trackerrors_trtP->Fill(ppt, drtrack);
970 if (ppt > ptcutforerrors)
971 pull_vs_tb_trt->Fill(time, pull);
972 if (ppt > ptcutforerrors)
973 errors_vs_tb_trt->Fill(time, dr);
974 if (ppt > ptcutforerrors)
975 trackerrors_vs_tb_trt->Fill(time, drtrack);
976 if (ppt > ptcutforerrors)
977 residual_vs_tb_trt->Fill(time, residual);
978
979 if ((int)det == -2)
980 {
981 tres_vs_ToT_ea->Fill(ToT, tresidual);
982 res_vs_ToT_ea->Fill(ToT, residual);
983 tres_vs_HT_ea->Fill(HT, tresidual);
984 res_vs_HT_ea->Fill(HT, residual);
985 tres_vs_SinOverP_ea->Fill(sinoverp, tresidual);
986 res_vs_SinOverP_ea->Fill(sinoverp, residual);
987 tres_vs_CosOverP_ea->Fill(cosoverp, tresidual);
988 res_vs_CosOverP_ea->Fill(cosoverp, residual);
989
990 pull_eaP->Fill(ppt, pull);
991 residual_eaP->Fill(ppt, residual);
992 tresidual_eaP->Fill(ppt, tresidual);
993 trackerrors_eaP->Fill(ppt, drtrack);
994 if (ppt > ptcutforerrors)
995 pull_vs_tb_ea->Fill(time, pull);
996 if (ppt > ptcutforerrors)
997 errors_vs_tb_ea->Fill(time, dr);
998 if (ppt > ptcutforerrors)
999 trackerrors_vs_tb_ea->Fill(time, drtrack);
1000 if (ppt > ptcutforerrors)
1001 residual_vs_tb_ea->Fill(time, residual);
1002 }
1003 if ((int)det == 2)
1004 {
1005 tres_vs_ToT_ec->Fill(ToT, tresidual);
1006 res_vs_ToT_ec->Fill(ToT, residual);
1007 tres_vs_HT_ec->Fill(HT, tresidual);
1008 res_vs_HT_ec->Fill(HT, residual);
1009 tres_vs_SinOverP_ec->Fill(sinoverp, tresidual);
1010 res_vs_SinOverP_ec->Fill(sinoverp, residual);
1011 tres_vs_CosOverP_ec->Fill(cosoverp, tresidual);
1012 res_vs_CosOverP_ec->Fill(cosoverp, residual);
1013
1014 pull_ecP->Fill(ppt, pull);
1015 residual_ecP->Fill(ppt, residual);
1016 tresidual_ecP->Fill(ppt, tresidual);
1017 trackerrors_ecP->Fill(ppt, drtrack);
1018 if (ppt > ptcutforerrors)
1019 pull_vs_tb_ec->Fill(time, pull);
1020 if (ppt > ptcutforerrors)
1021 errors_vs_tb_ec->Fill(time, dr);
1022 if (ppt > ptcutforerrors)
1023 trackerrors_vs_tb_ec->Fill(time, drtrack);
1024 if (ppt > ptcutforerrors)
1025 residual_vs_tb_ec->Fill(time, residual);
1026 }
1027 if ((int)det == -1)
1028 {
1029 tres_vs_ToT_ba->Fill(ToT, tresidual);
1030 res_vs_ToT_ba->Fill(ToT, residual);
1031 tres_vs_HT_ba->Fill(HT, tresidual);
1032 res_vs_HT_ba->Fill(HT, residual);
1033 tres_vs_SinOverP_ba->Fill(sinoverp, tresidual);
1034 res_vs_SinOverP_ba->Fill(sinoverp, residual);
1035 tres_vs_CosOverP_ba->Fill(cosoverp, tresidual);
1036 res_vs_CosOverP_ba->Fill(cosoverp, residual);
1037
1038 pull_baP->Fill(ppt, pull);
1039 residual_baP->Fill(ppt, residual);
1040 tresidual_baP->Fill(ppt, tresidual);
1041 trackerrors_baP->Fill(ppt, drtrack);
1042 if (ppt > ptcutforerrors)
1043 pull_vs_tb_ba->Fill(time, pull);
1044 if (ppt > ptcutforerrors)
1045 errors_vs_tb_ba->Fill(time, dr);
1046 if (ppt > ptcutforerrors)
1047 trackerrors_vs_tb_ba->Fill(time, drtrack);
1048 if (ppt > ptcutforerrors)
1049 residual_vs_tb_ba->Fill(time, residual);
1050 }
1051 if ((int)det == 1)
1052 {
1053 tres_vs_ToT_bc->Fill(ToT, tresidual);
1054 res_vs_ToT_bc->Fill(ToT, residual);
1055 tres_vs_HT_bc->Fill(HT, tresidual);
1056 res_vs_HT_bc->Fill(HT, residual);
1057 tres_vs_SinOverP_bc->Fill(sinoverp, tresidual);
1058 res_vs_SinOverP_bc->Fill(sinoverp, residual);
1059 tres_vs_CosOverP_bc->Fill(cosoverp, tresidual);
1060 res_vs_CosOverP_bc->Fill(cosoverp, residual);
1061
1062 pull_bcP->Fill(ppt, pull);
1063 residual_bcP->Fill(ppt, residual);
1064 tresidual_bcP->Fill(ppt, tresidual);
1065 trackerrors_bcP->Fill(ppt, drtrack);
1066 if (ppt > ptcutforerrors)
1067 pull_vs_tb_bc->Fill(time, pull);
1068 if (ppt > ptcutforerrors)
1069 errors_vs_tb_bc->Fill(time, dr);
1070 if (ppt > ptcutforerrors)
1071 trackerrors_vs_tb_bc->Fill(time, drtrack);
1072 if (ppt > ptcutforerrors)
1073 residual_vs_tb_bc->Fill(time, residual);
1074 }
1075 }
1076 else
1077 { // DUMP THE ARGON PART:
1078
1079 existAr = true;
1080
1081 resArhist_trt->Fill(t - t0, fabs(r) - fabs(rtrack));
1082 if ((int)det == -2)
1083 resArhist1->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
1084 if ((int)det == -1)
1085 resArhist2->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
1086 if ((int)det == 1)
1087 resArhist3->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
1088 if ((int)det == 2)
1089 resArhist4->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
1090 if (abs((int)det) == 1)
1091 resArhist_bar->Fill(t - t0 + ephase, fabs(r) - fabs(rtrack));
1092
1093 tresArhist_trt->Fill(t - t0 - ttrack, fabs(rtrack));
1094 if ((int)det == -2)
1095 tresArhist1->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
1096 if ((int)det == -1)
1097 tresArhist2->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
1098 if ((int)det == 1)
1099 tresArhist3->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
1100 if ((int)det == 2)
1101 tresArhist4->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
1102 if (abs((int)det) == 1)
1103 tresArhist_bar->Fill(t - t0 - ttrack + ephase, fabs(rtrack));
1104
1105 float residual = r - rtrack;
1106 // float timebin = int((t-t0 + ephase)/3.125 );
1107 float time = t - t0 + ephase;
1108 float pull = (r - rtrack) / TMath::Sqrt(dr * dr + drtrack * drtrack);
1109 float tresidual = t - t0 - ttrack + ephase;
1110
1111 ToT = ToT / 3.125;
1112 float sinoverp = 0.0;
1113 float cosoverp = 0.0;
1114 if (pt != 0.0)
1115 {
1116 sinoverp = 1000. * sin(theta) * sin(theta) / pt;
1117 cosoverp = 1000. * cos(theta) * cos(theta) / pt;
1118 }
1119 float ppt = pt / 1000.;
1120 if (ppt > 16)
1121 ppt = 15.5;
1122
1123 tresAr_vs_ToT_trt->Fill(ToT, tresidual);
1124 resAr_vs_ToT_trt->Fill(ToT, residual);
1125 tresAr_vs_HT_trt->Fill(HT, tresidual);
1126 resAr_vs_HT_trt->Fill(HT, residual);
1127 tresAr_vs_SinOverP_trt->Fill(sinoverp, tresidual);
1128 resAr_vs_SinOverP_trt->Fill(sinoverp, residual);
1129 tresAr_vs_CosOverP_trt->Fill(cosoverp, tresidual);
1130 resAr_vs_CosOverP_trt->Fill(cosoverp, residual);
1131
1132 pullAr_trtP->Fill(ppt, pull);
1133 residualAr_trtP->Fill(ppt, residual);
1134 tresidualAr_trtP->Fill(ppt, tresidual);
1135 trackerrorsAr_trtP->Fill(ppt, drtrack);
1136 if (ppt > ptcutforerrors)
1137 pullAr_vs_tb_trt->Fill(time, pull);
1138 if (ppt > ptcutforerrors)
1139 errors_vs_tb_trt->Fill(time, dr);
1140 if (ppt > ptcutforerrors)
1141 trackerrors_vs_tb_trt->Fill(time, drtrack);
1142 if (ppt > ptcutforerrors)
1143 residualAr_vs_tb_trt->Fill(time, residual);
1144
1145 if ((int)det == 2)
1146 {
1147 tresAr_vs_ToT_ea->Fill(ToT, tresidual);
1148 resAr_vs_ToT_ea->Fill(ToT, residual);
1149 tresAr_vs_HT_ea->Fill(HT, tresidual);
1150 resAr_vs_HT_ea->Fill(HT, residual);
1151 tresAr_vs_SinOverP_ea->Fill(sinoverp, tresidual);
1152 resAr_vs_SinOverP_ea->Fill(sinoverp, residual);
1153 tresAr_vs_CosOverP_ea->Fill(cosoverp, tresidual);
1154 resAr_vs_CosOverP_ea->Fill(cosoverp, residual);
1155
1156 pullAr_eaP->Fill(ppt, pull);
1157 residualAr_eaP->Fill(ppt, residual);
1158 tresidualAr_eaP->Fill(ppt, tresidual);
1159 trackerrorsAr_eaP->Fill(ppt, drtrack);
1160 if (ppt > ptcutforerrors)
1161 pullAr_vs_tb_ea->Fill(time, pull);
1162 if (ppt > ptcutforerrors)
1163 errorsAr_vs_tb_ea->Fill(time, dr);
1164 if (ppt > ptcutforerrors)
1165 trackerrorsAr_vs_tb_ea->Fill(time, drtrack);
1166 if (ppt > ptcutforerrors)
1167 residualAr_vs_tb_ea->Fill(time, residual);
1168 }
1169 else if ((int)det == -2)
1170 {
1171 tresAr_vs_ToT_ec->Fill(ToT, tresidual);
1172 resAr_vs_ToT_ec->Fill(ToT, residual);
1173 tresAr_vs_HT_ec->Fill(HT, tresidual);
1174 resAr_vs_HT_ec->Fill(HT, residual);
1175 tresAr_vs_SinOverP_ec->Fill(sinoverp, tresidual);
1176 resAr_vs_SinOverP_ec->Fill(sinoverp, residual);
1177 tresAr_vs_CosOverP_ec->Fill(cosoverp, tresidual);
1178 resAr_vs_CosOverP_ec->Fill(cosoverp, residual);
1179
1180 pullAr_ecP->Fill(ppt, pull);
1181 residualAr_ecP->Fill(ppt, residual);
1182 tresidualAr_ecP->Fill(ppt, tresidual);
1183 trackerrorsAr_ecP->Fill(ppt, drtrack);
1184 if (ppt > ptcutforerrors)
1185 pullAr_vs_tb_ec->Fill(time, pull);
1186 if (ppt > ptcutforerrors)
1187 errorsAr_vs_tb_ec->Fill(time, dr);
1188 if (ppt > ptcutforerrors)
1189 trackerrorsAr_vs_tb_ec->Fill(time, drtrack);
1190 if (ppt > ptcutforerrors)
1191 residualAr_vs_tb_ec->Fill(time, residual);
1192 }
1193 else if ((int)det == 1)
1194 {
1195 tresAr_vs_ToT_ba->Fill(ToT, tresidual);
1196 resAr_vs_ToT_ba->Fill(ToT, residual);
1197 tresAr_vs_HT_ba->Fill(HT, tresidual);
1198 resAr_vs_HT_ba->Fill(HT, residual);
1199 tresAr_vs_SinOverP_ba->Fill(sinoverp, tresidual);
1200 resAr_vs_SinOverP_ba->Fill(sinoverp, residual);
1201 tresAr_vs_CosOverP_ba->Fill(cosoverp, tresidual);
1202 resAr_vs_CosOverP_ba->Fill(cosoverp, residual);
1203
1204 pullAr_baP->Fill(ppt, pull);
1205 residualAr_baP->Fill(ppt, residual);
1206 tresidualAr_baP->Fill(ppt, tresidual);
1207 trackerrorsAr_baP->Fill(ppt, drtrack);
1208 if (ppt > ptcutforerrors)
1209 pullAr_vs_tb_ba->Fill(time, pull);
1210 if (ppt > ptcutforerrors)
1211 errorsAr_vs_tb_ba->Fill(time, dr);
1212 if (ppt > ptcutforerrors)
1213 trackerrorsAr_vs_tb_ba->Fill(time, drtrack);
1214 if (ppt > ptcutforerrors)
1215 residualAr_vs_tb_ba->Fill(time, residual);
1216 }
1217 else if ((int)det == -1)
1218 {
1219 tresAr_vs_ToT_bc->Fill(ToT, tresidual);
1220 resAr_vs_ToT_bc->Fill(ToT, residual);
1221 tresAr_vs_HT_bc->Fill(HT, tresidual);
1222 resAr_vs_HT_bc->Fill(HT, residual);
1223 tresAr_vs_SinOverP_bc->Fill(sinoverp, tresidual);
1224 resAr_vs_SinOverP_bc->Fill(sinoverp, residual);
1225 tresAr_vs_CosOverP_bc->Fill(cosoverp, tresidual);
1226 resAr_vs_CosOverP_bc->Fill(cosoverp, residual);
1227
1228 pullAr_bcP->Fill(ppt, pull);
1229 residualAr_bcP->Fill(ppt, residual);
1230 tresidualAr_bcP->Fill(ppt, tresidual);
1231 trackerrorsAr_bcP->Fill(ppt, drtrack);
1232 if (ppt > ptcutforerrors)
1233 pullAr_vs_tb_bc->Fill(time, pull);
1234 if (ppt > ptcutforerrors)
1235 errorsAr_vs_tb_bc->Fill(time, dr);
1236 if (ppt > ptcutforerrors)
1237 trackerrorsAr_vs_tb_bc->Fill(time, drtrack);
1238 if (ppt > ptcutforerrors)
1239 residualAr_vs_tb_bc->Fill(time, residual);
1240 }
1241 } // Close the argon part
1242
1243 } // ... if event selection (dr<1) -> Precision hits:
1244
1245 float ppt = pt / 1000.;
1246 if (ppt > 16)
1247 ppt = 15.5;
1248 float residual = r - rtrack;
1249 float tresidual = t - t0 - ttrack + ephase;
1250 float pull = (r - rtrack) / TMath::Sqrt(dr * dr + drtrack * drtrack);
1251
1252 if (isArgonStraw == 0)
1253 { // ONLY STORE THE XENON CIRCLES:
1254 pull_trt->Fill(ppt, pull);
1255 residual_trt->Fill(ppt, residual);
1256 tresidual_trt->Fill(ppt, residual);
1257 if ((int)det == -2)
1258 {
1259 pull_ea->Fill(ppt, pull);
1260 tresidual_ea->Fill(ppt, tresidual);
1261 }
1262 if ((int)det == 2)
1263 {
1264 pull_ec->Fill(ppt, pull);
1265 residual_ec->Fill(ppt, residual);
1266 tresidual_ec->Fill(ppt, tresidual);
1267 }
1268 if ((int)det == -1)
1269 {
1270 pull_ba->Fill(ppt, pull);
1271 residual_ba->Fill(ppt, residual);
1272 tresidual_ba->Fill(ppt, tresidual);
1273 }
1274 if ((int)det == 1)
1275 {
1276 pull_bc->Fill(ppt, pull);
1277 residual_bc->Fill(ppt, residual);
1278 tresidual_bc->Fill(ppt, tresidual);
1279 }
1280 }
1281 else
1282 {
1283 pullAr_trt->Fill(ppt, pull);
1284 residualAr_trt->Fill(ppt, residual);
1285 tresidualAr_trt->Fill(ppt, residual);
1286 if ((int)det == 2)
1287 {
1288 pullAr_ea->Fill(ppt, pull);
1289 tresidualAr_ea->Fill(ppt, tresidual);
1290 }
1291 if ((int)det == -2)
1292 {
1293 pullAr_ec->Fill(ppt, pull);
1294 residualAr_ec->Fill(ppt, residual);
1295 tresidualAr_ec->Fill(ppt, tresidual);
1296 }
1297 if ((int)det == 1)
1298 {
1299 pullAr_ba->Fill(ppt, pull);
1300 residualAr_ba->Fill(ppt, residual);
1301 tresidualAr_ba->Fill(ppt, tresidual);
1302 }
1303 if ((int)det == -1)
1304 {
1305 pullAr_bc->Fill(ppt, pull);
1306 residualAr_bc->Fill(ppt, residual);
1307 tresidualAr_bc->Fill(ppt, tresidual);
1308 }
1309 } // save all AR hits
1310
1311 } // END OF EVENT LOOP
1312 if (npt0 !=0){
1313 printf("REPLACED %d T0S WITH PRECISION VALUES (MEAN DIFF=%f)\n", npt0, diffpt0 / float(npt0));
1314 } else {
1315 printf("ERROR!! npt0 is zero!");
1316 }
1317
1318 if (ifiles == nfiles - 1)
1319 {
1320 // Since we access here in the last loop iteration, better to use RECREATE than UPDATE
1321 TFile *ttfile = new TFile("tracktuple.root", "RECREATE");
1322
1323 TDirectory *trtdir = ttfile->mkdir("TRT_all");
1324 TDirectory *binhist = trtdir->mkdir("reshists");
1325 TDirectory *detdir1 = trtdir->mkdir("Detector_-2");
1326 TDirectory *binhist1 = detdir1->mkdir("reshists1"); // NB! changed reshist to reshist1 for safety
1327 TDirectory *detdir2 = trtdir->mkdir("Detector_-1");
1328 TDirectory *binhist2 = detdir2->mkdir("reshists2");
1329 TDirectory *detdir3 = trtdir->mkdir("Detector_1");
1330 TDirectory *binhist3 = detdir3->mkdir("reshists3");
1331 TDirectory *detdir4 = trtdir->mkdir("Detector_2");
1332 TDirectory *binhist4 = detdir4->mkdir("reshists4");
1333 TDirectory *detdir5 = trtdir->mkdir("WholeBarrel_1");
1334 TDirectory *binhist5 = detdir5->mkdir("reshists5");
1335
1336 TDirectory *errordir = ttfile->mkdir("Errors");
1337 TDirectory *corrdir = ttfile->mkdir("Correction");
1338
1339 int npointsX = reshist1->GetNbinsX();
1340 int npointsY = treshist1->GetNbinsY();
1341 TH1D **hslizesX = new TH1D *[npointsX];
1342 TH1D **hslizesY = new TH1D *[npointsY];
1343
1344 TVectorD tbins(npointsX);
1345 TVectorD rbins(npointsY);
1346
1347 string chname;
1348 for (int i = 0; i < npointsX; i++)
1349 {
1350
1351 tbins(i) = reshist1->GetXaxis()->GetBinCenter(i + 1);
1352
1353 binhist->cd();
1354 chname = string(Form("res_tbin%i_trt", i));
1355 hslizesX[i] = reshist_trt->ProjectionY(chname.data(), i + 1, i + 1);
1356 binhist1->cd();
1357 chname = string(Form("res_tbin%i_-2", i));
1358 hslizesX[i] = reshist1->ProjectionY(chname.data(), i + 1, i + 1);
1359 binhist2->cd();
1360 chname = string(Form("res_tbin%i_-1", i));
1361 hslizesX[i] = reshist2->ProjectionY(chname.data(), i + 1, i + 1);
1362 binhist3->cd();
1363 chname = string(Form("res_tbin%i_1", i));
1364 hslizesX[i] = reshist3->ProjectionY(chname.data(), i + 1, i + 1);
1365 binhist4->cd();
1366 chname = string(Form("res_tbin%i_2", i));
1367 hslizesX[i] = reshist4->ProjectionY(chname.data(), i + 1, i + 1);
1368 binhist5->cd();
1369 chname = string(Form("res_tbin%i_bar", i));
1370 hslizesX[i] = reshist_bar->ProjectionY(chname.data(), i + 1, i + 1);
1371 }
1372
1373 for (int i = 0; i < npointsY; i++)
1374 {
1375
1376 rbins(i) = treshist1->GetYaxis()->GetBinCenter(i + 1);
1377
1378 binhist->cd();
1379 chname = string(Form("tres_rbin%i_trt", i));
1380 hslizesY[i] = treshist_trt->ProjectionX(chname.data(), i + 1, i + 1);
1381 binhist1->cd();
1382 chname = string(Form("tres_rbin%i_-2", i));
1383 hslizesY[i] = treshist1->ProjectionX(chname.data(), i + 1, i + 1);
1384 binhist2->cd();
1385 chname = string(Form("tres_rbin%i_-1", i));
1386 hslizesY[i] = treshist1->ProjectionX(chname.data(), i + 1, i + 1);
1387 binhist3->cd();
1388 chname = string(Form("tres_rbin%i_1", i));
1389 hslizesY[i] = treshist3->ProjectionX(chname.data(), i + 1, i + 1);
1390 binhist4->cd();
1391 chname = string(Form("tres_rbin%i_2", i));
1392 hslizesY[i] = treshist4->ProjectionX(chname.data(), i + 1, i + 1);
1393 binhist5->cd();
1394 chname = string(Form("tres_rbin%i_bar", i));
1395 hslizesY[i] = treshist_bar->ProjectionX(chname.data(), i + 1, i + 1);
1396 }
1397
1398 ttfile->cd();
1399 tbins.Write("tbins");
1400 rbins.Write("rbins");
1401
1402 trtdir->cd();
1403 residual_ref->Write();
1404 timeresidual_ref->Write();
1405 rt_ref->Write();
1406 rt_ref2->Write();
1407
1408 binhist->cd();
1409 reshist_trt->Write();
1410 binhist5->cd();
1411 reshist_bar->Write();
1412 binhist1->cd();
1413 reshist1->Write();
1414 binhist2->cd();
1415 reshist2->Write();
1416 binhist3->cd();
1417 reshist3->Write();
1418 binhist4->cd();
1419 reshist4->Write();
1420 binhist->cd();
1421 treshist_trt->Write();
1422 binhist5->cd();
1423 treshist_bar->Write();
1424 binhist1->cd();
1425 treshist1->Write();
1426 binhist2->cd();
1427 treshist2->Write();
1428 binhist3->cd();
1429 treshist3->Write();
1430 binhist4->cd();
1431 treshist4->Write();
1432
1433 // Save stuff for errors:
1434 errordir->cd();
1435
1436 pull_trt->Write();
1437 pull_ba->Write();
1438 pull_bc->Write();
1439 pull_ea->Write();
1440 pull_ec->Write();
1441
1442 residual_trt->Write();
1443 residual_ba->Write();
1444 residual_bc->Write();
1445 residual_ea->Write();
1446 residual_ec->Write();
1447
1448 tresidual_trt->Write();
1449 tresidual_ba->Write();
1450 tresidual_bc->Write();
1451 tresidual_ea->Write();
1452 tresidual_ec->Write();
1453
1454 pull_trtP->Write();
1455 pull_baP->Write();
1456 pull_bcP->Write();
1457 pull_eaP->Write();
1458 pull_ecP->Write();
1459
1460 residual_trtP->Write();
1461 residual_baP->Write();
1462 residual_bcP->Write();
1463 residual_eaP->Write();
1464 residual_ecP->Write();
1465
1466 tresidual_trtP->Write();
1467 tresidual_baP->Write();
1468 tresidual_bcP->Write();
1469 tresidual_eaP->Write();
1470 tresidual_ecP->Write();
1471
1472 pull_vs_tb_trt->Write();
1473 pull_vs_tb_trt->Write();
1474 pull_vs_tb_ba->Write();
1475 pull_vs_tb_bc->Write();
1476 pull_vs_tb_ea->Write();
1477 pull_vs_tb_ec->Write();
1478
1479 errors_vs_tb_trt->Write();
1480 errors_vs_tb_ba->Write();
1481 errors_vs_tb_bc->Write();
1482 errors_vs_tb_ea->Write();
1483 errors_vs_tb_ec->Write();
1484
1485 residual_vs_tb_trt->Write();
1486 residual_vs_tb_ba->Write();
1487 residual_vs_tb_bc->Write();
1488 residual_vs_tb_ea->Write();
1489 residual_vs_tb_ec->Write();
1490
1491 trackerrors_trtP->Write();
1492 trackerrors_baP->Write();
1493 trackerrors_bcP->Write();
1494 trackerrors_eaP->Write();
1495 trackerrors_ecP->Write();
1496
1497 trackerrors_vs_tb_trt->Write();
1498 trackerrors_vs_tb_ba->Write();
1499 trackerrors_vs_tb_bc->Write();
1500 trackerrors_vs_tb_ea->Write();
1501 trackerrors_vs_tb_ec->Write();
1502
1503 // Save HT correction and ToT
1504 corrdir->cd();
1505
1506 tres_vs_ToT_trt->Write();
1507 tres_vs_ToT_ba->Write();
1508 tres_vs_ToT_bc->Write();
1509 tres_vs_ToT_ea->Write();
1510 tres_vs_ToT_ec->Write();
1511
1512 res_vs_ToT_trt->Write();
1513 res_vs_ToT_ba->Write();
1514 res_vs_ToT_bc->Write();
1515 res_vs_ToT_ea->Write();
1516 res_vs_ToT_ec->Write();
1517
1518 tres_vs_HT_trt->Write();
1519 tres_vs_HT_ba->Write();
1520 tres_vs_HT_bc->Write();
1521 tres_vs_HT_ea->Write();
1522 tres_vs_HT_ec->Write();
1523
1524 res_vs_HT_trt->Write();
1525 res_vs_HT_ba->Write();
1526 res_vs_HT_bc->Write();
1527 res_vs_HT_ea->Write();
1528 res_vs_HT_ec->Write();
1529
1530 tres_vs_SinOverP_trt->Write();
1531 tres_vs_SinOverP_ba->Write();
1532 tres_vs_SinOverP_bc->Write();
1533 tres_vs_SinOverP_ea->Write();
1534 tres_vs_SinOverP_ec->Write();
1535
1536 res_vs_SinOverP_trt->Write();
1537 res_vs_SinOverP_ba->Write();
1538 res_vs_SinOverP_bc->Write();
1539 res_vs_SinOverP_ea->Write();
1540 res_vs_SinOverP_ec->Write();
1541
1542 tres_vs_CosOverP_trt->Write();
1543 tres_vs_CosOverP_ba->Write();
1544 tres_vs_CosOverP_bc->Write();
1545 tres_vs_CosOverP_ea->Write();
1546 tres_vs_CosOverP_ec->Write();
1547
1548 res_vs_CosOverP_trt->Write();
1549 res_vs_CosOverP_ba->Write();
1550 res_vs_CosOverP_bc->Write();
1551 res_vs_CosOverP_ea->Write();
1552 res_vs_CosOverP_ec->Write();
1553
1554 // Save the ARGON Part:
1555
1556 if (existAr)
1557 {
1558 TDirectory *trtArdir = ttfile->mkdir("TRT_Ar_all");
1559 TDirectory *binArhist = trtArdir->mkdir("reshists_Ar"); // NB! changed directory names from reshists
1560 TDirectory *detArdir1 = trtArdir->mkdir("Detector_Ar_-2");
1561 TDirectory *binArhist1 = detArdir1->mkdir("reshists_Ar1");
1562 TDirectory *detArdir2 = trtArdir->mkdir("Detector_Ar_-1");
1563 TDirectory *binArhist2 = detArdir2->mkdir("reshists_Ar2");
1564 TDirectory *detArdir3 = trtArdir->mkdir("Detector_Ar_1");
1565 TDirectory *binArhist3 = detArdir3->mkdir("reshists_Ar3");
1566 TDirectory *detArdir4 = trtArdir->mkdir("Detector_Ar_2");
1567 TDirectory *binArhist4 = detArdir4->mkdir("reshists_Ar4");
1568 TDirectory *detArdir5 = trtArdir->mkdir("WholeBarrel_Ar_1");
1569 TDirectory *binArhist5 = detArdir5->mkdir("reshists_Ar5");
1570
1571 TDirectory *errordirAr = ttfile->mkdir("ErrorsAr");
1572 TDirectory *corrdirAr = ttfile->mkdir("CorrectionAr");
1573 trtArdir->cd();
1574
1575 binArhist->cd();
1576 resArhist_trt->Write();
1577 binArhist5->cd();
1578 resArhist_bar->Write();
1579 binArhist1->cd();
1580 resArhist1->Write();
1581 binArhist2->cd();
1582 resArhist2->Write();
1583 binArhist3->cd();
1584 resArhist3->Write();
1585 binArhist4->cd();
1586 resArhist4->Write();
1587 binArhist->cd();
1588 tresArhist_trt->Write();
1589 binArhist5->cd();
1590 tresArhist_bar->Write();
1591 binArhist1->cd();
1592 tresArhist1->Write();
1593 binArhist2->cd();
1594 tresArhist2->Write();
1595 binArhist3->cd();
1596 tresArhist3->Write();
1597 binArhist4->cd();
1598 tresArhist4->Write();
1599
1600 TH1D **hslizesArX = new TH1D *[npointsX];
1601 TH1D **hslizesArY = new TH1D *[npointsY];
1602 TVectorD tbinsAr(npointsX);
1603 TVectorD rbinsAr(npointsY);
1604 for (int i = 0; i < npointsX; i++)
1605 {
1606 tbinsAr(i) = resArhist1->GetXaxis()->GetBinCenter(i + 1);
1607 binArhist->cd();
1608 chname = string(Form("res_tbin%i_trt", i));
1609 hslizesArX[i] = resArhist_trt->ProjectionY(chname.data(), i + 1, i + 1);
1610 binArhist1->cd();
1611 chname = string(Form("res_tbin%i_-2", i));
1612 hslizesArX[i] = resArhist1->ProjectionY(chname.data(), i + 1, i + 1);
1613 binArhist2->cd();
1614 chname = string(Form("res_tbin%i_-1", i));
1615 hslizesArX[i] = resArhist2->ProjectionY(chname.data(), i + 1, i + 1);
1616 binArhist3->cd();
1617 chname = string(Form("res_tbin%i_1", i));
1618 hslizesArX[i] = resArhist3->ProjectionY(chname.data(), i + 1, i + 1);
1619 binArhist4->cd();
1620 chname = string(Form("res_tbin%i_2", i));
1621 hslizesArX[i] = resArhist4->ProjectionY(chname.data(), i + 1, i + 1);
1622 binArhist5->cd();
1623 chname = string(Form("res_tbin%i_bar", i));
1624 hslizesArX[i] = resArhist_bar->ProjectionY(chname.data(), i + 1, i + 1);
1625 }
1626
1627 for (int i = 0; i < npointsY; i++)
1628 {
1629 rbinsAr(i) = tresArhist1->GetYaxis()->GetBinCenter(i + 1);
1630 binArhist->cd();
1631 chname = string(Form("tres_rbin%i_trt", i));
1632 hslizesArY[i] = tresArhist_trt->ProjectionX(chname.data(), i + 1, i + 1);
1633 binArhist1->cd();
1634 chname = string(Form("tres_rbin%i_-2", i));
1635 hslizesArY[i] = tresArhist1->ProjectionX(chname.data(), i + 1, i + 1);
1636 binArhist2->cd();
1637 chname = string(Form("tres_rbin%i_-1", i));
1638 hslizesArY[i] = tresArhist1->ProjectionX(chname.data(), i + 1, i + 1);
1639 binArhist3->cd();
1640 chname = string(Form("tres_rbin%i_1", i));
1641 hslizesArY[i] = tresArhist3->ProjectionX(chname.data(), i + 1, i + 1);
1642 binArhist4->cd();
1643 chname = string(Form("tres_rbin%i_2", i));
1644 hslizesArY[i] = tresArhist4->ProjectionX(chname.data(), i + 1, i + 1);
1645 binArhist5->cd();
1646 chname = string(Form("tres_rbin%i_bar", i));
1647 hslizesArY[i] = tresArhist_bar->ProjectionX(chname.data(), i + 1, i + 1);
1648 }
1649
1650 // Save stuff for errors:
1651 errordirAr->cd();
1652
1653 pullAr_trt->Write();
1654 pullAr_ba->Write();
1655 pullAr_bc->Write();
1656 pullAr_ea->Write();
1657 pullAr_ec->Write();
1658
1659 residualAr_trt->Write();
1660 residualAr_ba->Write();
1661 residualAr_bc->Write();
1662 residualAr_ea->Write();
1663 residualAr_ec->Write();
1664
1665 tresidualAr_trt->Write();
1666 tresidualAr_ba->Write();
1667 tresidualAr_bc->Write();
1668 tresidualAr_ea->Write();
1669 tresidualAr_ec->Write();
1670
1671 pullAr_trtP->Write();
1672 pullAr_baP->Write();
1673 pullAr_bcP->Write();
1674 pullAr_eaP->Write();
1675 pullAr_ecP->Write();
1676
1677 residualAr_trtP->Write();
1678 residualAr_baP->Write();
1679 residualAr_bcP->Write();
1680 residualAr_eaP->Write();
1681 residualAr_ecP->Write();
1682
1683 tresidualAr_trtP->Write();
1684 tresidualAr_baP->Write();
1685 tresidualAr_bcP->Write();
1686 tresidualAr_eaP->Write();
1687 tresidualAr_ecP->Write();
1688
1689 pullAr_vs_tb_trt->Write();
1690 pullAr_vs_tb_trt->Write();
1691 pullAr_vs_tb_ba->Write();
1692 pullAr_vs_tb_bc->Write();
1693 pullAr_vs_tb_ea->Write();
1694 pullAr_vs_tb_ec->Write();
1695
1696 errorsAr_vs_tb_trt->Write();
1697 errorsAr_vs_tb_ba->Write();
1698 errorsAr_vs_tb_bc->Write();
1699 errorsAr_vs_tb_ea->Write();
1700 errorsAr_vs_tb_ec->Write();
1701
1702 residualAr_vs_tb_trt->Write();
1703 residualAr_vs_tb_ba->Write();
1704 residualAr_vs_tb_bc->Write();
1705 residualAr_vs_tb_ea->Write();
1706 residualAr_vs_tb_ec->Write();
1707
1708 trackerrorsAr_trtP->Write();
1709 trackerrorsAr_baP->Write();
1710 trackerrorsAr_bcP->Write();
1711 trackerrorsAr_eaP->Write();
1712 trackerrorsAr_ecP->Write();
1713
1714 trackerrorsAr_vs_tb_trt->Write();
1715 trackerrorsAr_vs_tb_ba->Write();
1716 trackerrorsAr_vs_tb_bc->Write();
1717 trackerrorsAr_vs_tb_ea->Write();
1718 trackerrorsAr_vs_tb_ec->Write();
1719
1720 // Save HT correction and ToT
1721 corrdirAr->cd();
1722
1723 tresAr_vs_ToT_trt->Write();
1724 tresAr_vs_ToT_ba->Write();
1725 tresAr_vs_ToT_bc->Write();
1726 tresAr_vs_ToT_ea->Write();
1727 tresAr_vs_ToT_ec->Write();
1728
1729 resAr_vs_ToT_trt->Write();
1730 resAr_vs_ToT_ba->Write();
1731 resAr_vs_ToT_bc->Write();
1732 resAr_vs_ToT_ea->Write();
1733 resAr_vs_ToT_ec->Write();
1734
1735 tresAr_vs_HT_trt->Write();
1736 tresAr_vs_HT_ba->Write();
1737 tresAr_vs_HT_bc->Write();
1738 tresAr_vs_HT_ea->Write();
1739 tresAr_vs_HT_ec->Write();
1740
1741 resAr_vs_HT_trt->Write();
1742 resAr_vs_HT_ba->Write();
1743 resAr_vs_HT_bc->Write();
1744 resAr_vs_HT_ea->Write();
1745 resAr_vs_HT_ec->Write();
1746
1747 tresAr_vs_SinOverP_trt->Write();
1748 tresAr_vs_SinOverP_ba->Write();
1749 tresAr_vs_SinOverP_bc->Write();
1750 tresAr_vs_SinOverP_ea->Write();
1751 tresAr_vs_SinOverP_ec->Write();
1752
1753 resAr_vs_SinOverP_trt->Write();
1754 resAr_vs_SinOverP_ba->Write();
1755 resAr_vs_SinOverP_bc->Write();
1756 resAr_vs_SinOverP_ea->Write();
1757 resAr_vs_SinOverP_ec->Write();
1758
1759 tresAr_vs_CosOverP_trt->Write();
1760 tresAr_vs_CosOverP_ba->Write();
1761 tresAr_vs_CosOverP_bc->Write();
1762 tresAr_vs_CosOverP_ea->Write();
1763 tresAr_vs_CosOverP_ec->Write();
1764
1765 resAr_vs_CosOverP_trt->Write();
1766 resAr_vs_CosOverP_ba->Write();
1767 resAr_vs_CosOverP_bc->Write();
1768 resAr_vs_CosOverP_ea->Write();
1769 resAr_vs_CosOverP_ec->Write();
1770 } // Close argon part
1771
1772 ttfile->Write();
1773 ttfile->Close();
1774 }
1775
1776 //*********************************************************************
1777 // SCANNING BINARY FILE
1778 }
1779 else
1780 {
1781 //*********************************************************************
1782
1783 cout << "SCANNING BINARY HISTOGRAM FILE " << argv[ifiles + 2] << endl;
1784
1785 nhistsadd = 0;
1786 nhists = 0;
1787
1788 ifstream *ifile = new ifstream(argv[ifiles + 2], ios::in | ios::binary);
1789
1790 // histogram loop
1791 while (true)
1792 {
1793
1794 // reading histogram
1795 //dangerous cast, reproduces previous c-style cast
1796 ifile->read(reinterpret_cast<char *>(&npop), sizeof(int));
1797 if (ifile->eof())
1798 break;
1799 int *chist = new int[(2 * npop)];
1800 if (npop > 0){
1801 //dangerous cast, reproduces previous c-style cast
1802 ifile->read(reinterpret_cast<char *>(chist), sizeof(int) * 2 * npop);
1803 }
1804 //dangerous cast, reproduces previous c-style cast
1805 ifile->read(reinterpret_cast<char *>(&sid), sizeof(int));
1806
1807 if (histmap.find(sid) == histmap.end())
1808 { // create histogram if seen for the first time
1809 // cppcheck-suppress stlFindInsert
1810 histmap[sid] = new CalHist();
1811 nhistsadd++;
1812 }
1813
1814 for (int ipop = 0; ipop < 2 * npop; ipop = ipop + 2)
1815 {
1816 histmap[sid]->IncreaseBin(static_cast<short>(chist[ipop]), static_cast<unsigned short>(chist[ipop + 1])); // increase the bins
1817 if (chist[ipop] < 100)
1818 {
1819 ntres += chist[ipop + 1];
1820 } // increase hit counters
1821 else if (chist[ipop] >= 100 && chist[ipop] < 200)
1822 nres += chist[ipop + 1];
1823 else
1824 nrt += chist[ipop + 1];
1825 }
1826
1827 nhists++;
1828
1829 if (nhists % 10000 == 0)
1830 {
1831 process_mem_usage(vm, rss);
1832 printf("%i HISTOGRAMS READ! %i HISTOGRAMS ADDED! %i BINS ADDED! MAXVALUE: %i VM: %0f RSS: %0f\n", nhists, nhistsadd, nbinsadd, maxvalue, vm, rss);
1833 }
1834
1835 std::destroy_at(chist);
1836
1837 if (nhists % 1000000000 == 0)
1838 {
1839 process_mem_usage(vm, rss);
1840 printf("%i HISTOGRAMS READ! %i HISTOGRAMS ADDED! %i BINS ADDED! MAXVALUE: %i VM: %0f RSS: %0f\n", nhists, nhistsadd, nbinsadd, maxvalue, vm, rss);
1841 dump_hists(&histmap, ntbins, nrbins, ntres, argv[1], 0);
1842 for (map<int, CalHist *>::iterator it = histmap.begin(); it != histmap.end(); ++it)
1843 {
1844 delete it->second;
1845 }
1846 histmap.clear();
1847 }
1848 }
1849
1850 delete ifile;
1851
1852 cout << nhists << " HISTOGRAMS SACNNED" << endl;
1853 cout << nhistsadd << " HISTOGRAMS ADDED" << endl;
1854
1855 } // end of if (root or binary histograms)
1856 } // end of file loop
1857
1858 cout << endl;
1859 cout << "TOT HISTOGRAMS............. " << nhists << endl;
1860 if (nhits!=0){
1861 cout << "TOT HITS IN TIMERES HIST .. " << ntres << " (" << 100 * float(ntres) / float(nhits) << "%)" << endl;
1862 cout << "TOT HITS IN RRES HIST ..... " << nres << " (" << 100 * float(nres) / float(nhits) << "%)" << endl;
1863 cout << "TOT HITS IN RT HIST ....... " << nrt << " (" << 100 * float(nrt) / float(nhits) << "%)" << endl;
1864 } else {
1865 cout << "Number of hits was zero !!\n";
1866 }
1867 cout << endl;
1868
1869 ofilestat << "TRESHITS " << ntres << endl;
1870 ofilestat.close();
1871 dump_hists(&histmap, ntbins, nrbins, ntres, argv[1], 0);
1872 dump_tracktuple(&trackmap);
1873
1874 return 0;
1875}
Scalar phi() const
phi method
Scalar theta() const
theta method
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
static Double_t t0
map< string, calibpars > readoldpars()
void process_mem_usage(double &vm_usage, double &resident_set)
int Simple1dHist(float min, float max, int nbins, float value)
int dump_hists(map< int, CalHist * > *histmap, int ntbins, int nrbins, int, char *fname, int fileno)
int dump_tracktuple(map< string, trackdata > *trackmap)
int Simple2dHist(float minx, float maxx, int nbinsx, float miny, float maxy, int nbinsy, float valuex, float valuey)
int r
Definition globals.cxx:22
time(flags, cells_name, *args, **kw)
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
int run(int argc, char *argv[])

◆ process_mem_usage()

void process_mem_usage ( double & vm_usage,
double & resident_set )

Definition at line 42 of file TRTCalib_bhadd.cxx.

43{
44 using std::ifstream;
45 using std::ios_base;
46 using std::string;
47
48 vm_usage = 0.0;
49 resident_set = 0.0;
50
51 // 'file' stat seems to give the most reliable results
52 //
53 ifstream stat_stream("/proc/self/stat", ios_base::in);
54
55 // dummy vars for leading entries in stat that we don't care about
56 //
57 string pid, comm, state, ppid, pgrp, session, tty_nr;
58 string tpgid, flags, minflt, cminflt, majflt, cmajflt;
59 string utime, stime, cutime, cstime, priority, nice;
60 string O, itrealvalue, starttime;
61
62 // the two fields we want
63 //
64 unsigned long vsize;
65 long rss;
66
67 stat_stream >> pid >> comm >> state >> ppid >> pgrp >> session >> tty_nr >> tpgid >> flags >> minflt >> cminflt >> majflt >> cmajflt >> utime >> stime >> cutime >> cstime >> priority >> nice >> O >> itrealvalue >> starttime >> vsize >> rss; // don't care about the rest
68
69 long page_size_kb = sysconf(_SC_PAGE_SIZE) / 1024; // in case x86-64 is configured to use 2MB pages
70 vm_usage = vsize / 1024.0;
71 resident_set = rss * page_size_kb;
72}
std::string stime()
return the current data and time

◆ readoldpars()

map< string, calibpars > readoldpars ( )

Definition at line 243 of file TRTCalib_bhadd.cxx.

244{
245 string line;
246 ifstream oldconstfile("precision_constants.txt");
247 int det, lay, mod, stl, stw;
248 float t0, dt0;
249 map<string, calibpars> oldparsmap;
250 char key[100];
251
252 if (oldconstfile.is_open())
253 {
254 while (!oldconstfile.eof())
255 {
256 getline(oldconstfile, line);
257 if (line.size() > 25 && line.size() < 55)
258 {
259 sscanf(line.data(), "%d %d %d %d %d : %e %e", &det, &lay, &mod, &stl, &stw, &t0, &dt0);
260 if (det >= -3 && lay > -1 && mod > -1 && stl > -1 && stw > -1)
261 {
262 sprintf(key, "_%i_%i_%i_%i_%i", det, lay, mod, stl, stw);
263 oldparsmap[string(key)].t0 = t0;
264 oldparsmap[string(key)].dt0 = dt0;
265 }
266 }
267 }
268 oldconstfile.close();
269 }
270 else{
271 // Output to std::cout since for online calibration this file is not used.
272 std::cout << "Unable to open precision t0 file!\n";
273 }
274
275
276 return oldparsmap;
277}

◆ Simple1dHist()

int Simple1dHist ( float min,
float max,
int nbins,
float value )

Definition at line 25 of file TRTCalib_bhadd.cxx.

26{
27 if ((value < min) or (value > max))
28 return -1;
29 int binno = (int)(nbins * ((value - min) / (max - min)));
30 return binno;
31}
#define min(a, b)
Definition cfImp.cxx:40
#define max(a, b)
Definition cfImp.cxx:41

◆ Simple2dHist()

int Simple2dHist ( float minx,
float maxx,
int nbinsx,
float miny,
float maxy,
int nbinsy,
float valuex,
float valuey )

Definition at line 33 of file TRTCalib_bhadd.cxx.

34{
35 if (valuex < minx or valuex > maxx or valuey < miny or valuey > maxy)
36 return -1;
37 int binnox = (int)(nbinsx * ((valuex - minx) / (maxx - minx)));
38 int binnoy = (int)(nbinsy * ((valuey - miny) / (maxy - miny)));
39 return binnoy * nbinsx + binnox;
40}