All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
Compute_dEdxProcessor2021.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <sstream>
4 #include <fstream>
5 #include <utility>
6 
7 #define _USE_MATH_DEFINES
8 #include <cmath>
9 #include <math.h>
10 
11 #include <marlin/Global.h>
12 #include <marlin/Processor.h>
13 #include <marlin/AIDAProcessor.h>
14 
15 #include <lcio.h>
16 #include <EVENT/LCCollection.h>
17 #include <EVENT/Track.h>
18 #include <EVENT/TrackerHit.h>
19 #include <IMPL/TrackImpl.h>
20 
21 #include <DD4hep/Detector.h>
22 #include <DDRec/DetectorData.h>
23 #include <DD4hep/DD4hepUnits.h>
24 
25 #include <UTIL/BitField64.h>
26 #include "UTIL/LCTrackerConf.h"
27 
28 // MarlinUtil
29 #include "LCGeometryTypes.h"
30 #include "SimpleHelix.h"
31 #include "HelixClass.h"
32 
33 #include "TCanvas.h"
34 #include "TImage.h"
35 #include "TStyle.h"
36 
38 
40 
42  : Processor("Compute_dEdxProcessor2021") {
43 
44  // Processor description
45  _description = "dE/dx calculation using truncation method, angular correction applied after smearing" ;
46 
47  registerInputCollection(LCIO::TRACK,
48  "LDCTrackCollection",
49  "LDC track collection name",
51  std::string("MarlinTrkTracks"));
52 
53  registerProcessorParameter( "Write_dEdx",
54  "If set, the calculated dEdx value will be attached to its corresponding track (default: true).",
55  _writedEdx,
56  bool(true));
57 
58  registerProcessorParameter( "EnergyLossErrorTPC",
59  "Fractional error of dEdx in the TPC",
61  float(0.054));
62 
63  registerProcessorParameter( "LowerTruncationFraction",
64  "Lower truncation fraction, default: 0.08 (8%)",
66  float(0.08));
67 
68  registerProcessorParameter( "UpperTruncationFraction",
69  "Upper truncation fraction, default: 0.3 (30%)",
71  float(0.3));
72 
73  registerProcessorParameter( "isSmearing",
74  "Flag for dEdx Smearing",
76  bool(0));
77 
78  registerProcessorParameter( "smearingFactor",
79  "Smearing factor for dEdx smearing",
81  float(0.035));
82 
83  std::vector<float> errexp = {-0.34, -0.45};
84  registerProcessorParameter( "dEdxErrorScalingExponents",
85  "scaling exponents of the dE/dx error for path length and number of hits, respectively",
86  _errexp,
87  errexp);
88 
89  registerProcessorParameter( "dxStrategy",
90  "ID of the track length calculation: 1 - hit-to-hit distance; 2 - hit-to-hit path length of projected hits; 3 - path over hit row",
92  int(1));
93 
94  registerProcessorParameter( "StrategyCompHist",
95  "If set, an AIDA root file will be generated with Bethe-Bloch histograms (dE/dx over momentum) for each dx strategy (default: false).",
97  bool(false));
98 
99  registerProcessorParameter( "StrategyCompHistWeight",
100  "If set, the Bethe-Bloch histograms will have a sqrt(nHits) weighting (default: false).",
102  bool(false));
103 
104  registerProcessorParameter( "StrategyCompHistFiles",
105  "File names of the generated dx strategy comparison histograms (if chosen). The respective strategy number and '.png' is added.",
107  std::string("dEdx_Histo_Strategy"));
108 
109  registerProcessorParameter( "doAngularCorr_dEdx",
110  "If set, the calculated dEdx value is corrected as function of a pol3(lambda)",
112  bool(false));
113 
114  std::vector<float> _newpar = { 0.948185,
115  0.000520502,
116  3.33231e-5,
117  -1.51052e-7};
118  registerProcessorParameter( "AngularCorrectionParameters",
119  "parameter for new angular correction dedx= uncorrected_dedx / f, with f= pol3(lambda)",
120  _par,
121  _newpar);
122 
123 
124 
125 }
126 
128  streamlog_out(DEBUG) << " init called " << std::endl ;
129 
130  // it's usually a good idea to
131  printParameters();
132 
133  //get TPC radii, pad height and B field
134  dd4hep::Detector& lcdd = dd4hep::Detector::getInstance();
135  dd4hep::DetElement tpcDE = lcdd.detector("TPC") ;
136  const dd4hep::rec::FixedPadSizeTPCData* tpc = tpcDE.extension<dd4hep::rec::FixedPadSizeTPCData>() ;
137  _TPC_inner = tpc->rMinReadout/dd4hep::mm;
138  _TPC_outer = tpc->rMaxReadout/dd4hep::mm;
139  _TPC_padHeight = tpc->padHeight/dd4hep::mm;
140  double bfieldV[3] ;
141  lcdd.field().magneticField( { 0., 0., 0. } , bfieldV ) ;
142  _bField = bfieldV[2]/dd4hep::tesla ;
143 
144  engine = new std::default_random_engine();
145 
146  streamlog_out(MESSAGE) << "TPC inner radius = " << _TPC_inner << " mm TPC outer radius = " << _TPC_outer << " mm"<< std::endl;
147  streamlog_out(MESSAGE) << "TPC pad height = " << _TPC_padHeight << " mm B field = " << _bField << " T" <<std::endl;
148  streamlog_out(DEBUG) << "dd4hep::mm = " << dd4hep::mm << " dd4hep::cm = " << dd4hep::cm << std::endl;
149 
150  if (_StratCompHist)
151  {
152  // prepare the histogram bins - this could also be implemented as processor parameters
153  const int nBinsX=1000, nBinsY=300;
154  double minBinX=-1, maxBinX=2; // log10(min/max momentum / GeV)
155  double minBinY= 0, maxBinY=1e-6;
156  double histbinsX[nBinsX+1], histbinsY[nBinsY+1];
157  for (int i=0; i<nBinsX+1; i++) histbinsX[i] = pow( 10, minBinX + (maxBinX-minBinX)*i/(float)(nBinsX) );
158  for (int i=0; i<nBinsY+1; i++) histbinsY[i] = minBinY + (maxBinY-minBinY)*i/(float)(nBinsY);
159 
160  // create the histograms
161  AIDAProcessor::tree(this);
162  _BBHist_Strategy1 = new TH2D("BBHist_Strategy1","Bethe-Bloch curve for dx strategy 1: hit-to-hit distance",nBinsX,histbinsX,nBinsY,histbinsY);
163  _BBHist_Strategy2 = new TH2D("BBHist_Strategy2","Bethe-Bloch curve for dx strategy 2: hit-to-hit path length of projected hits",nBinsX,histbinsX,nBinsY,histbinsY);
164  _BBHist_Strategy3 = new TH2D("BBHist_Strategy3","Bethe-Bloch curve for dx strategy 3: path over hit row",nBinsX,histbinsX,nBinsY,histbinsY);
165  }
166 }
167 
169 }
170 
172  _LDCCol = evt->getCollection( _LDCTrackCollection ) ;
173  int nTrkCand = _LDCCol->getNumberOfElements();
174 
175  for (int iTRK=0;iTRK<nTrkCand;++iTRK) {
176 
177  TrackImpl * trkCand = (TrackImpl*) _LDCCol->getElementAt( iTRK );
178  TrackerHitVec trkHits = trkCand->getTrackerHits();
179 
180  //calculate dEdx here when the subtrack is the TPC track!
181  //silicon dE/dx can also be calculated. is it necessary??
182  float dedx=0.0;
183  float unum=0.0;
184 
185  //get tpc hits only
186  TrackerHitVec tpcHitVec;
187  for(unsigned int sdet=0;sdet<trkHits.size(); sdet++){
188  if(trkHits[sdet]->getType() == 0){ //it is very suspicious condition...
189  //check whether this hit is in TPC
190  float x=trkHits[sdet]->getPosition()[0];
191  float y=trkHits[sdet]->getPosition()[1];
192  if(sqrt(x*x+y*y)>=_TPC_inner && sqrt(x*x+y*y)<=_TPC_outer) tpcHitVec.push_back(trkHits[sdet]); // add hit
193  else streamlog_out(DEBUG) << "hit not in TPC: " << x << " " << y << " " << trkHits[sdet]->getPosition()[2] << " " << sqrt(x*x+y*y) << std::endl;
194  }
195  }
196 
197  std::pair<double,double> dummy = CalculateEnergyLoss(tpcHitVec, trkCand);
198  dedx = dummy.first;
199  unum = dummy.second*_energyLossErrorTPC;
200 
201  if(_angularcorrdEdx==true) {
202  float trklambda = trkCand->getTanLambda();
203  dedx = get_Corrected_dEdx(dedx,trklambda);
204 
205  }
206 
207  //fill values
208  if (_writedEdx)
209  {
210  trkCand->setdEdx(dedx);
211  trkCand->setdEdxError(unum);
212  }
213  }
214 }
215 
217 }
218 
220 
221  if (_StratCompHist)
222  {
223  TCanvas *can = new TCanvas;
224  TImage *img = TImage::Create();
225  gStyle->SetPalette(kBird);
226  can->SetGrid();
227  can->SetLogx();
228  _BBHist_Strategy1->SetXTitle("momentum / GeV");
229  _BBHist_Strategy1->SetYTitle("dE/dx / (GeV/mm)");
230  _BBHist_Strategy2->SetXTitle("momentum / GeV");
231  _BBHist_Strategy2->SetYTitle("dE/dx / (GeV/mm)");
232  _BBHist_Strategy3->SetXTitle("momentum / GeV");
233  _BBHist_Strategy3->SetYTitle("dE/dx / (GeV/mm)");
234 
235  std::string StratCompHistFiles2 = _StratCompHistFiles;
236  _BBHist_Strategy1->Draw("colz");
237  img->FromPad(can);
238  _StratCompHistFiles.append("1.png");
239  img->WriteImage(_StratCompHistFiles.c_str());
240  _BBHist_Strategy2->Draw("colz");
241  img->FromPad(can);
242  _StratCompHistFiles.replace(_StratCompHistFiles.end()-5,_StratCompHistFiles.end()-4,"2");
243  img->WriteImage(_StratCompHistFiles.c_str());
244  _BBHist_Strategy3->Draw("colz");
245  img->FromPad(can);
246  _StratCompHistFiles.replace(_StratCompHistFiles.end()-5,_StratCompHistFiles.end()-4,"3");
247  img->WriteImage(_StratCompHistFiles.c_str());
248 
249  delete can;
250  delete img;
251  }
252 
253  delete engine;
254 }
255 
256 
257 std::pair<double,double> Compute_dEdxProcessor2021::CalculateEnergyLoss(TrackerHitVec& hitVec, Track* trk){
258  int tH=hitVec.size();
259 
260  // fill the track into a SimpleHelix and also into a HelixClass to make use of their specific geometric utilities
261  const float* refPoint = trk->getReferencePoint();
262  LCVector3D refPointV3D = LCVector3D(refPoint[0],refPoint[1],refPoint[2]); // MarlinUtil::LCVector3D is CLHEP::Hep3Vector
263  SimpleHelix trkHelix = SimpleHelix(trk->getD0(),trk->getPhi(),trk->getOmega(),trk->getZ0(),trk->getTanLambda(),refPointV3D);
264  HelixClass trkHelixC;
265  trkHelixC.Initialize_Canonical(trk->getPhi(),trk->getD0(),trk->getZ0(),trk->getOmega(),trk->getTanLambda(),_bField);
266 
267  // set up some variables
268  std::vector<double> dEdxVec1, dEdxVec2, dEdxVec3; // one vector per dx strategy
269  double dx=0, dy=0, dz=0, length=0, totlength=0;
270  double dxy=0, r=0, theta=0, l=0;
271  double path=0, totpath=0, cpath=0;
272  double over=0, totover=0, cover=0;
273  double hitX =0, hitY =0, hitZ =0;
274  double hit0X=0, hit0Y=0, hit0Z=0;
275 
276  for(int iH=0; iH<tH; iH++){
277  TrackerHit * hit = hitVec[iH];
278  hitX = hit ->getPosition()[0];
279  hitY = hit ->getPosition()[1];
280  hitZ = hit ->getPosition()[2];
281 
282  if (iH>0)
283  {
284  if (_dxStrategy==1 || _dxStrategy==2 || _StratCompHist)
285  {
286  TrackerHit * hit0 = hitVec[iH-1];
287  hit0X = hit0->getPosition()[0];
288  hit0Y = hit0->getPosition()[1];
289  hit0Z = hit0->getPosition()[2];
290 
291  if (_dxStrategy==1 || _StratCompHist)
292  { // strategy 1: dx = distance between 2 hits
293  dx=hitX-hit0X;
294  dy=hitY-hit0Y;
295  dz=hitZ-hit0Z;
296 
297  //cal. length of curvature
298  dxy=sqrt(dx*dx+dy*dy);
299  r=1.0/trk->getOmega();
300  theta=2.0*asin(0.5*dxy/r);
301  l=r*theta;
302  //total helix length
303  length=sqrt(l*l+dz*dz);
304 
305  //sum of flight length
306  totlength += length;
307  dEdxVec1.push_back(hit->getEDep()/length);
308  }
309 
310  if (_dxStrategy==2 || _StratCompHist)
311  { // strategy 2: dx = path difference on the track via SimpleHelix
312  double path1 = trkHelix.getPathAt(LCVector3D(hitX ,hitY ,hitZ ));
313  double path2 = trkHelix.getPathAt(LCVector3D(hit0X,hit0Y,hit0Z));
314  path = fabs(path2-path1);
315  totpath += path;
316  dEdxVec2.push_back(hit->getEDep()/path);
317  cpath++;
318  }
319  }
320  }
321 
322  if (_dxStrategy==3 || _StratCompHist)
323  { // strategy 3: dx = path of track over the hit's row via SimpleHelix and HelixClass
324  // since this strategy does not require a distance between two hits, it can be used for every single hit
325  double hitRadius = sqrt(hitX*hitX + hitY*hitY);
326 
327  // create cylinders of the row's inner and outer edge, check if the track crosses both, then get the path difference of these crossing points
328  double cyl1_radius = hitRadius - _TPC_padHeight*.5; // assuming the radial hit position is always in the center of the row
329  double cyl2_radius = hitRadius + _TPC_padHeight*.5;
330  float trkRefP[3];
331  trkRefP[0] = (trkHelixC.getReferencePoint())[0];
332  trkRefP[1] = (trkHelixC.getReferencePoint())[1];
333  trkRefP[2] = (trkHelixC.getReferencePoint())[2];
334  float cross1[6], cross2[6];
335  float time1 = trkHelixC.getPointOnCircle(cyl1_radius, trkRefP, cross1);
336  float time2 = trkHelixC.getPointOnCircle(cyl2_radius, trkRefP, cross2);
337 
338  bool rejected = false;
339  if (fabs(time1)<1e+10 && fabs(time2)<1e+10) { // getPointOnCircle returns -1e+20 when the helix doesn't cross the cylinder
340  if (sqrt(pow(cross1[0]-cross1[3],2) + pow(cross1[1]-cross1[4],2) + pow(cross1[2]-cross1[5],2)) < 2*_TPC_padHeight) rejected = true; // curler protection
341  if (sqrt(pow(cross2[0]-cross2[3],2) + pow(cross2[1]-cross2[4],2) + pow(cross2[2]-cross2[5],2)) < 2*_TPC_padHeight) rejected = true;
342  }
343  else rejected = true;
344 
345  if (!rejected)
346  {
347  float* cross_1 = cross1; // check which crossing point is the correct one by comparing to the hit position
348  double dist1 = sqrt(pow(cross1[0]-hitX,2) + pow(cross1[1]-hitY,2) + pow(cross1[2]-hitZ,2));
349  double dist2 = sqrt(pow(cross1[3]-hitX,2) + pow(cross1[4]-hitY,2) + pow(cross1[5]-hitZ,2));
350  if (dist2 < dist1) cross_1 = &cross1[3];
351  float* cross_2 = cross2;
352  dist1 = sqrt(pow(cross2[0]-hitX,2) + pow(cross2[1]-hitY,2) + pow(cross2[2]-hitZ,2));
353  dist2 = sqrt(pow(cross2[3]-hitX,2) + pow(cross2[4]-hitY,2) + pow(cross2[5]-hitZ,2));
354  if (dist2 < dist1) cross_2 = &cross2[3];
355 
356  // again correct for curvature
357  dxy = sqrt(pow(cross_1[0]-cross_2[0],2) + pow(cross_1[1]-cross_2[1],2));
358  r = 1./trk->getOmega();
359  l = r*2.*asin(0.5*dxy/r);
360  over= sqrt(l*l + pow(cross_1[2]-cross_2[2],2));
361 
362  totover += over;
363  cover++;
364  dEdxVec3.push_back(hit->getEDep()/over);
365  }
366  }
367 
368  streamlog_out(DEBUG3) << "1: " << length << " 2: " << path << " 3: " << over << " EDep: " << hit->getEDep() << std::endl;
369  if (length/over>1.5 || over/length>1.5) streamlog_out(DEBUG2) << " " << hitX << " " << hitY << " " << hitZ << " " << sqrt(hitX*hitX + hitY*hitY) << std::endl;
370  }
371  streamlog_out(DEBUG4) << "1: " << totlength << " 2: " << totpath << " of " << cpath << " 3: " << totover << " of " << cover << std::endl;
372 
373 
374  int nTruncate=0; // number of remaining track hits after truncation
375  double normdedx=0, totway=0;
376  // double trkcos = sqrt(trk->getTanLambda()*trk->getTanLambda()/(1.0+trk->getTanLambda()*trk->getTanLambda()));
377 
378  // if comparison plots should be made, the dE/dx for all three strategies is calculated, filled in histograms,
379  // and the one of the chosen strategy is also transferred to be returned
380  if (_StratCompHist)
381  {
382  double dedx1=0, dedx2=0, dedx3=0;
383  int nTruncate1=0, nTruncate2=0, nTruncate3=0;
384 
385  std::sort(dEdxVec1.begin(),dEdxVec1.end());
386  std::sort(dEdxVec2.begin(),dEdxVec2.end());
387  std::sort(dEdxVec3.begin(),dEdxVec3.end());
388 
389  for (int idE=dEdxVec1.size()*_lowerTrunFrac; idE<dEdxVec1.size()*(1-_upperTrunFrac); idE++) {dedx1 += dEdxVec1[idE]; nTruncate1++;}
390  for (int idE=dEdxVec2.size()*_lowerTrunFrac; idE<dEdxVec2.size()*(1-_upperTrunFrac); idE++) {dedx2 += dEdxVec2[idE]; nTruncate2++;}
391  for (int idE=dEdxVec3.size()*_lowerTrunFrac; idE<dEdxVec3.size()*(1-_upperTrunFrac); idE++) {dedx3 += dEdxVec3[idE]; nTruncate3++;}
392 
393  //cal. truncated mean
394  dedx1 = dedx1/(float)(nTruncate1);
395  dedx2 = dedx2/(float)(nTruncate2);
396  dedx3 = dedx3/(float)(nTruncate3);
397  streamlog_out(DEBUG4) << " check dedx: " << dedx1 << " " << dedx2 << " " << dedx3 << std::endl;
398 
399  double normdedx1 = dedx1;//getNormalization(dedx1, (float)nTruncate1, trkcos);
400  double normdedx2 = dedx2;//getNormalization(dedx2, (float)nTruncate2, trkcos);
401  double normdedx3 = dedx3;//getNormalization(dedx3, (float)nTruncate3, trkcos);
402  if(_isSmearing)
403  {
404  normdedx1 += getSmearing(normdedx1);
405  normdedx2 += getSmearing(normdedx2);
406  normdedx3 += getSmearing(normdedx3);
407  }
408 
409  double lmom = sqrt(pow(trkHelixC.getMomentum()[0],2) + pow(trkHelixC.getMomentum()[1],2) + pow(trkHelixC.getMomentum()[2],2));
410 
411  if (_StratCompHistWeight) {
412  _BBHist_Strategy1->Fill(lmom,normdedx1,sqrt(nTruncate1));
413  _BBHist_Strategy2->Fill(lmom,normdedx2,sqrt(nTruncate2));
414  _BBHist_Strategy3->Fill(lmom,normdedx3,sqrt(nTruncate3));
415  }
416  else {
417  _BBHist_Strategy1->Fill(lmom,normdedx1);
418  _BBHist_Strategy2->Fill(lmom,normdedx2);
419  _BBHist_Strategy3->Fill(lmom,normdedx3);
420  }
421 
422  switch (_dxStrategy)
423  {
424  case 1: normdedx = normdedx1; nTruncate = nTruncate1; totway = totlength; break;
425  case 2: normdedx = normdedx2; nTruncate = nTruncate2; totway = totpath; break;
426  case 3: normdedx = normdedx3; nTruncate = nTruncate3; totway = totover; break;
427  }
428  }
429  else // to if (_StratCompHist)
430  {
431  double dedx=0;
432  std::vector<double> dEdxVec;
433  switch (_dxStrategy)
434  {
435  case 1: dEdxVec = dEdxVec1; totway = totlength; break;
436  case 2: dEdxVec = dEdxVec2; totway = totpath; break;
437  case 3: dEdxVec = dEdxVec3; totway = totover; break;
438  }
439 
440  std::sort(dEdxVec.begin(),dEdxVec.end());
441  for (int idE=dEdxVec.size()*_lowerTrunFrac; idE<dEdxVec.size()*(1-_upperTrunFrac); idE++) {dedx += dEdxVec[idE]; nTruncate++;}
442  dedx = dedx/(float)(nTruncate);
443  normdedx = dedx;//getNormalization(dedx, (float)nTruncate, trkcos);
444  if(_isSmearing) normdedx += getSmearing(normdedx);
445  }
446 
447 
448  std::pair<double,double> ret;
449  if(std::isnan(normdedx))
450  {
451  ret.first =0;
452  ret.second=0;
453  streamlog_out(DEBUG4) << "NAN!" << std::endl;
454  }
455  else
456  {
457  ret.first =normdedx;
458  ret.second=normdedx*std::pow(0.001*totway, _errexp[0])*std::pow(nTruncate, _errexp[1]); //todo: what is gas pressure in TPC?
459  }
460 
461  return ret;
462 }
463 
464 
465 //create smearing factor
467  double z=sqrt( -2.0*std::log(dist(*engine)) ) * std::sin( 2.0*M_PI*dist(*engine) ); // 3.141592
468 
469  //std::cout << dist(*engine) << " " << z << std::endl;
470 
471  return 0.0 + dEdx * _smearingFactor * z;
472 }
473 
474 double Compute_dEdxProcessor2021::get_Corrected_dEdx(float dedx, float trklambda){
475 
476  // // DEFAULT CORRECTION APPLIED FOR THE MC2020 production obtained using versions v02-02-01 and previous
477  // //cal. hit dep.
478  //double f1 = 1 + std::exp(-nHit/_ncorrpar);
479  // //cal. polar angle dep.
480  // // double c=1.0/sqrt(1.0-trkcos*trkcos);
481  // // double f2=1.0/(1.0-0.08887*std::log(c));
482 
483  // //cal. polar angle dep. 20160702
484  // //double c = std::acos(trkcos);
485  // //if(c>3.141592/2.0) c= 3.141592-c;
486  // //double f2 = 1.0/std::pow(c, 0.0703);
487 
488  // //change polar angle dep. 20170901
489  // double f2 = _acorrpar[0] / (_acorrpar[0] + _acorrpar[1] * trkcos * trkcos);
490  // return dedx/(f1*f2);
491 
492  // ****************************************
493  // NEW CORRECTION
494  //Parametrization by A. Irles, U. Einhaus 2021/04
495  //Using single particle mc2020 samples v02-02-01
496  //Analyzed with dEdxAnalyser, fitting histogram NormLambdaFullAll_1
497  float lambda = fabs(atan(trklambda)*180./M_PI);
498  double f3 = 1 / (_par[0] + _par[1] * lambda + _par[2] * pow(lambda,2) + _par[3] * pow(lambda,3) );
499 
500  return dedx*f3;
501 
502 }
503 
504 
static const float mm
virtual void processRunHeader(LCRunHeader *run)
virtual void check(LCEvent *evt)
Compute dE/dx Processor2021 This processor calculates the dE/dx for every track.
std::uniform_real_distribution dist
double get_Corrected_dEdx(float dedx, float trklambda)
Compute_dEdxProcessor2021 aCompute_dEdxProcessor2021
std::default_random_engine * engine
static const float e
virtual void processEvent(LCEvent *evt)
static const float cm
std::pair< double, double > CalculateEnergyLoss(TrackerHitVec &hitVec, Track *trk)