3 #include <EVENT/LCCollection.h>
4 #include <EVENT/MCParticle.h>
5 #include <IMPL/TrackImpl.h>
7 #include <EVENT/LCRelation.h>
8 #include <UTIL/LCRelationNavigator.h>
9 #include <IMPL/TrackStateImpl.h>
10 #include <IMPL/TrackerHitImpl.h>
13 #include "marlin/VerbosityLevels.h"
15 #include <DD4hep/Detector.h>
16 #include "DD4hep/DD4hepUnits.h"
17 #include "DDRec/DDGear.h"
18 #include "DDRec/DetectorData.h"
19 #include <DD4hep/DetType.h>
20 #include <DDRec/Vector3D.h>
22 #include "marlin/Global.h"
23 #include <MarlinTrk/IMarlinTrack.h>
24 #include <MarlinTrk/IMarlinTrkSystem.h>
25 #include "MarlinTrk/MarlinTrkUtils.h"
34 dE(_dE), dx(_dx), dEdx(_dE/_dx) {}
37 dE(orig.Get_dE()), dx(orig.Get_dx()), dEdx(orig.Get_dE()/orig.Get_dx()) {}
43 m_trackCollName(
""), m_trkHitCollNames(),
44 surfMap(NULL), trkSystem(NULL), _bField(0),
46 lastRunHeaderProcessed(-1),
48 lastTP(std::chrono::high_resolution_clock::now()),
49 newTP(std::chrono::high_resolution_clock::now())
53 _description =
"SiTracker_dEdxProcessor calculates dE/dx for planar silicon trackers" ;
58 registerInputCollection( LCIO::TRACK,
59 "TrackCollectionName" ,
60 "Name of the input Track collection" ,
62 std::string(
"SiTracks")
67 defaultTrkHitCollections.push_back(std::string(
"ITrackerHits"));
68 defaultTrkHitCollections.push_back(std::string(
"ITrackerEndcapHits"));
69 defaultTrkHitCollections.push_back(std::string(
"OTrackerHits"));
70 defaultTrkHitCollections.push_back(std::string(
"OTrackerEndcapHits"));
71 defaultTrkHitCollections.push_back(std::string(
"VXDTrackerHits"));
72 defaultTrkHitCollections.push_back(std::string(
"VXDEndcapTrackerHits"));
74 registerProcessorParameter(
"TrkHitCollections" ,
75 "Tracker hit collections that will be analysed",
77 defaultTrkHitCollections ) ;
80 for (
unsigned ibit=0; ibit<
sizeof(int)*CHAR_BIT; ibit++) {
81 elementMask += 1 << ibit;
84 registerProcessorParameter(
"CheatSensorThicknesses" ,
85 "Shall we use the sensitive thicknesses from parameters?",
89 FloatVec sensThicknessCheatVals;
90 for (
unsigned i=0; i<defaultTrkHitCollections.size(); i++) {
91 sensThicknessCheatVals.push_back(-1);
94 registerProcessorParameter(
"SensorThicknessCheatValues" ,
95 "Sensor thicknesses to use instead of automatic values from DD4hep (if CheatSensorThicknesses==true).",
97 sensThicknessCheatVals ) ;
102 registerProcessorParameter(
"dEdxEstimator" ,
103 "Type of estimator for dEdx.",
105 std::string(
"median") ) ;
118 streamlog_out(DEBUG) <<
" init called " << std::endl ;
145 streamlog_out(ERROR) <<
"Unknown dE/dx evaluation method " <<
m_dEdxEstimator <<
". Exiting.\n";
149 dd4hep::Detector& theDetector = dd4hep::Detector::getInstance();
151 dd4hep::rec::SurfaceManager& surfMan = *theDetector.extension< dd4hep::rec::SurfaceManager >() ;
152 surfMap = surfMan.map(
"tracker" ) ;
165 const double pos[3]={0,0,0};
166 double bFieldVec[3]={0,0,0};
167 theDetector.field().magneticField(pos,bFieldVec);
168 _bField = bFieldVec[2]/dd4hep::tesla;
172 trkSystem = MarlinTrk::Factory::createMarlinTrkSystem(
"DDKalTest" , marlin::Global::GEAR ,
"" ) ;
174 if(
trkSystem == 0 )
throw EVENT::Exception( std::string(
" Cannot initialize MarlinTrkSystem of Type: ") + std::string(
"DDKalTest") );
176 trkSystem->setOption( MarlinTrk::IMarlinTrkSystem::CFG::useQMS,
true );
178 trkSystem->setOption( MarlinTrk::IMarlinTrkSystem::CFG::usedEdx,
false) ;
179 trkSystem->setOption( MarlinTrk::IMarlinTrkSystem::CFG::useSmoothing,
true) ;
183 gROOT->ProcessLine(
"#include <vector>");
187 for (
unsigned i=0; i<
nTimers; i++) {
188 timers.push_back(std::chrono::duration<double>(std::chrono::duration_values<double>::zero()));
191 streamlog_out(DEBUG) <<
" init done " << std::endl ;
211 if(evt->getEventNumber()%100 == 0) {
212 streamlog_out(MESSAGE) <<
" processing event: " << evt->getEventNumber()
213 <<
" in run: " << evt->getRunNumber() << std::endl ;
221 lastTP = std::chrono::high_resolution_clock::now();
223 LCCollection* tracks = NULL;
227 catch(EVENT::DataNotAvailableException &dataex) {
228 streamlog_out(MESSAGE) <<
"Collection " <<
m_trackCollName <<
" not found. Skipping event #" << evt->getEventNumber() <<
".\n";
241 streamlog_out(WARNING) <<
"None of the requested collections found in event #" << evt->getEventNumber() <<
". Skipping event.\n";
249 int nTracks = tracks->getNumberOfElements() ;
251 for (
int i = 0; i < nTracks; i++)
253 streamlog_out(DEBUG5) <<
"Processing track #" << i <<
".\n";
254 TrackImpl * track =
dynamic_cast<TrackImpl*
>( tracks->getElementAt(i) );
258 EVENT::TrackerHitVec trackhits = track->getTrackerHits();
260 MarlinTrk::IMarlinTrack* marlin_trk =
trkSystem->createTrack();
261 for(
auto it : trackhits ){
262 marlin_trk->addHit(it);
265 const TrackStateImpl *trackState =
dynamic_cast<const TrackStateImpl*
>(track->getTrackState(TrackState::AtFirstHit));
267 streamlog_out(WARNING) <<
"Cannot get track state for track #" << i
268 <<
" in event " << evt->getEventNumber() << std::endl;
269 streamlog_out(WARNING) <<
"Skipping track.\n";
273 marlin_trk->initialise( *trackState,
_bField, MarlinTrk::IMarlinTrack::forward ) ;
279 for(
unsigned int ihit = 0; ihit < trackhits.size(); ihit++) {
284 IMPL::TrackStateImpl ts;
287 marlin_trk->extrapolate(hitpos, ts, chi2, ndf);
293 float tanLambda = ts.getTanLambda();
294 float sinTheta = 1. / sqrt(1.+pow(tanLambda,2));
295 float phi = ts.getPhi();
296 float trackDirX = cos(phi)*sinTheta;
297 float trackDirY = sin(phi)*sinTheta;
298 float trackDirZ = tanLambda*sinTheta;
304 unsigned long cellid = trackhits[ihit]->getCellID0();
305 dd4hep::rec::SurfaceMap::const_iterator surface =
surfMap->find(cellid);
306 if (surface ==
surfMap->end()) {
307 streamlog_out(ERROR) <<
"Cannot find the surface corresponding to track hit ID " << cellid
308 <<
" in event " << evt->getEventNumber() <<
"!\n";
315 double norm = sqrt(trackDir.dot(trackDir)*surfaceNormal.dot(surfaceNormal));
316 if (norm < FLT_MIN)
continue;
317 double cosAngle = fabs(trackDir.dot(surfaceNormal)) / norm ;
320 if (thickness < 0.) {
321 streamlog_out(ERROR) <<
"Could not find hit collection corresponding to hit CellID " << cellid
322 <<
", hit ID " << trackhits.at(ihit)->id() <<
" .\n";
323 streamlog_out(ERROR) <<
"Event #" << evt->getEventNumber() <<
".\n";
326 if (thickness < 0.00001) {
327 streamlog_out(ERROR) <<
"GetThickness returned zero!\n";
331 double effThickness = thickness / cosAngle;
333 dEdxHitVec.push_back(
dEdxPoint(trackhits[ihit]->getEDep(), effThickness) );
340 if(dEdxHitVec.size() == 0)
continue;
342 double dEdx, dEdxError;
343 dEdx =
dEdxEval(dEdxHitVec, dEdxError);
346 track->setdEdx(dEdx);
347 track->setdEdxError(dEdxError);
365 for (
unsigned i=0; i<
timers.size(); i++) {
366 streamlog_out(MESSAGE) <<
"Total time in timer #" << i <<
": " <<
timers.at(i).count() <<
" s\n";
390 const double truncLo,
const double truncHi) {
392 const unsigned n = hitVec.size();
393 const unsigned iStart =
static_cast<unsigned>(floor(n*truncLo + 0.5));
394 const unsigned iEnd =
static_cast<unsigned>(floor(n*(1-truncHi) + 0.5));
396 if(iEnd-iStart == 0) {
400 if(iEnd-iStart == 1) {
401 dEdxError = hitVec.at(iStart).Get_dEdx() ;
402 return hitVec.at(iStart).Get_dEdx() ;
405 sort(hitVec.begin(), hitVec.end(),
dEdxOrder);
408 double thickness = 0.;
410 for (
unsigned i=iStart; i<iEnd; i++) {
411 eDepSum += hitVec.at(i).Get_dE();
412 thickness += hitVec.at(i).Get_dx();
413 mu2dEdx += pow(hitVec.at(i).Get_dE(), 2) / hitVec.at(i).Get_dx();
416 mu2dEdx /= thickness;
418 double track_dEdx = eDepSum / thickness;
419 dEdxError = sqrt((mu2dEdx - pow(track_dEdx, 2))/(iEnd-iStart));
435 const unsigned n = hitVec.size();
441 dEdxError = hitVec.at(0).Get_dEdx();
442 return hitVec.at(0).Get_dEdx();
445 sort(hitVec.begin(), hitVec.end(),
dEdxOrder);
448 median = hitVec.at(n/2).Get_dEdx();
451 median = (hitVec.at(n/2-1).Get_dEdx() + hitVec.at(n/2).Get_dEdx()) / 2;
472 const unsigned n = hitVec.size();
478 dEdxError = hitVec.at(0).Get_dEdx();
479 return hitVec.at(0).Get_dEdx();
486 for (
unsigned i=0; i<n; i++) {
487 double inverse = 1/hitVec.at(i).Get_dEdx();
489 mu2sum += pow( inverse, 2 );
492 double mu2 = mu2sum / n;
493 double mu1 = mu1sum / n;
494 double sigma = sqrt( (mu2 - pow(mu1, 2)) / n );
496 double dEdx = 1 / mu1;
497 dEdxError = sigma / pow(mu1, 2) ;
506 const unsigned n = hitVec.size();
512 dEdxError = hitVec.at(0).Get_dEdx();
513 return hitVec.at(0).Get_dEdx();
520 for (
unsigned i=0; i<n; i++) {
521 double sqinverse = pow(hitVec.at(i).Get_dEdx(), -2);
523 mu2sum += pow( sqinverse, 2 );
526 double mu2 = mu2sum / n;
527 double mu1 = mu1sum / n;
528 double sigma = sqrt( (mu2 - pow(mu1, 2)) / n );
530 double dEdx = 1 / sqrt(mu1);
531 dEdxError = sigma * pow(dEdx, 3) / 2 ;
541 const unsigned n = hitVec.size();
547 dEdxError = hitVec.at(0).Get_dEdx();
548 return hitVec.at(0).Get_dEdx();
556 for (
unsigned i=0; i<n; i++) {
557 double inverse = 1/hitVec.at(i).Get_dEdx();
558 double wgt = hitVec.at(i).Get_dx();
559 mu1sum += wgt*inverse;
560 mu2sum += wgt*pow( inverse, 2 );
564 double mu2 = mu2sum / wgtsum;
565 double mu1 = mu1sum / wgtsum;
566 double sigma = sqrt( (mu2 - pow(mu1, 2)) / n );
568 double dEdx = 1 / mu1;
569 dEdxError = sigma * pow(dEdx, 2) ;
579 const unsigned n = hitVec.size();
585 dEdxError = hitVec.at(0).Get_dEdx();
586 return hitVec.at(0).Get_dEdx();
594 for (
unsigned i=0; i<n; i++) {
595 double sqinverse = pow(hitVec.at(i).Get_dEdx(), -2);
596 double wgt = hitVec.at(i).Get_dx();
597 mu1sum += wgt*sqinverse;
598 mu2sum += wgt*pow( sqinverse, 2 );
602 double mu2 = mu2sum / wgtsum;
603 double mu1 = mu1sum / wgtsum;
604 double sigma = sqrt( (mu2 - pow(mu1, 2)) / n );
606 double dEdx = 1 / sqrt(mu1);
607 dEdxError = sigma * pow(dEdx, 3) / 2 ;
const dd4hep::rec::SurfaceMap * surfMap
int lastRunHeaderProcessed
static double dEdxHarmonic2(dEdxVec, double &dEdxError)
static double truncFractionUp
StringVec m_trkHitCollNames
CLHEP::Hep3Vector Vector3D
static double dEdxHarmonic(dEdxVec, double &dEdxError)
virtual void end()
Called after data processing for clean up.
virtual void check(LCEvent *evt)
FloatVec m_sensThicknessCheatVals
static double dEdxTruncMean(dEdxVec, double &dEdxError)
std::vector< std::chrono::duration< double > > timers
SiTracker_dEdxProcessor aSiTracker_dEdxProcessor
virtual void processRunHeader(LCRunHeader *run)
Called at the end of every run.
SiTracker_dEdxProcessor()
static double dEdxGeneralTruncMean(dEdxVec, double &dEdxError, const double truncLo=0, const double truncHi=0)
int ReadCollections(EVENT::LCEvent *)
static double truncFractionLo
virtual void processEvent(LCEvent *evt)
Called for every event - the working horse.
LayerFinder * layerFinder
std::string m_dEdxEstimator
bool m_cheatSensorThicknesses
static double dEdxWgtHarmonic2(dEdxVec, double &dEdxError)
static double dEdxWgtHarmonic(dEdxVec, double &dEdxError)
virtual ~SiTracker_dEdxProcessor()
MarlinTrk::IMarlinTrkSystem * trkSystem
SiTracker_dEdxProcessor for Marlin.
static double dEdxMean(dEdxVec, double &dEdxError)
std::vector< dEdxPoint > dEdxVec
std::chrono::high_resolution_clock::time_point lastTP
double SensitiveThickness(lcio::TrackerHitPlane *)
std::string m_trackCollName
static bool dEdxOrder(dEdxPoint p1, dEdxPoint p2)
static double dEdxMedian(dEdxVec, double &dEdxError)
dEdxPoint(const double _dE, const double _dx)
std::vector< std::string > StringVec
virtual void init()
Called at the begin of the job before anything is read.