7 #define _USE_MATH_DEFINES
11 #include <marlin/Global.h>
12 #include <marlin/Processor.h>
13 #include <marlin/AIDAProcessor.h>
16 #include <EVENT/LCCollection.h>
17 #include <EVENT/Track.h>
18 #include <EVENT/TrackerHit.h>
19 #include <IMPL/TrackImpl.h>
21 #include <DD4hep/Detector.h>
22 #include <DDRec/DetectorData.h>
23 #include <DD4hep/DD4hepUnits.h>
25 #include <UTIL/BitField64.h>
26 #include "UTIL/LCTrackerConf.h"
29 #include "LCGeometryTypes.h"
30 #include "SimpleHelix.h"
31 #include "HelixClass.h"
42 : Processor(
"Compute_dEdxProcessor2021") {
45 _description =
"dE/dx calculation using truncation method, angular correction applied after smearing" ;
47 registerInputCollection(LCIO::TRACK,
49 "LDC track collection name",
51 std::string(
"MarlinTrkTracks"));
53 registerProcessorParameter(
"Write_dEdx",
54 "If set, the calculated dEdx value will be attached to its corresponding track (default: true).",
58 registerProcessorParameter(
"EnergyLossErrorTPC",
59 "Fractional error of dEdx in the TPC",
63 registerProcessorParameter(
"LowerTruncationFraction",
64 "Lower truncation fraction, default: 0.08 (8%)",
68 registerProcessorParameter(
"UpperTruncationFraction",
69 "Upper truncation fraction, default: 0.3 (30%)",
73 registerProcessorParameter(
"isSmearing",
74 "Flag for dEdx Smearing",
78 registerProcessorParameter(
"smearingFactor",
79 "Smearing factor for dEdx smearing",
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",
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",
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).",
99 registerProcessorParameter(
"StrategyCompHistWeight",
100 "If set, the Bethe-Bloch histograms will have a sqrt(nHits) weighting (default: false).",
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"));
109 registerProcessorParameter(
"doAngularCorr_dEdx",
110 "If set, the calculated dEdx value is corrected as function of a pol3(lambda)",
114 std::vector<float> _newpar = { 0.948185,
118 registerProcessorParameter(
"AngularCorrectionParameters",
119 "parameter for new angular correction dedx= uncorrected_dedx / f, with f= pol3(lambda)",
128 streamlog_out(DEBUG) <<
" init called " << std::endl ;
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>() ;
141 lcdd.field().magneticField( { 0., 0., 0. } , bfieldV ) ;
142 _bField = bfieldV[2]/dd4hep::tesla ;
144 engine =
new std::default_random_engine();
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;
153 const int nBinsX=1000, nBinsY=300;
154 double minBinX=-1, maxBinX=2;
155 double minBinY= 0, maxBinY=1
e-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);
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);
173 int nTrkCand =
_LDCCol->getNumberOfElements();
175 for (
int iTRK=0;iTRK<nTrkCand;++iTRK) {
177 TrackImpl * trkCand = (TrackImpl*)
_LDCCol->getElementAt( iTRK );
178 TrackerHitVec trkHits = trkCand->getTrackerHits();
186 TrackerHitVec tpcHitVec;
187 for(
unsigned int sdet=0;sdet<trkHits.size(); sdet++){
188 if(trkHits[sdet]->getType() == 0){
190 float x=trkHits[sdet]->getPosition()[0];
191 float y=trkHits[sdet]->getPosition()[1];
193 else streamlog_out(DEBUG) <<
"hit not in TPC: " << x <<
" " << y <<
" " << trkHits[sdet]->getPosition()[2] <<
" " << sqrt(x*x+y*y) << std::endl;
202 float trklambda = trkCand->getTanLambda();
210 trkCand->setdEdx(dedx);
211 trkCand->setdEdxError(unum);
223 TCanvas *can =
new TCanvas;
224 TImage *img = TImage::Create();
225 gStyle->SetPalette(kBird);
258 int tH=hitVec.size();
261 const float* refPoint = trk->getReferencePoint();
262 LCVector3D refPointV3D = LCVector3D(refPoint[0],refPoint[1],refPoint[2]);
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);
268 std::vector<double> dEdxVec1, dEdxVec2, dEdxVec3;
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;
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];
286 TrackerHit * hit0 = hitVec[iH-1];
287 hit0X = hit0->getPosition()[0];
288 hit0Y = hit0->getPosition()[1];
289 hit0Z = hit0->getPosition()[2];
298 dxy=sqrt(dx*dx+dy*dy);
299 r=1.0/trk->getOmega();
300 theta=2.0*asin(0.5*dxy/r);
303 length=sqrt(l*l+dz*dz);
307 dEdxVec1.push_back(hit->getEDep()/length);
312 double path1 = trkHelix.getPathAt(LCVector3D(hitX ,hitY ,hitZ ));
313 double path2 = trkHelix.getPathAt(LCVector3D(hit0X,hit0Y,hit0Z));
314 path = fabs(path2-path1);
316 dEdxVec2.push_back(hit->getEDep()/path);
325 double hitRadius = sqrt(hitX*hitX + hitY*hitY);
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);
338 bool rejected =
false;
339 if (fabs(time1)<1
e+10 && fabs(time2)<1
e+10) {
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;
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;
343 else rejected =
true;
347 float* cross_1 = cross1;
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];
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));
364 dEdxVec3.push_back(hit->getEDep()/over);
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;
371 streamlog_out(DEBUG4) <<
"1: " << totlength <<
" 2: " << totpath <<
" of " << cpath <<
" 3: " << totover <<
" of " << cover << std::endl;
375 double normdedx=0, totway=0;
382 double dedx1=0, dedx2=0, dedx3=0;
383 int nTruncate1=0, nTruncate2=0, nTruncate3=0;
385 std::sort(dEdxVec1.begin(),dEdxVec1.end());
386 std::sort(dEdxVec2.begin(),dEdxVec2.end());
387 std::sort(dEdxVec3.begin(),dEdxVec3.end());
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;
399 double normdedx1 = dedx1;
400 double normdedx2 = dedx2;
401 double normdedx3 = dedx3;
409 double lmom = sqrt(pow(trkHelixC.getMomentum()[0],2) + pow(trkHelixC.getMomentum()[1],2) + pow(trkHelixC.getMomentum()[2],2));
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;
432 std::vector<double> dEdxVec;
435 case 1: dEdxVec = dEdxVec1; totway = totlength;
break;
436 case 2: dEdxVec = dEdxVec2; totway = totpath;
break;
437 case 3: dEdxVec = dEdxVec3; totway = totover;
break;
440 std::sort(dEdxVec.begin(),dEdxVec.end());
442 dedx = dedx/(float)(nTruncate);
448 std::pair<double,double> ret;
449 if(std::isnan(normdedx))
453 streamlog_out(DEBUG4) <<
"NAN!" << std::endl;
458 ret.second=normdedx*std::pow(0.001*totway,
_errexp[0])*std::pow(nTruncate,
_errexp[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) );
virtual void processRunHeader(LCRunHeader *run)
std::string _LDCTrackCollection
virtual void check(LCEvent *evt)
Compute dE/dx Processor2021 This processor calculates the dE/dx for every track.
std::uniform_real_distribution dist
bool _StratCompHistWeight
double get_Corrected_dEdx(float dedx, float trklambda)
Compute_dEdxProcessor2021 aCompute_dEdxProcessor2021
std::default_random_engine * engine
std::vector< float > _errexp
float _energyLossErrorTPC
double getSmearing(double dEdx)
virtual void processEvent(LCEvent *evt)
std::vector< float > _par
std::string _StratCompHistFiles
Compute_dEdxProcessor2021()
std::pair< double, double > CalculateEnergyLoss(TrackerHitVec &hitVec, Track *trk)