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 <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 293 of file TRTCalib_bhadd.cxx.

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

◆ dump_tracktuple()

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

Definition at line 278 of file TRTCalib_bhadd.cxx.

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

◆ main()

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

Definition at line 320 of file TRTCalib_bhadd.cxx.

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

◆ process_mem_usage()

void process_mem_usage ( double & vm_usage,
double & resident_set )

Definition at line 41 of file TRTCalib_bhadd.cxx.

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

◆ readoldpars()

map< string, calibpars > readoldpars ( )

Definition at line 242 of file TRTCalib_bhadd.cxx.

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

◆ Simple1dHist()

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

Definition at line 24 of file TRTCalib_bhadd.cxx.

25{
26 if ((value < min) or (value > max))
27 return -1;
28 int binno = (int)(nbins * ((value - min) / (max - min)));
29 return binno;
30}
#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 32 of file TRTCalib_bhadd.cxx.

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