File indexing completed on 2025-02-23 09:22:33
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030 #include "Par02FastSimModelTracker.hh"
0031
0032 #include "Par02EventInformation.hh"
0033 #include "Par02Output.hh"
0034 #include "Par02PrimaryParticleInformation.hh"
0035 #include "Par02Smearer.hh"
0036
0037 #include "G4AnalysisManager.hh"
0038 #include "G4Electron.hh"
0039 #include "G4Event.hh"
0040 #include "G4FieldTrack.hh"
0041 #include "G4FieldTrackUpdator.hh"
0042 #include "G4Gamma.hh"
0043 #include "G4PathFinder.hh"
0044 #include "G4Positron.hh"
0045 #include "G4RunManager.hh"
0046 #include "G4SystemOfUnits.hh"
0047 #include "G4Track.hh"
0048 #include "Randomize.hh"
0049
0050
0051
0052 Par02FastSimModelTracker::Par02FastSimModelTracker(
0053 G4String aModelName, G4Region* aEnvelope, Par02DetectorParametrisation::Parametrisation aType)
0054 : G4VFastSimulationModel(aModelName, aEnvelope),
0055 fCalculateParametrisation(),
0056 fParametrisation(aType)
0057 {}
0058
0059
0060
0061 Par02FastSimModelTracker::Par02FastSimModelTracker(G4String aModelName, G4Region* aEnvelope)
0062 : G4VFastSimulationModel(aModelName, aEnvelope),
0063 fCalculateParametrisation(),
0064 fParametrisation(Par02DetectorParametrisation::eCMS)
0065 {}
0066
0067
0068
0069 Par02FastSimModelTracker::Par02FastSimModelTracker(G4String aModelName)
0070 : G4VFastSimulationModel(aModelName),
0071 fCalculateParametrisation(),
0072 fParametrisation(Par02DetectorParametrisation::eCMS)
0073 {}
0074
0075
0076
0077 Par02FastSimModelTracker::~Par02FastSimModelTracker() = default;
0078
0079
0080
0081 G4bool Par02FastSimModelTracker::IsApplicable(const G4ParticleDefinition& aParticleType)
0082 {
0083 return aParticleType.GetPDGCharge() != 0;
0084 }
0085
0086
0087
0088 G4bool Par02FastSimModelTracker::ModelTrigger(const G4FastTrack& )
0089 {
0090 return true;
0091 }
0092
0093
0094
0095 void Par02FastSimModelTracker::DoIt(const G4FastTrack& aFastTrack, G4FastStep& aFastStep)
0096 {
0097 G4cout << " ________Tracker model triggered _________" << G4endl;
0098
0099
0100
0101
0102 G4Track track = *aFastTrack.GetPrimaryTrack();
0103 G4FieldTrack aFieldTrack('0');
0104 G4FieldTrackUpdator::Update(&aFieldTrack, &track);
0105
0106 G4double retSafety = -1.0;
0107 ELimited retStepLimited;
0108 G4FieldTrack endTrack('a');
0109 G4double currentMinimumStep = 10.0 * m;
0110
0111 G4PathFinder* fPathFinder = G4PathFinder::GetInstance();
0112
0113 fPathFinder->ComputeStep(aFieldTrack, currentMinimumStep, 0,
0114 aFastTrack.GetPrimaryTrack()->GetCurrentStepNumber(), retSafety,
0115 retStepLimited, endTrack, aFastTrack.GetPrimaryTrack()->GetVolume());
0116
0117
0118
0119 aFastStep.ProposePrimaryTrackFinalPosition(endTrack.GetPosition());
0120
0121
0122 G4ThreeVector Porg = aFastTrack.GetPrimaryTrack()->GetMomentum();
0123 if (!aFastTrack.GetPrimaryTrack()->GetParentID()) {
0124 auto info = (Par02EventInformation*)G4EventManager::GetEventManager()->GetUserInformation();
0125 if (info->GetDoSmearing()) {
0126
0127 G4double res = fCalculateParametrisation->GetResolution(
0128 Par02DetectorParametrisation::eTRACKER, fParametrisation, Porg.mag());
0129 G4double eff = fCalculateParametrisation->GetEfficiency(
0130 Par02DetectorParametrisation::eTRACKER, fParametrisation, Porg.mag());
0131 G4ThreeVector Psm;
0132 Psm = Par02Smearer::Instance()->SmearMomentum(aFastTrack.GetPrimaryTrack(), res);
0133 Par02Output::Instance()->FillHistogram(0, ((Psm.mag() / MeV) / (Porg.mag() / MeV)));
0134
0135 ((Par02PrimaryParticleInformation*)(const_cast<G4PrimaryParticle*>(
0136 aFastTrack.GetPrimaryTrack()
0137 ->GetDynamicParticle()
0138 ->GetPrimaryParticle())
0139 ->GetUserInformation()))
0140 ->SetTrackerMomentum(Psm);
0141 ((Par02PrimaryParticleInformation*)(const_cast<G4PrimaryParticle*>(
0142 aFastTrack.GetPrimaryTrack()
0143 ->GetDynamicParticle()
0144 ->GetPrimaryParticle())
0145 ->GetUserInformation()))
0146 ->SetTrackerResolution(res);
0147 ((Par02PrimaryParticleInformation*)(const_cast<G4PrimaryParticle*>(
0148 aFastTrack.GetPrimaryTrack()
0149 ->GetDynamicParticle()
0150 ->GetPrimaryParticle())
0151 ->GetUserInformation()))
0152 ->SetTrackerEfficiency(eff);
0153 }
0154 else {
0155
0156 ((Par02PrimaryParticleInformation*)(const_cast<G4PrimaryParticle*>(
0157 aFastTrack.GetPrimaryTrack()
0158 ->GetDynamicParticle()
0159 ->GetPrimaryParticle())
0160 ->GetUserInformation()))
0161 ->SetTrackerMomentum(Porg);
0162 }
0163 }
0164 }
0165
0166