Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:21:27

0001 //
0002 // =======================================================================
0003 //
0004 // Class MonopoleTransportation
0005 //
0006 // Created:  3 May 2010, J. Apostolakis, B. Bozsogi
0007 //                       G4MonopoleTransportation class for
0008 //                       Geant4 extended example "monopole"
0009 //
0010 // Adopted for CMSSW by V.Ivanchenko 30 April 2018
0011 // from Geant4 global tag geant4-10-04-ref-03
0012 //
0013 // =======================================================================
0014 //
0015 // Class description:
0016 //
0017 // G4MonopoleTransportation is a process responsible for the transportation of
0018 // magnetic monopoles, i.e. the geometrical propagation encountering the
0019 // geometrical sub-volumes of the detectors.
0020 // It is also tasked with part of updating the "safety".
0021 // For monopoles, uses a different equation of motion and ignores energy
0022 // conservation.
0023 //
0024 
0025 #ifndef SimG4Core_PhysicsLists_MonopoleTransportation_h
0026 #define SimG4Core_PhysicsLists_MonopoleTransportation_h 1
0027 
0028 #include "G4VProcess.hh"
0029 #include "G4FieldManager.hh"
0030 
0031 #include "G4Navigator.hh"
0032 #include "G4TransportationManager.hh"
0033 #include "G4PropagatorInField.hh"
0034 #include "G4Track.hh"
0035 #include "G4Step.hh"
0036 #include "G4ParticleChangeForTransport.hh"
0037 #include <CLHEP/Units/SystemOfUnits.h>
0038 
0039 #include <memory>
0040 
0041 class G4SafetyHelper;
0042 class Monopole;
0043 class CMSFieldManager;
0044 
0045 class MonopoleTransportation : public G4VProcess {
0046 public:
0047   MonopoleTransportation(const Monopole* p, G4int verbosityLevel = 1);
0048   ~MonopoleTransportation() override;
0049 
0050   G4double AlongStepGetPhysicalInteractionLength(const G4Track& track,
0051                                                  G4double previousStepSize,
0052                                                  G4double currentMinimumStep,
0053                                                  G4double& currentSafety,
0054                                                  G4GPILSelection* selection) override;
0055 
0056   G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step& stepData) override;
0057 
0058   G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step& stepData) override;
0059   // Responsible for the relocation.
0060 
0061   G4double PostStepGetPhysicalInteractionLength(const G4Track&,
0062                                                 G4double previousStepSize,
0063                                                 G4ForceCondition* pForceCond) override;
0064   // Forces the PostStepDoIt action to be called,
0065   // but does not limit the step.
0066 
0067   G4PropagatorInField* GetPropagatorInField();
0068   void SetPropagatorInField(G4PropagatorInField* pFieldPropagator);
0069   // Access/set the assistant class that Propagate in a Field.
0070 
0071   inline G4double GetThresholdWarningEnergy() const;
0072   inline G4double GetThresholdImportantEnergy() const;
0073   inline G4int GetThresholdTrials() const;
0074 
0075   inline void SetThresholdWarningEnergy(G4double newEnWarn);
0076   inline void SetThresholdImportantEnergy(G4double newEnImp);
0077   inline void SetThresholdTrials(G4int newMaxTrials);
0078 
0079   // Get/Set parameters for killing loopers:
0080   //   Above 'important' energy a 'looping' particle in field will
0081   //   *NOT* be abandoned, except after fThresholdTrials attempts.
0082   // Below Warning energy, no verbosity for looping particles is issued
0083 
0084   inline G4double GetMaxEnergyKilled() const;
0085   inline G4double GetSumEnergyKilled() const;
0086   void ResetKilledStatistics(G4int report = 1);
0087   // Statistics for tracks killed (currently due to looping in field)
0088 
0089   inline void EnableShortStepOptimisation(G4bool optimise = true);
0090   // Whether short steps < safety will avoid to call Navigator (if field=0)
0091 
0092   G4double AtRestGetPhysicalInteractionLength(const G4Track&, G4ForceCondition*) override;
0093 
0094   G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override;
0095   // No operation in  AtRestDoIt.
0096 
0097   void StartTracking(G4Track* aTrack) override;
0098   // Reset state for new (potentially resumed) track
0099 
0100 protected:
0101   G4bool DoesGlobalFieldExist();
0102   // Checks whether a field exists for the "global" field manager.
0103 
0104 private:
0105   const Monopole* fParticleDef;
0106 
0107   CMSFieldManager* fieldMgrCMS;
0108 
0109   G4Navigator* fLinearNavigator;
0110   G4PropagatorInField* fFieldPropagator;
0111   // The Propagators used to transport the particle
0112 
0113   G4ThreeVector fTransportEndPosition;
0114   G4ThreeVector fTransportEndMomentumDir;
0115   G4double fTransportEndKineticEnergy;
0116   G4ThreeVector fTransportEndSpin;
0117   G4bool fMomentumChanged;
0118 
0119   G4bool fEndGlobalTimeComputed;
0120   G4double fCandidateEndGlobalTime;
0121   // The particle's state after this Step, Store for DoIt
0122 
0123   G4bool fParticleIsLooping;
0124 
0125   G4TouchableHandle fCurrentTouchableHandle;
0126 
0127   G4bool fGeometryLimitedStep;
0128   // Flag to determine whether a boundary was reached.
0129 
0130   G4ThreeVector fPreviousSftOrigin;
0131   G4double fPreviousSafety;
0132   // Remember last safety origin & value.
0133 
0134   G4ParticleChangeForTransport fParticleChange;
0135   // New ParticleChange
0136 
0137   G4double endpointDistance;
0138 
0139   // Thresholds for looping particles:
0140   //
0141   G4double fThreshold_Warning_Energy;    //  Warn above this energy
0142   G4double fThreshold_Important_Energy;  //  Hesitate above this
0143   G4int fThresholdTrials;                //    for this no of trials
0144   // Above 'important' energy a 'looping' particle in field will
0145   //   *NOT* be abandoned, except after fThresholdTrials attempts.
0146 
0147   // Counter for steps in which particle reports 'looping',
0148   //   if it is above 'Important' Energy
0149   G4int fNoLooperTrials;
0150   // Statistics for tracks abandoned
0151   G4double fSumEnergyKilled;
0152   G4double fMaxEnergyKilled;
0153 
0154   // Whether to avoid calling G4Navigator for short step ( < safety)
0155   //   If using it, the safety estimate for endpoint will likely be smaller.
0156   G4bool fShortStepOptimisation;
0157 
0158   G4SafetyHelper* fpSafetyHelper;  // To pass it the safety value obtained
0159 };
0160 
0161 inline void MonopoleTransportation::SetPropagatorInField(G4PropagatorInField* pFieldPropagator) {
0162   fFieldPropagator = pFieldPropagator;
0163 }
0164 
0165 inline G4PropagatorInField* MonopoleTransportation::GetPropagatorInField() { return fFieldPropagator; }
0166 
0167 inline G4bool MonopoleTransportation::DoesGlobalFieldExist() {
0168   G4TransportationManager* transportMgr = G4TransportationManager::GetTransportationManager();
0169   return transportMgr->GetFieldManager()->DoesFieldExist();
0170 }
0171 
0172 inline G4double MonopoleTransportation::GetThresholdWarningEnergy() const { return fThreshold_Warning_Energy; }
0173 
0174 inline G4double MonopoleTransportation::GetThresholdImportantEnergy() const { return fThreshold_Important_Energy; }
0175 
0176 inline G4int MonopoleTransportation::GetThresholdTrials() const { return fThresholdTrials; }
0177 
0178 inline void MonopoleTransportation::SetThresholdWarningEnergy(G4double newEnWarn) {
0179   fThreshold_Warning_Energy = newEnWarn;
0180 }
0181 
0182 inline void MonopoleTransportation::SetThresholdImportantEnergy(G4double newEnImp) {
0183   fThreshold_Important_Energy = newEnImp;
0184 }
0185 
0186 inline void MonopoleTransportation::SetThresholdTrials(G4int newMaxTrials) { fThresholdTrials = newMaxTrials; }
0187 
0188 // Get parameters for killing loopers:
0189 //   Above 'important' energy a 'looping' particle in field will
0190 //   *NOT* be abandoned, except after fThresholdTrials attempts.
0191 // Below Warning energy, no verbosity for looping particles is issued
0192 
0193 inline G4double MonopoleTransportation::GetMaxEnergyKilled() const { return fMaxEnergyKilled; }
0194 
0195 inline G4double MonopoleTransportation::GetSumEnergyKilled() const { return fSumEnergyKilled; }
0196 
0197 inline void MonopoleTransportation::EnableShortStepOptimisation(G4bool optimiseShortStep) {
0198   fShortStepOptimisation = optimiseShortStep;
0199 }
0200 
0201 #endif