10 #include <tracking/trackExtrapolateG4e/TrackExtrapolateG4e.h>
13 #include <ecl/geometry/ECLGeometryPar.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <framework/datastore/StoreObjPtr.h>
16 #include <framework/geometry/BFieldManager.h>
17 #include <framework/logging/Logger.h>
18 #include <genfit/Exception.h>
19 #include <klm/dataobjects/bklm/BKLMElementNumbers.h>
20 #include <klm/dataobjects/bklm/BKLMStatus.h>
21 #include <klm/bklm/geometry/GeometryPar.h>
22 #include <klm/bklm/geometry/Module.h>
23 #include <klm/dataobjects/KLMChannelIndex.h>
24 #include <klm/dataobjects/KLMMuidLikelihood.h>
25 #include <klm/dataobjects/KLMMuidHit.h>
26 #include <klm/dataobjects/eklm/EKLMAlignmentHit.h>
27 #include <klm/muid/MuidBuilder.h>
28 #include <klm/muid/MuidElementNumbers.h>
29 #include <mdst/dataobjects/ECLCluster.h>
30 #include <mdst/dataobjects/KLMCluster.h>
31 #include <mdst/dataobjects/Track.h>
32 #include <simulation/kernel/ExtCylSurfaceTarget.h>
33 #include <simulation/kernel/ExtManager.h>
36 #include <CLHEP/Matrix/Vector.h>
37 #include <CLHEP/Units/PhysicalConstants.h>
38 #include <CLHEP/Units/SystemOfUnits.h>
41 #include <G4ErrorFreeTrajState.hh>
42 #include <G4ErrorMatrix.hh>
43 #include <G4ErrorPropagatorData.hh>
44 #include <G4ErrorSymMatrix.hh>
45 #include <G4ParticleTable.hh>
46 #include <G4PhysicalVolumeStore.hh>
47 #include <G4Point3D.hh>
48 #include <G4StateManager.hh>
50 #include <G4StepPoint.hh>
51 #include <G4ThreeVector.hh>
52 #include <G4TouchableHandle.hh>
54 #include <G4UImanager.hh>
55 #include <G4VPhysicalVolume.hh>
56 #include <G4VTouchable.hh>
63 #define TWOPI (2.0*M_PI)
64 #define PI_8 (0.125*M_PI)
66 #define DEPTH_SCINT 11
80 m_ExtInitialized(false),
81 m_MuidInitialized(false),
85 m_MaxDistSqInVariances(0.0),
86 m_MaxKLMTrackClusterDistance(0.0),
87 m_MaxECLTrackClusterDistance(0.0),
91 m_HypothesesExt(nullptr),
92 m_HypothesesMuid(nullptr),
93 m_DefaultHypotheses(nullptr),
95 m_BKLMVolumes(nullptr),
97 m_TargetMuid(nullptr),
103 m_BarrelHalfLength(0.0),
104 m_OutermostActiveBarrelLayer(0),
105 m_BarrelScintVariance(0.0),
108 m_EndcapMiddleZ(0.0),
109 m_EndcapHalfLength(0.0),
110 m_OutermostActiveForwardEndcapLayer(0),
111 m_OutermostActiveBackwardEndcapLayer(0),
112 m_EndcapScintVariance(0.0),
113 m_eklmTransformData(nullptr)
139 std::vector<Const::ChargedStable>& hypotheses)
153 m_MinPt = std::max(0.0, minPt) * CLHEP::GeV;
154 m_MinKE = std::max(0.0, minKE) * CLHEP::GeV;
169 B2FATAL(
"Coil geometry data are not available.");
179 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetExt);
182 B2FATAL(
"Beam pipe geometry data are not available.");
183 double beampipeRadius =
m_BeamPipeGeo->getParameter(
"Lv2OutBe.R2") * CLHEP::cm;
185 beampipeRadius =
m_BeamPipeGeo->getParameter(
"Lv2OutBe.R2") * CLHEP::cm;
192 double maxKLMTrackClusterDistance,
double maxECLTrackClusterDistance,
193 double minPt,
double minKE,
bool addHitsToRecoTrack,
194 std::vector<Const::ChargedStable>& hypotheses)
234 m_MinPt = std::max(0.0, minPt) * CLHEP::GeV;
235 m_MinKE = std::max(0.0, minKE) * CLHEP::GeV;
250 B2FATAL(
"Coil geometry data are not available.");
260 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetExt);
263 B2FATAL(
"Beam pipe geometry data are not available.");
264 double beampipeRadius =
m_BeamPipeGeo->getParameter(
"Lv2OutBe.R2") * CLHEP::cm;
266 beampipeRadius =
m_BeamPipeGeo->getParameter(
"Lv2OutBe.R2") * CLHEP::cm;
281 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetMuid);
300 width = module->getPhiStripWidth();
302 width = module->getZStripWidth();
309 double phi = M_PI_4 * (sector - 1);
317 bklmGeometry->
getActiveMiddleRadius(klmLayer.getSection(), klmLayer.getSector(), klmLayer.getLayer());
337 B2DEBUG(20, (byMuid ?
"muid" :
"ext"));
340 B2FATAL(
"KLM channel status data are not available.");
342 B2FATAL(
"KLM strip efficiency data are not available.");
344 B2FATAL(
"KLM likelihood parameters are not available.");
354 for (
int pdg : muidPdgCodes)
363 if (G4StateManager::GetStateManager()->GetCurrentState() == G4State_Idle) {
364 G4StateManager::GetStateManager()->SetNewState(G4State_GeomClosed);
367 G4ThreeVector directionAtIP, positionG4e, momentumG4e;
368 G4ErrorTrajErr covG4e(5);
373 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetMuid);
374 std::vector<std::pair<ECLCluster*, G4ThreeVector> > eclClusterInfo(
m_eclClusters.getEntries());
377 eclClusterInfo[c].second = G4ThreeVector(
m_eclClusters[c]->getClusterPosition().X(),
381 std::vector<std::pair<KLMCluster*, G4ThreeVector> > klmClusterInfo(
m_klmClusters.getEntries());
384 klmClusterInfo[c].second = G4ThreeVector(
m_klmClusters[c]->getClusterPosition().X(),
389 std::vector<std::map<const Track*, double> > bklmHitUsed(
m_bklmHit2ds.getEntries());
392 int pdgCode = hypothesis.getPDGCode();
395 G4ErrorFreeTrajState g4eState(
"g4e_mu+", G4ThreeVector(), G4ThreeVector());
397 swim(extState, g4eState, &eclClusterInfo, &klmClusterInfo, &bklmHitUsed);
401 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetExt);
404 int pdgCode = hypothesis.getPDGCode();
406 G4ErrorFreeTrajState g4eState(
"g4e_mu+", G4ThreeVector(), G4ThreeVector());
408 swim(extState, g4eState);
444 const G4ThreeVector& position,
445 const G4ThreeVector& momentum,
446 const G4ErrorSymMatrix& covariance)
448 bool isCosmic =
false;
453 extMgr->
Initialize(
"Ext",
"default", 0.0, 0.25,
false, 0, std::vector<std::string>());
458 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/stepLength 250 mm");
459 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/magField 0.001");
460 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/energyLoss 0.05");
466 if (G4StateManager::GetStateManager()->GetCurrentState() == G4State_Idle) {
467 G4StateManager::GetStateManager()->SetNewState(G4State_GeomClosed);
470 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetExt);
475 G4ThreeVector positionG4e = position * CLHEP::cm;
476 G4ThreeVector momentumG4e = momentum * CLHEP::GeV;
479 momentumG4e = -momentumG4e;
480 G4ErrorSymMatrix covarianceG4e(5, 0);
482 G4String nameG4e(
"g4e_" + G4ParticleTable::GetParticleTable()->FindParticle(pdgCode)->GetParticleName());
483 G4ErrorFreeTrajState g4eState(nameG4e, positionG4e, momentumG4e, covarianceG4e);
484 ExtState extState = {
nullptr, pdgCode, isCosmic, tof, 0.0,
485 momentumG4e.unit(), 0.0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
false
487 swim(extState, g4eState);
493 const G4ThreeVector& position,
494 const G4ThreeVector& momentum,
495 const G4ErrorSymMatrix& covariance)
502 extMgr->
Initialize(
"Muid",
"default", 0.0, 0.25,
false, 0, std::vector<std::string>());
507 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/stepLength 250 mm");
508 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/magField 0.001");
509 G4UImanager::GetUIpointer()->ApplyCommand(
"/geant4e/limits/energyLoss 0.05");
515 if (G4StateManager::GetStateManager()->GetCurrentState() == G4State_Idle) {
516 G4StateManager::GetStateManager()->SetNewState(G4State_GeomClosed);
519 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(
m_TargetMuid);
524 G4ThreeVector positionG4e = position * CLHEP::cm;
525 G4ThreeVector momentumG4e = momentum * CLHEP::GeV;
527 momentumG4e = -momentumG4e;
528 G4ErrorSymMatrix covarianceG4e(5, 0);
530 G4String nameG4e(
"g4e_" + G4ParticleTable::GetParticleTable()->FindParticle(pdgCode)->GetParticleName());
531 G4ErrorFreeTrajState g4eState(nameG4e, positionG4e, momentumG4e, covarianceG4e);
532 ExtState extState = {
nullptr, pdgCode, isCosmic, tof, 0.0,
533 momentumG4e.unit(), 0.0, 0, 0, 0, -1, -1, -1, -1, 0, 0,
false
535 swim(extState, g4eState,
nullptr,
nullptr,
nullptr);
540 const std::vector<std::pair<ECLCluster*, G4ThreeVector> >* eclClusterInfo,
541 const std::vector<std::pair<KLMCluster*, G4ThreeVector> >* klmClusterInfo,
542 std::vector<std::map<const Track*, double> >* bklmHitUsed)
546 if (g4eState.GetMomentum().perp() <=
m_MinPt)
550 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.
pdgCode);
551 double mass = particle->GetPDGMass();
552 double minPSq = (mass +
m_MinKE) * (mass +
m_MinKE) - mass * mass;
554 std::vector<double> eclClusterDistance;
555 std::vector<ExtHit> eclHit1, eclHit2, eclHit3;
556 if (eclClusterInfo !=
nullptr) {
557 eclClusterDistance.resize(eclClusterInfo->size(), 1.0E10);
558 ExtHit tempExtHit(extState.
pdgCode, Const::EDetector::ECL, 0, EXT_FIRST,
560 G4ThreeVector(), G4ThreeVector(), G4ErrorSymMatrix(6));
561 eclHit1.resize(eclClusterInfo->size(), tempExtHit);
562 eclHit2.resize(eclClusterInfo->size(), tempExtHit);
563 eclHit3.resize(eclClusterInfo->size(), tempExtHit);
566 std::vector<TrackClusterSeparation> klmHit;
567 if (klmClusterInfo !=
nullptr) {
568 klmHit.resize(klmClusterInfo->size());
572 if (extState.
track !=
nullptr)
574 G4ErrorMode propagationMode = (extState.
isCosmic ? G4ErrorMode_PropBackwards : G4ErrorMode_PropForwards);
578 G4Track* track = g4eState.GetG4Track();
579 const G4Step* step = track->GetStep();
580 const G4StepPoint* preStepPoint = step->GetPreStepPoint();
581 const G4StepPoint* postStepPoint = step->GetPostStepPoint();
582 G4TouchableHandle preTouch = preStepPoint->GetTouchableHandle();
583 G4VPhysicalVolume* pVol = preTouch->GetVolume();
584 const G4int preStatus = preStepPoint->GetStepStatus();
585 const G4int postStatus = postStepPoint->GetStepStatus();
586 G4ThreeVector pos = track->GetPosition();
587 G4ThreeVector mom = track->GetMomentum();
590 if (step->GetStepLength() > 0.0) {
591 double dt = step->GetDeltaTime();
592 double dl = step->GetStepLength() / track->GetMaterial()->GetRadlen();
593 if (preStatus == fGeomBoundary) {
596 createExtHit(EXT_ENTER, extState, g4eState, preStepPoint, preTouch);
607 if (postStatus == fGeomBoundary) {
610 createExtHit(EXT_EXIT, extState, g4eState, postStepPoint, preTouch);
614 if (
createMuidHit(extState, g4eState, klmMuidLikelihood, bklmHitUsed)) {
618 if (eclClusterInfo !=
nullptr) {
619 for (
unsigned int c = 0; c < eclClusterInfo->size(); ++c) {
620 G4ThreeVector eclPos((*eclClusterInfo)[c].second);
621 G4ThreeVector prePos(preStepPoint->GetPosition());
622 G4ThreeVector diff(prePos - eclPos);
623 double distance = diff.mag();
626 if (distance < eclClusterDistance[c]) {
627 eclClusterDistance[c] = distance;
628 G4ErrorSymMatrix covariance(6, 0);
630 eclHit3[c].update(EXT_ECLNEAR, extState.
tof, pos / CLHEP::cm, mom / CLHEP::GeV, covariance);
633 if (eclHit1[c].getStatus() == EXT_FIRST) {
634 if (pos.mag2() >= eclPos.mag2()) {
635 double r = eclPos.mag();
636 double preD = prePos.mag() - r;
637 double postD = pos.mag() - r;
638 double f = postD / (postD - preD);
639 G4ThreeVector midPos = pos + (prePos - pos) * f;
640 double tof = extState.
tof + dt * f * (extState.
isCosmic ? +1 : -1);
641 G4ErrorSymMatrix covariance(6, 0);
643 eclHit1[c].update(EXT_ECLCROSS, tof, midPos / CLHEP::cm, mom / CLHEP::GeV, covariance);
648 if (eclHit2[c].getStatus() == EXT_FIRST) {
649 G4ThreeVector delta(pos - prePos);
650 G4ThreeVector perp(eclPos.cross(delta));
651 double perpMag2 = perp.mag2();
652 if (perpMag2 > 1.0E-10) {
653 double dist = std::fabs(diff * perp) / std::sqrt(perpMag2);
655 double f = eclPos * (prePos.cross(perp)) / perpMag2;
656 if ((f > -0.5) && (f <= 1.0)) {
657 G4ThreeVector midPos(prePos + f * delta);
658 double length = extState.
length + dl * (1.0 - f) * (extState.
isCosmic ? +1 : -1);
659 G4ErrorSymMatrix covariance(6, 0);
661 eclHit2[c].update(EXT_ECLDL, length, midPos / CLHEP::cm, mom / CLHEP::GeV, covariance);
668 if (klmClusterInfo !=
nullptr) {
669 for (
unsigned int c = 0; c < klmClusterInfo->size(); ++c) {
670 G4ThreeVector klmPos = (*klmClusterInfo)[c].second;
671 G4ThreeVector separation = klmPos - pos;
672 double distance = separation.mag();
673 if (distance < klmHit[c].getDistance()) {
674 klmHit[c].setDistance(distance);
675 klmHit[c].setTrackClusterAngle(mom.angle(separation));
676 klmHit[c].setTrackClusterSeparationAngle(mom.angle(klmPos));
677 klmHit[c].setTrackRotationAngle(extState.
directionAtIP.angle(mom));
678 klmHit[c].setTrackClusterInitialSeparationAngle(extState.
directionAtIP.angle(klmPos));
684 if (errCode || (mom.mag2() < minPSq)) {
701 if (eclClusterInfo !=
nullptr) {
702 for (
unsigned int c = 0; c < eclClusterInfo->size(); ++c) {
703 if (eclHit1[c].getStatus() != EXT_FIRST) {
705 (*eclClusterInfo)[c].first->addRelationTo(h);
708 if (eclHit2[c].getStatus() != EXT_FIRST) {
710 (*eclClusterInfo)[c].first->addRelationTo(h);
713 if (eclHit3[c].getStatus() != EXT_FIRST) {
715 (*eclClusterInfo)[c].first->addRelationTo(h);
721 if (klmClusterInfo !=
nullptr) {
725 unsigned int closestCluster = 0;
726 for (
unsigned int c = 0; c < klmClusterInfo->size(); ++c) {
727 if (klmHit[c].getDistance() > 1.0E9) {
731 (*klmClusterInfo)[c].first->addRelationTo(h);
733 if (klmHit[c].getDistance() < minDistance) {
735 minDistance = klmHit[c].getDistance();
741 extState.
track->
addRelationTo((*klmClusterInfo)[closestCluster].first, 1. / (minDistance / CLHEP::cm));
752 if (g4eState.GetMomentum().perp() <=
m_MinPt)
756 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.
pdgCode);
757 double mass = particle->GetPDGMass();
758 double minPSq = (mass +
m_MinKE) * (mass +
m_MinKE) - mass * mass;
759 G4ErrorMode propagationMode = (extState.
isCosmic ? G4ErrorMode_PropBackwards : G4ErrorMode_PropForwards);
763 G4Track* track = g4eState.GetG4Track();
764 const G4Step* step = track->GetStep();
765 const G4StepPoint* preStepPoint = step->GetPreStepPoint();
766 const G4StepPoint* postStepPoint = step->GetPostStepPoint();
767 G4TouchableHandle preTouch = preStepPoint->GetTouchableHandle();
768 G4VPhysicalVolume* pVol = preTouch->GetVolume();
769 const G4int preStatus = preStepPoint->GetStepStatus();
770 const G4int postStatus = postStepPoint->GetStepStatus();
771 G4ThreeVector pos = track->GetPosition();
772 G4ThreeVector mom = track->GetMomentum();
776 if (preStatus == fUndefined) {
778 createExtHit(EXT_FIRST, extState, g4eState, preStepPoint, preTouch);
782 if (step->GetStepLength() > 0.0) {
783 double dt = step->GetDeltaTime();
784 double dl = step->GetStepLength() / track->GetMaterial()->GetRadlen();
785 if (preStatus == fGeomBoundary) {
787 createExtHit(EXT_ENTER, extState, g4eState, preStepPoint, preTouch);
798 if (postStatus == fGeomBoundary) {
800 createExtHit(EXT_EXIT, extState, g4eState, postStepPoint, preTouch);
805 if (errCode || (mom.mag2() < minPSq)) {
807 createExtHit(EXT_STOP, extState, g4eState, postStepPoint, preTouch);
814 createExtHit(EXT_ESCAPE, extState, g4eState, postStepPoint, preTouch);
831 G4PhysicalVolumeStore* pvStore = G4PhysicalVolumeStore::GetInstance();
832 if (pvStore->size() == 0) {
833 B2FATAL(
"No geometry defined. Please create the geometry first.");
837 m_EnterExit =
new std::map<G4VPhysicalVolume*, enum VolTypes>;
839 for (std::vector<G4VPhysicalVolume*>::iterator iVol = pvStore->begin();
840 iVol != pvStore->end(); ++iVol) {
841 const G4String name = (*iVol)->GetName();
844 if (name.find(
"TOPModule") != std::string::npos) {
848 else if (name.find(
"_TOPPrism_") != std::string::npos ||
849 name.find(
"_TOPBarSegment") != std::string::npos ||
850 name.find(
"_TOPMirrorSegment") != std::string::npos) {
854 else if (name.find(
"TOPBarSegment1Glue") != std::string::npos ||
855 name.find(
"TOPBarSegment2Glue") != std::string::npos ||
856 name.find(
"TOPMirrorSegmentGlue") != std::string::npos) {
859 }
else if (name ==
"ARICH.AerogelSupportPlate") {
861 }
else if (name ==
"ARICH.AerogelImgPlate") {
863 }
else if (name.find(
"ARICH.HAPDWindow") != std::string::npos) {
867 else if (name.find(
"lv_barrel_crystal_") != std::string::npos ||
868 name.find(
"lv_forward_crystal_") != std::string::npos ||
869 name.find(
"lv_backward_crystal_") != std::string::npos) {
874 else if (name.compare(0, 5,
"BKLM.") == 0) {
875 if (name.find(
"GasPhysical") != std::string::npos) {
877 }
else if (name.find(
"ScintActiveType") != std::string::npos) {
879 }
else if ((name.find(
"ScintType") != std::string::npos) ||
880 (name.find(
"ElectrodePhysical") != std::string::npos)) {
885 else if (name.compare(0, 14,
"StripSensitive") == 0) {
897 detID = Const::EDetector::invalidDetector;
900 G4VPhysicalVolume* pv = touch->GetVolume(0);
901 std::map<G4VPhysicalVolume*, enum VolTypes>::iterator it =
m_EnterExit->find(pv);
905 switch (it->second) {
907 detID = Const::EDetector::CDC;
908 copyID = pv->GetCopyNo();
911 detID = Const::EDetector::TOP;
912 copyID = -(pv->GetCopyNo());
915 detID = Const::EDetector::TOP;
916 if (touch->GetHistoryDepth() >= 1)
917 copyID = touch->GetVolume(1)->GetCopyNo();
920 detID = Const::EDetector::TOP;
921 if (touch->GetHistoryDepth() >= 2)
922 copyID = touch->GetVolume(2)->GetCopyNo();
925 detID = Const::EDetector::ARICH;
929 detID = Const::EDetector::ARICH;
933 detID = Const::EDetector::ARICH;
934 if (touch->GetHistoryDepth() >= 2)
935 copyID = touch->GetVolume(2)->GetCopyNo();
938 detID = Const::EDetector::ECL;
942 detID = Const::EDetector::BKLM;
943 if (touch->GetHistoryDepth() == DEPTH_RPC) {
945 int layer = touch->GetCopyNumber(4);
946 int sector = touch->GetCopyNumber(6);
947 int section = touch->GetCopyNumber(7);
952 detID = Const::EDetector::BKLM;
953 if (touch->GetHistoryDepth() == DEPTH_SCINT) {
954 int strip = touch->GetCopyNumber(1);
955 int plane = (touch->GetCopyNumber(2) == BKLM_INNER) ?
958 int layer = touch->GetCopyNumber(6);
959 int sector = touch->GetCopyNumber(8);
960 int section = touch->GetCopyNumber(9);
962 section, sector, layer, plane, strip);
967 detID = Const::EDetector::EKLM;
969 touch->GetVolume(7)->GetCopyNo(),
970 touch->GetVolume(6)->GetCopyNo(),
971 touch->GetVolume(5)->GetCopyNo(),
972 touch->GetVolume(4)->GetCopyNo(),
973 touch->GetVolume(1)->GetCopyNo());
981 ExtState extState = {&b2track, pdgCode,
false, 0.0, 0.0,
982 G4ThreeVector(0, 0, 1), 0.0, 0, 0, 0, -1, -1, -1, -1, 0, 0,
false
985 if (recoTrack ==
nullptr) {
986 B2WARNING(
"Track without associated RecoTrack: skipping extrapolation for this track.");
993 B2WARNING(
"RecoTrack fit failed for cardinal representation: skipping extrapolation for this track.");
1003 TVector3 firstPosition, firstMomentum, lastPosition, lastMomentum;
1004 TMatrixDSym firstCov(6), lastCov(6);
1007 trackRep->
getPosMomCov(firstState, firstPosition, firstMomentum, firstCov);
1009 trackRep->
getPosMomCov(lastState, lastPosition, lastMomentum, lastCov);
1011 extState.
tof = lastState.getTime();
1012 if (lastPosition.Mag2() < firstPosition.Mag2()) {
1013 firstPosition = lastPosition;
1014 firstMomentum = -lastMomentum;
1016 trackRep->
getPosMomCov(firstState, lastPosition, lastMomentum, lastCov);
1017 lastMomentum *= -1.0;
1019 extState.
tof = firstState.getTime();
1022 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.
pdgCode);
1024 double pSq = lastMomentum.Mag2();
1025 double mass = particle->GetPDGMass() / CLHEP::GeV;
1026 extState.
tof *= std::sqrt((pSq + mass * mass) / (pSq + lastState.getMass() * lastState.getMass()));
1029 extState.
directionAtIP.set(firstMomentum.Unit().X(), firstMomentum.Unit().Y(), firstMomentum.Unit().Z());
1031 double radius = (firstMomentum.Perp() * CLHEP::GeV / CLHEP::eV) / (CLHEP::c_light * charge *
m_MagneticField);
1033 double centerX = firstPosition.X() + radius * std::cos(centerPhi);
1034 double centerY = firstPosition.Y() + radius * std::sin(centerPhi);
1035 double pocaPhi = atan2(charge * centerY, charge * centerX) + M_PI;
1045 lastCov[0][0] = 5.0E-5;
1046 lastCov[1][1] = 1.0E-7;
1047 lastCov[2][2] = 5.0E-4;
1048 lastCov[3][3] = 3.5E-3;
1049 lastCov[4][4] = 3.5E-3;
1050 lastMomentum = lastMomentum.Unit() * 10.0;
1053 G4ThreeVector posG4e(lastPosition.X() * CLHEP::cm, lastPosition.Y() * CLHEP::cm,
1054 lastPosition.Z() * CLHEP::cm);
1055 G4ThreeVector momG4e(lastMomentum.X() * CLHEP::GeV, lastMomentum.Y() * CLHEP::GeV,
1056 lastMomentum.Z() * CLHEP::GeV);
1057 G4ErrorSymMatrix covG4e(5, 0);
1059 g4eState.SetData(
"g4e_" + particle->GetParticleName(), posG4e, momG4e);
1060 g4eState.SetParameters(posG4e, momG4e);
1061 g4eState.SetError(covG4e);
1063 B2WARNING(
"genfit::MeasuredStateOnPlane() exception: skipping extrapolation for this track. initial momentum = ("
1064 << firstMomentum.X() <<
"," << firstMomentum.Y() <<
"," << firstMomentum.Z() <<
")");
1084 G4ErrorFreeTrajParam param = g4eState.GetParameters();
1085 double p = 1.0 / (param.GetInvP() * CLHEP::GeV);
1087 double lambda = param.GetLambda();
1088 double sinLambda = std::sin(lambda);
1089 double cosLambda = std::cos(lambda);
1090 double phi = param.GetPhi();
1091 double sinPhi = std::sin(phi);
1092 double cosPhi = std::cos(phi);
1096 G4ErrorMatrix jacobian(6, 5, 0);
1098 jacobian(4, 1) = -pSq * cosLambda * cosPhi;
1099 jacobian(5, 1) = -pSq * cosLambda * sinPhi;
1100 jacobian(6, 1) = -pSq * sinLambda;
1102 jacobian(4, 2) = -p * sinLambda * cosPhi;
1103 jacobian(5, 2) = -p * sinLambda * sinPhi;
1104 jacobian(6, 2) = p * cosLambda;
1106 jacobian(4, 3) = -p * cosLambda * sinPhi;
1107 jacobian(5, 3) = p * cosLambda * cosPhi;
1109 jacobian(1, 4) = -sinPhi;
1110 jacobian(2, 4) = cosPhi;
1112 jacobian(1, 5) = -sinLambda * cosPhi;
1113 jacobian(2, 5) = -sinLambda * sinPhi;
1114 jacobian(3, 5) = cosLambda;
1116 G4ErrorTrajErr g4eCov = g4eState.GetError();
1117 covariance.assign(g4eCov.similarity(jacobian));
1133 G4ErrorSymMatrix temp(6, 0);
1134 for (
int k = 0; k < 6; ++k) {
1135 for (
int j = k; j < 6; ++j) {
1136 temp[j][k] = covariance[j][k];
1140 double pInvSq = 1.0 / momentum.Mag2();
1141 double pInv = std::sqrt(pInvSq);
1142 double pPerpInv = 1.0 / momentum.Perp();
1143 double sinLambda = momentum.CosTheta();
1144 double cosLambda = std::sqrt(1.0 - sinLambda * sinLambda);
1145 double phi = momentum.Phi();
1146 double cosPhi = std::cos(phi);
1147 double sinPhi = std::sin(phi);
1150 G4ErrorMatrix jacobian(5, 6, 0);
1152 jacobian(1, 4) = -pInvSq * cosLambda * cosPhi;
1153 jacobian(1, 5) = -pInvSq * cosLambda * sinPhi;
1154 jacobian(1, 6) = -pInvSq * sinLambda;
1156 jacobian(2, 4) = -pInv * sinLambda * cosPhi;
1157 jacobian(2, 5) = -pInv * sinLambda * sinPhi;
1158 jacobian(2, 6) = pInv * cosLambda;
1160 jacobian(3, 4) = -pPerpInv * sinPhi;
1161 jacobian(3, 5) = pPerpInv * cosPhi;
1163 jacobian(4, 1) = -sinPhi;
1164 jacobian(4, 2) = cosPhi;
1166 jacobian(5, 1) = -sinLambda * cosPhi;
1167 jacobian(5, 2) = -sinLambda * sinPhi;
1168 jacobian(5, 3) = cosLambda;
1170 covG4e = temp.similarity(jacobian);
1175 G4ErrorTrajErr& covG4e)
1187 G4ErrorSymMatrix temp(covariance);
1189 double pInvSq = 1.0 / momentum.mag2();
1190 double pInv = std::sqrt(pInvSq);
1191 double pPerpInv = 1.0 / momentum.perp();
1192 double sinLambda = momentum.cosTheta();
1193 double cosLambda = std::sqrt(1.0 - sinLambda * sinLambda);
1194 double phi = momentum.phi();
1195 double cosPhi = std::cos(phi);
1196 double sinPhi = std::sin(phi);
1199 G4ErrorMatrix jacobian(5, 6, 0);
1201 jacobian(1, 4) = -pInvSq * cosLambda * cosPhi;
1202 jacobian(1, 5) = -pInvSq * cosLambda * sinPhi;
1203 jacobian(1, 6) = -pInvSq * sinLambda;
1205 jacobian(2, 4) = -pInv * sinLambda * cosPhi;
1206 jacobian(2, 5) = -pInv * sinLambda * sinPhi;
1207 jacobian(2, 6) = pInv * cosLambda;
1209 jacobian(3, 4) = -pPerpInv * sinPhi;
1210 jacobian(3, 5) = pPerpInv * cosPhi;
1212 jacobian(4, 1) = -sinPhi;
1213 jacobian(4, 2) = cosPhi;
1215 jacobian(5, 1) = -sinLambda * cosPhi;
1216 jacobian(5, 2) = -sinLambda * sinPhi;
1217 jacobian(5, 3) = cosLambda;
1218 covG4e = temp.similarity(jacobian);
1224 const G4ErrorFreeTrajState& g4eState,
1225 const G4StepPoint* stepPoint,
const G4TouchableHandle& touch)
1230 G4ThreeVector pos(stepPoint->GetPosition() / CLHEP::cm);
1231 G4ThreeVector mom(stepPoint->GetMomentum() / CLHEP::GeV);
1234 G4ErrorSymMatrix covariance(6, 0);
1238 pos, mom, covariance);
1240 if (extState.
track !=
nullptr)
1248 std::vector<std::map<const Track*, double> >* bklmHitUsed)
1252 intersection.
hit = -1;
1253 intersection.
chi2 = -1.0;
1254 intersection.
position = g4eState.GetPosition() / CLHEP::cm;
1255 intersection.
momentum = g4eState.GetMomentum() / CLHEP::GeV;
1256 G4ThreeVector prePos = g4eState.GetG4Track()->GetStep()->GetPreStepPoint()->GetPosition() / CLHEP::cm;
1257 G4ThreeVector oldPosition(prePos.x(), prePos.y(), prePos.z());
1258 double r = intersection.
position.perp();
1267 (*bklmHitUsed)[intersection.
hit].insert(std::pair<const Track*, double>(extState.
track, intersection.
chi2));
1269 float layerBarrelEfficiency = 1.;
1273 intersection.
sector + 1, intersection.
layer + 1, plane, 1);
1287 intersection.
chi2 = -1.0;
1292 g4eState.GetG4Track()->GetVolume());
1298 int sector = intersection.
sector + 1;
1299 int layer = intersection.
layer + 1;
1302 const CLHEP::Hep3Vector localPosition = m->globalToLocal(intersection.
position);
1303 int zStrip = m->getZStripNumber(localPosition);
1304 int phiStrip = m->getPhiStripNumber(localPosition);
1305 if (zStrip >= 0 && phiStrip >= 0) {
1306 uint16_t channel1, channel2;
1308 section, sector, layer,
1311 section, sector, layer,
1318 B2ERROR(
"No KLM channel status data."
1319 <<
LogVar(
"Section", section)
1320 <<
LogVar(
"Sector", sector)
1321 <<
LogVar(
"Layer", layer)
1322 <<
LogVar(
"Z strip", zStrip)
1323 <<
LogVar(
"Phi strip", phiStrip));
1330 float layerBarrelEfficiency = 1.;
1334 intersection.
sector + 1, intersection.
layer + 1, plane, 1);
1355 float layerEndcapEfficiency = 1.;
1359 intersection.
sector + 1, intersection.
layer + 1, plane, 1);
1373 intersection.
chi2 = -1.0;
1377 int result, strip1, strip2;
1379 intersection.
position, &strip1, &strip2);
1381 uint16_t channel1, channel2;
1389 B2ERROR(
"Incomplete KLM channel status data.");
1395 float layerEndcapEfficiency = 1.;
1399 intersection.
sector + 1, intersection.
layer + 1, plane, 1);
1414 if (intersection.
chi2 >= 0.0) {
1421 intersection.
layer, tpos,
1422 tposAtHitPlane, extState.
tof, intersection.
time, intersection.
chi2);
1424 G4Point3D newPos(intersection.
position.x() * CLHEP::cm,
1425 intersection.
position.y() * CLHEP::cm,
1426 intersection.
position.z() * CLHEP::cm);
1427 g4eState.SetPosition(newPos);
1428 G4Vector3D newMom(intersection.
momentum.x() * CLHEP::GeV,
1429 intersection.
momentum.y() * CLHEP::GeV,
1430 intersection.
momentum.z() * CLHEP::GeV);
1431 g4eState.SetMomentum(newMom);
1432 G4ErrorTrajErr covG4e;
1434 g4eState.SetError(covG4e);
1435 extState.
chi2 += intersection.
chi2;
1451 double phi = intersection.
position.phi();
1454 if (phi > TWOPI - PI_8)
1456 int sector = (int)((phi + PI_8) / M_PI_4);
1469 intersection.
layer = layer;
1470 intersection.
sector = sector;
1486 double oldZ = std::fabs(oldPosition.z() -
m_OffsetZ);
1491 for (
int layer = extState.
firstEndcapLayer; layer <= outermostLayer; ++layer) {
1500 intersection.
layer = layer;
1502 isForward ? 2 : 1, intersection.
position) - 1;
1512 G4ThreeVector extPos0(intersection.
position);
1513 double diffBestMagSq = 1.0E60;
1515 int matchingLayer = intersection.
layer + 1;
1519 if (hit->getLayer() != matchingLayer)
1521 if (hit->isOutOfTime())
1525 G4ThreeVector diff(hit->getGlobalPositionX() - intersection.
position.x(),
1526 hit->getGlobalPositionY() - intersection.
position.y(),
1527 hit->getGlobalPositionZ() - intersection.
position.z());
1528 double dn = diff * n;
1529 if (std::fabs(dn) < 2.0) {
1532 if (diff.mag2() < diffBestMagSq) {
1533 diffBestMagSq = diff.mag2();
1539 if (std::fabs(dn) > 50.0)
1541 int sector = hit->getSector() - 1;
1542 int dSector = abs(intersection.
sector - sector);
1551 dn = diff * nHit + dn2;
1552 if (std::fabs(dn) > 1.0)
1555 G4ThreeVector extDir(intersection.
momentum.unit());
1556 double extDirA = extDir * nHit;
1557 if (std::fabs(extDirA) < 1.0E-6)
1559 G4ThreeVector projection = extDir * (dn2 / extDirA);
1560 if (projection.mag() > 15.0)
1562 diff += projection - nHit * dn;
1563 if (diff.mag2() < diffBestMagSq) {
1564 diffBestMagSq = diff.mag2();
1566 extPos0 = intersection.
position - projection;
1573 intersection.
isForward = (hit->getSection() == 1);
1574 intersection.
sector = hit->getSector() - 1;
1575 intersection.
time = hit->getTime();
1578 int nStrips = hit->getPhiStripMax() - hit->getPhiStripMin() + 1;
1579 double dn = nStrips - 1.5;
1580 double factor = std::pow((0.9 + 0.4 * dn * dn), 1.5) * 0.60;
1582 nStrips = hit->getZStripMax() - hit->getZStripMin() + 1;
1584 factor = std::pow((0.9 + 0.4 * dn * dn), 1.5) * 0.55;
1587 G4ThreeVector hitPos(hit->getGlobalPositionX(), hit->getGlobalPositionY(), hit->getGlobalPositionZ());
1589 if (intersection.
chi2 >= 0.0) {
1590 intersection.
hit = bestHit;
1591 hit->isOnTrack(
true);
1592 if (track !=
nullptr) {
1593 track->addRelationTo(hit);
1601 return intersection.
chi2 >= 0.0;
1607 double diffBestMagSq = 1.0E60;
1609 int matchingLayer = intersection.
layer + 1;
1610 int matchingEndcap = (intersection.
isForward ? 2 : 1);
1611 G4ThreeVector n(0.0, 0.0, (intersection.
isForward ? 1.0 : -1.0));
1614 if (hit->getLayer() != matchingLayer)
1616 if (hit->getSection() != matchingEndcap)
1622 G4ThreeVector diff(hit->getPositionX() - intersection.
position.x(),
1623 hit->getPositionY() - intersection.
position.y(),
1624 hit->getPositionZ() - intersection.
position.z());
1625 double dn = diff * n;
1626 if (std::fabs(dn) > 2.0)
1629 if (diff.mag2() < diffBestMagSq) {
1630 diffBestMagSq = diff.mag2();
1637 intersection.
hit = bestHit;
1638 intersection.
isForward = (hit->getSection() == 2);
1639 intersection.
sector = hit->getSector() - 1;
1640 intersection.
time = hit->getTime();
1642 G4ThreeVector hitPos(hit->getPositionX(), hit->getPositionY(), hit->getPositionZ());
1644 if (intersection.
chi2 >= 0.0) {
1647 if (track !=
nullptr) {
1648 track->addRelationTo(hit);
1658 return intersection.
chi2 >= 0.0;
1663 const G4ThreeVector& hitPos,
const G4ThreeVector& extPos0)
1683 G4ThreeVector extPos(extPos0);
1684 G4ThreeVector extMom(intersection.
momentum);
1685 G4ThreeVector extDir(extMom.unit());
1686 G4ThreeVector diffPos(hitPos - extPos);
1687 G4ErrorSymMatrix extCov(intersection.
covariance);
1689 G4ErrorMatrix extPar(6, 1);
1690 extPar[0][0] = extPos.x();
1691 extPar[1][0] = extPos.y();
1692 extPar[2][0] = extPos.z();
1693 extPar[3][0] = extMom.x();
1694 extPar[4][0] = extMom.y();
1695 extPar[5][0] = extMom.z();
1702 nC = G4ThreeVector(0.0, 0.0, 1.0);
1704 double out = (intersection.
isForward ? 1.0 : -1.0);
1705 nA = G4ThreeVector(0.0, 0.0, out);
1706 nB = G4ThreeVector(out, 0.0, 0.0);
1707 nC = G4ThreeVector(0.0, out, 0.0);
1710 double extDirA = extDir * nA;
1711 if (std::fabs(extDirA) < 1.0E-6)
1713 double extDirBA = extDir * nB / extDirA;
1714 double extDirCA = extDir * nC / extDirA;
1717 G4ThreeVector move = extDir * ((diffPos * nA) / extDirA);
1722 G4ErrorMatrix jacobian(2, 6);
1723 jacobian[0][0] = nB.x() - nA.x() * extDirBA;
1724 jacobian[0][1] = nB.y() - nA.y() * extDirBA;
1725 jacobian[0][2] = nB.z() - nA.z() * extDirBA;
1726 jacobian[1][0] = nC.x() - nA.x() * extDirCA;
1727 jacobian[1][1] = nC.y() - nA.y() * extDirCA;
1728 jacobian[1][2] = nC.z() - nA.z() * extDirCA;
1730 G4ErrorMatrix residual(2, 1);
1731 residual[0][0] = diffPos.x() * jacobian[0][0] + diffPos.y() * jacobian[0][1] + diffPos.z() * jacobian[0][2];
1732 residual[1][0] = diffPos.x() * jacobian[1][0] + diffPos.y() * jacobian[1][1] + diffPos.z() * jacobian[1][2];
1734 G4ErrorSymMatrix hitCov(2, 0);
1735 hitCov[0][0] = localVariance[0];
1736 hitCov[1][1] = localVariance[1];
1739 hitCov[0][0] *= 10.0;
1740 hitCov[1][1] *= 10.0;
1744 G4ErrorSymMatrix correction(extCov.similarity(jacobian) + hitCov);
1751 correction.invert(fail);
1757 intersection.
chi2 = (correction.similarityT(residual))[0][0];
1759 G4ErrorMatrix gain((extCov * jacobian.T()) * correction);
1760 G4ErrorSymMatrix HRH(correction.similarityT(jacobian));
1761 extCov -= HRH.similarity(extCov);
1762 extPar += gain * residual;
1763 extPos.set(extPar[0][0], extPar[1][0], extPar[2][0]);
1764 extMom.set(extPar[3][0], extPar[4][0], extPar[5][0]);
1766 correction = hitCov - extCov.similarity(jacobian);
1767 correction.invert(fail);
1770 diffPos = hitPos - extPos;
1771 residual[0][0] = diffPos.x() * jacobian[0][0] + diffPos.y() * jacobian[0][1] + diffPos.z() * jacobian[0][2];
1772 residual[1][0] = diffPos.x() * jacobian[1][0] + diffPos.y() * jacobian[1][1] + diffPos.z() * jacobian[1][2];
1773 intersection.
chi2 = (correction.similarityT(residual))[0][0];
1780 intersection.
position = extPos + extDir * (((intersection.
position - extPos) * nA) / extDirA);
1806 if (outcome != MuidElementNumbers::c_NotReached) {
1808 int charge = klmMuidLikelihood->
getCharge();
1810 std::map<int, double> mapPdgPDF;
1811 for (
int pdg : signedPdgVector) {
1814 B2FATAL(
"Something went wrong: PDF for PDG code " << pdg <<
" not found!");
1815 double pdf = (search->second)->getPDF(klmMuidLikelihood);
1817 mapPdgPDF.insert(std::pair<int, double>(std::abs(pdg), pdf));
1819 if (denom < 1.0E-20)
1822 for (
auto const& [pdg, pdf] : mapPdgPDF) {
1823 klmMuidLikelihood->
setPDFValue(pdf, std::abs(pdg));
1825 klmMuidLikelihood->
setLogL(std::log(pdf), std::abs(pdg));
static KLMChannelNumber channelNumber(int section, int sector, int layer, int plane, int strip)
Get channel number.
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
@ c_ForwardSection
Forward.
@ c_BackwardSection
Backward.
static constexpr int getMaximalSectorNumber()
Get maximal sector number (1-based).
static KLMModuleNumber moduleNumber(int section, int sector, int layer, bool fatalError=true)
Get module number.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number (0-based).
Store one BKLM strip hit as a ROOT object.
static void setMaximalStrip(int &module, int strip)
Set maximal strip number.
static const ChargedStable muon
muon particle
EDetector
Enum for identifying the detector components (detector and subdetector).
static const ChargedStable electron
electron particle
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
int ECLVolumeToCellID(const G4VTouchable *)
Get Cell Id (LEP: new way)
This dataobject is used only for EKLM alignment.
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
static const EKLMElementNumbers & Instance()
Instantiation.
int stripNumber(int section, int layer, int sector, int plane, int strip) const
Get strip number.
static constexpr int getMaximalLayerNumber()
Get maximal layer number.
@ c_ForwardSection
Forward.
@ c_BackwardSection
Backward.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
double getOuterR() const
Get outer radius.
double getInnerR() const
Get inner radius.
double getZ() const
Get Z coordinate.
double getLength() const
Get length.
double getWidth() const
Get width.
const StripGeometry * getStripGeometry() const
Get strip geometry data.
double getLayerShiftZ() const
Get Z distance between two layers.
const ElementPosition * getLayerPosition() const
Get position data for layers.
const ElementPosition * getSectionPosition() const
Get position data for sections.
Class for 2d hits handling.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Store one Ext hit as a ROOT object.
@ c_IndexLevelLayer
Layer.
ChannelStatus
Channel status.
@ c_Dead
Dead channel (no signal).
@ c_Unknown
Unknown status (no data).
KLMChannelNumber channelNumberBKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for BKLM.
static const KLMElementNumbers & Instance()
Instantiation.
KLMChannelNumber channelNumberEKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for EKLM.
Store one muon-identification hit in the KLM as a ROOT object.
Class to store the likelihoods from KLM with additional informations related to the extrapolation.
void setExtLayer(int layer)
Set the outermost EKLM layer crossed in the extrapolation.
void setExtLayerPattern(unsigned int pattern)
Set the pattern of the layers crossed in the extrapolation.
void setLogL(double logL, int pdg)
Set the log-likelihood.
void setExtEKLMEfficiencyValue(int layer, float efficiency)
Set the efficiency of a given EKLM layer.
void setExtBKLMEfficiencyValue(int layer, float efficiency)
Set the efficiency of a given BKLM layer.
void setChiSquared(double chiSquared)
Set the chi-squared of the extrapolation.
void setHitLayerPattern(unsigned int pattern)
Set the pattern of the layers actually crossed by the track.
void setBarrelExtLayer(int layer)
Set the outermost BKLM layer crossed in the extrapolation.
void setIsForward(bool isForward)
Set if this extrapolation is in forward or backward B/EKLM.
void setJunkPDFValue(bool flag)
Set the junk flag (1 if junk, 0 if not).
void setPDGCode(int pdg)
Set the PDG code of the particle hypothesis used during the extrapolation.
void setDegreesOfFreedom(int dof)
Set the number of degrees of freedom (= 2 times the number of KLM hits) for the chi-square computatio...
int getCharge() const
Get the charge of the particle hypothesis used during the extrapolation.
void setPDFValue(double pdfValue, int pdg)
Set the normalized PDF.
void setEndcapExtLayer(int layer)
Set the outermost EKLM layer crossed in the extrapolation.
void setBarrelHitLayer(int layer)
Set the outermost BKLM layer actually crossed by the track.
void setOutcome(unsigned int outcome)
Set the outcome of this extrapolation.
void setHitLayer(int layer)
Set the outermost KLM layer actually crossed by the track.
void setEndcapHitLayer(int layer)
Set the outermost EKLM layer actually crossed by the track.
Build the Muid likelihoods starting from the hit pattern and the transverse scattering in the KLM.
static unsigned int calculateExtrapolationOutcome(bool isForward, bool escaped, int lastBarrelLayer, int lastEndcapLayer)
Calculate the track extrapolation outcome.
static std::vector< int > getPDGVector()
Get a vector with all the hypothesis PDG codes used for Muid.
This is the Reconstruction Event-Data Model Track.
bool addBKLMHit(const UsedBKLMHit *bklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a bklm hit with the given information to the reco track.
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
bool addEKLMHit(const UsedEKLMHit *eklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds an eklm hit with the given information to the reco track.
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromLastHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the last hit in a fit useful for extrapolation of measuremen...
static void registerRequiredRelations(StoreArray< RecoTrack > &recoTracks, std::string const &pxdHitsStoreArrayName="", std::string const &svdHitsStoreArrayName="", std::string const &cdcHitsStoreArrayName="", std::string const &bklmHitsStoreArrayName="", std::string const &eklmHitsStoreArrayName="", std::string const &recoHitInformationStoreArrayName="")
Convenience method which registers all relations required to fully use a RecoTrack.
genfit::AbsTrackRep * getCardinalRepresentation() const
Get a pointer to the cardinal track representation. You are not allowed to modify or delete it!
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromFirstHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the first hit in a fit useful for extrapolation of measureme...
unsigned int getNumberOfTotalHits() const
Return the number of cdc + svd + pxd + bklm + eklm hits.
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
Defines a closed cylinder for the geant4e "target", the surface that encloses the volume within which...
virtual G4double GetDistanceFromPoint(const G4ThreeVector &point, const G4ThreeVector &direc) const
Get the distance from a point to the cylinder along direc.
It is the main interface for the user to define the setup and start the propagation.
static ExtManager * GetManager()
Get pointer to the instance of this singleton class (create if needed)
void RunTermination()
Terminate a run and set state to G4ErrorState_Init.
void EventTermination(G4ErrorMode)
Terminate an event and set state to G4ErrorState_Init.
void InitTrackPropagation(G4ErrorMode)
Initialize for propagation of a track and set state to G4ErrorState_Propagating.
void Initialize(const char[], const std::string &, double, double, bool, int, const std::vector< std::string > &)
Initialize Geant4 and Geant4e.
G4int PropagateOneStep(G4ErrorTrajState *currentTS, G4ErrorMode mode=G4ErrorMode_PropForwards)
Propagate a track by one step.
G4ErrorPropagator * GetPropagator() const
Get the propagator.
Store one Track-KLMCluster separation as a ROOT object.
Class that bundles various TrackFitResults.
static const double T
[tesla]
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
double getScintHalfWidth(void) const
Get the height of the entire volume of a scintillator strip (including TiO2 coating)
double getGap1InnerRadius(void) const
Get the radius of the inner tangent circle of gap 0 (innermost)
double getOuterRadius(void) const
Get the radius of the inscribed circle of the outer polygon.
double getOffsetZ(void) const
Get the global shift along a of the entire BKLM.
int getNSector(void) const
Get the number of sectors of the BKLM.
double getHalfLength(void) const
Get the half-length along z of the BKLM.
double getActiveMiddleRadius(int section, int sector, int layer) const
Get the radial midpoint of the detector module's active volume of specified layer.
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Class to store variables with their name which were sent to the logging service.
Abstract base class for a track representation.
double getPDGCharge() const
Get the charge of the particle of the pdg code.
int getPDG() const
Get the pdg code.
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
#StateOnPlane with additional covariance matrix.
B2Vector3< double > B2Vector3D
typedef for common usage with double
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
ExtHitStatus
Define state of extrapolation for each recorded hit.
@ VOLTYPE_ARICH2
ARICH Img plate.
@ VOLTYPE_TOP2
TOP quartz.
@ VOLTYPE_ARICH3
ARICH HAPD window.
@ VOLTYPE_BKLM2
BKLM scintillator.
@ VOLTYPE_ARICH1
ARICH aerogel.
@ VOLTYPE_TOP1
TOP container.
Abstract base class for different kinds of events.
Data structure to define extrapolation state.
int lastEndcapExtLayer
MUID: outermost endcap layer crossed by the extrapolated track.
double chi2
MUID: accumulated chi-squared of all in-plane transverse deviations between extrapolation and matchin...
G4ThreeVector directionAtIP
MUID: initial direction of track, used for KLID.
int nPoint
MUID: accumulated number of points with matching 2D hits.
int lastEndcapHitLayer
MUID: outermost endcap layer with a matching hit.
bool isCosmic
True for back-propagation of a cosmic ray.
int firstEndcapLayer
MUID: outermost barrel layer encountered by the extrapolated track in the prior steps.
int extLayerPattern
MUID: accumulated bit pattern of layers crossed by the extrapolated track.
int pdgCode
Particle hypothesis that is being extrapolated.
double length
Length from start of extrapolation (rad lengths), updated during extrapolation.
bool escaped
MUID: flag to indicate that the extrapolated track escaped from the KLM.
int firstBarrelLayer
MUID: outermost barrel layer encountered by the extrapolated track in the prior steps.
int lastBarrelHitLayer
MUID: outermost barrel layer with a matching hit.
int hitLayerPattern
MUID: accumulated bit pattern of layers with matching hits.
int lastBarrelExtLayer
MUID: outermost barrel layer crossed by the extrapolated track.
double tof
Time of flight from IP (ns), updated during extrapolation.
const Track * track
Pointer to the reconstructed track.
intersection of muid-extrapolated track with a KLM layer
int sector
sector number (0..7 for barrel, 0..3 for endcap) of this point
G4ThreeVector momentum
extrapolated-track momentum (GeV/c) at this intersection
double chi2
chi-squared value of transverse deviation between extrapolated and measured hit positions
bool inBarrel
flag to indicate if this point is in the barrel (true) or endcap (false)
int hit
index in {B,E}KLMHit2ds of matching hit
G4ThreeVector positionAtHitPlane
extrapolated-track position (cm) projected to the 2D hit's midplane
G4ThreeVector position
extrapolated-track global position (cm) of this intersection
G4ErrorSymMatrix covariance
extrapolated-track phase-space covariance matrix at this intersection
bool isForward
flag to indicate if this point is in the forward (true) or backward (false) end
double time
time (ns) of matching BKLMHit2d
int layer
layer number (0..14 for barrel, 0..13 for endcap) of this point