Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-30 09:17:35

0001 //==========================================================================
0002 //  AIDA Detector description implementation 
0003 //--------------------------------------------------------------------------
0004 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
0005 // All rights reserved.
0006 //
0007 // For the licensing terms see $DD4hepINSTALL/LICENSE.
0008 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
0009 //
0010 // Author     : M.Frank
0011 //
0012 //==========================================================================
0013 
0014 // Framework include files
0015 #include <DDG4/Geant4UserLimits.h>
0016 #include <DDG4/Geant4Particle.h>
0017 #include <DD4hep/InstanceCount.h>
0018 #include <DD4hep/DD4hepUnits.h>
0019 #include <DD4hep/Primitives.h>
0020 #include <DD4hep/Printout.h>
0021 
0022 // Geant 4 include files
0023 #include <G4Track.hh>
0024 #include <G4ParticleDefinition.hh>
0025 #include <CLHEP/Units/SystemOfUnits.h>
0026 
0027 using namespace dd4hep::sim;
0028 
0029 namespace  {
0030   bool user_limit_debug = false;
0031 }
0032 
0033 /// Allow for debugging user limits (very verbose)
0034 bool Geant4UserLimits::enable_debug(bool value)   {
0035   bool tmp = user_limit_debug;
0036   user_limit_debug = value;
0037   return tmp;
0038 }
0039 
0040 /// Access value according to track
0041 double Geant4UserLimits::Handler::value(const G4Track& track) const    {
0042   auto* def = track.GetParticleDefinition();
0043   if ( !particleLimits.empty() )  {
0044     auto i = particleLimits.find(track.GetDefinition());
0045     if ( i != particleLimits.end() )  {
0046       if ( user_limit_debug )   {
0047         dd4hep::printout(dd4hep::INFO,"Geant4UserLimits", "Apply explicit limit %f to track: %s",
0048                          def->GetParticleName().c_str());
0049       }
0050       return (*i).second;
0051     }
0052   }
0053   if ( user_limit_debug )   {
0054     dd4hep::printout(dd4hep::INFO,"Geant4UserLimits", "Apply default limit %f to track: %s",
0055                      def->GetParticleName().c_str());
0056   }
0057   return defaultValue;
0058 }
0059 
0060 /// Set the handler value(s)
0061 void Geant4UserLimits::Handler::set(const std::string& particles, double val)   {
0062   if ( particles == "*" || particles == ".(.*)" )   {
0063     defaultValue = val;
0064     return;
0065   }
0066   auto defs = Geant4ParticleHandle::g4DefinitionsRegEx(particles);
0067   for( auto* d : defs )  {
0068     particleLimits[d] = val;
0069   }
0070 }
0071 
0072 /// Initializing Constructor
0073 Geant4UserLimits::Geant4UserLimits(LimitSet limitset)
0074   : G4UserLimits(limitset.name()), limits(limitset)
0075 {
0076   this->update(limitset);
0077   InstanceCount::increment(this);
0078 }
0079 
0080 /// Standard destructor
0081 Geant4UserLimits::~Geant4UserLimits()  {
0082   InstanceCount::decrement(this);
0083 }
0084 
0085 /// Update the object
0086 void Geant4UserLimits::update(LimitSet limitset)    {
0087   const auto& lim = limitset.limits();
0088   /// Set defaults
0089   maxStepLength.defaultValue  = fMaxStep;
0090   maxTrackLength.defaultValue = fMaxTrack;
0091   maxTime.defaultValue        = fMaxTime;
0092   minEKine.defaultValue       = fMinEkine;
0093   minRange.defaultValue       = fMinRange;
0094   /// Overwrite with values if present:
0095   for(const Limit& l : lim)   {
0096     if (l.name == "step_length_max")
0097       maxStepLength.set(l.particles, l.value*CLHEP::mm/dd4hep::mm);
0098     else if (l.name == "track_length_max")
0099       maxTrackLength.set(l.particles, l.value*CLHEP::mm/dd4hep::mm);
0100     else if (l.name == "time_max")
0101       maxTime.set(l.particles, l.value*CLHEP::ns/dd4hep::ns);
0102     else if (l.name == "ekin_min")
0103       minEKine.set(l.particles, l.value*CLHEP::MeV/dd4hep::MeV);
0104     else if (l.name == "range_min")
0105       minRange.set(l.particles, l.value);
0106     else
0107       except("Geant4UserLimits", "Unknown Geant4 user limit: %s ", l.toString().c_str());
0108   }
0109 }
0110 
0111 /// Setters may not be called!
0112 void Geant4UserLimits::SetMaxAllowedStep(G4double /* ustepMax */)  {
0113   dd4hep::notImplemented(std::string(__PRETTY_FUNCTION__)+" May not be called!");
0114 }
0115 
0116 void Geant4UserLimits::SetUserMaxTrackLength(G4double /* utrakMax */)  {
0117   dd4hep::notImplemented(std::string(__PRETTY_FUNCTION__)+" May not be called!");
0118 }
0119 
0120 void Geant4UserLimits::SetUserMaxTime(G4double /* utimeMax */)  {
0121   dd4hep::notImplemented(std::string(__PRETTY_FUNCTION__)+" May not be called!");
0122 }
0123 
0124 void Geant4UserLimits::SetUserMinEkine(G4double /* uekinMin */)  {
0125   dd4hep::notImplemented(std::string(__PRETTY_FUNCTION__)+" May not be called!");
0126 }
0127 
0128 void Geant4UserLimits::SetUserMinRange(G4double /* urangMin */)  {
0129   dd4hep::notImplemented(std::string(__PRETTY_FUNCTION__)+" May not be called!");
0130 }
0131