Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-10 07:50:32

0001 #pragma once
0002 
0003 /**
0004 U4Process.h : Cherry Picking from cfg4/CProcessManager
0005 ================================================================
0006 
0007 Process identification relies on the processes using 
0008 their default names, for example from "g4-cls G4OpRayleigh"
0009 the default process name is "OpRayleigh". 
0010 
0011 **/
0012 
0013 #include <cassert>
0014 #include <cstring>
0015 #include <string>
0016 #include <sstream>
0017 
0018 #include "G4ProcessManager.hh"
0019 #include "G4VDiscreteProcess.hh"
0020 #include "G4OpticalPhoton.hh"
0021 #include "G4Track.hh"
0022 #include "G4Step.hh"
0023 
0024 enum {
0025    U4Process_Unknown        = 0, 
0026    U4Process_Transportation = 1, 
0027    U4Process_Scintillation  = 2,
0028    U4Process_OpAbsorption   = 3,
0029    U4Process_OpRayleigh     = 4,
0030    U4Process_OpBoundary     = 5
0031 }; 
0032 
0033 struct U4Process
0034 {
0035     static constexpr const char* Transportation_ = "Transportation" ; 
0036     static constexpr const char* Scintillation_ = "Scintillation" ; 
0037     static constexpr const char* OpAbsorption_ = "OpAbsorption" ; 
0038     static constexpr const char* OpRayleigh_ = "OpRayleigh" ; 
0039     static constexpr const char* OpBoundary_ = "OpBoundary" ; 
0040 
0041     template<typename T> static T* Get() ; 
0042 
0043     static const char* Name(const G4VProcess* proc); 
0044     static unsigned ProcessType(const char* name); 
0045     static bool IsKnownProcess(unsigned type); 
0046     static bool IsNormalProcess(unsigned type); 
0047     static bool IsPeculiarProcess(unsigned type); 
0048 
0049     static G4ProcessManager* GetManager(); 
0050     static std::string Desc(); 
0051 
0052     static void ClearNumberOfInteractionLengthLeft(const G4Track& aTrack, const G4Step& aStep); 
0053 
0054 };
0055 
0056 /**
0057 U4Process::Get
0058 ----------------
0059 
0060 Find a process of the OpticalPhoton process manager 
0061 with the template type. If none of the processes 
0062 are able to be dynamically cast to the template 
0063 type then nullptr is returned. Usage::
0064 
0065      G4OpRayleigh* p = U4Process::Get<G4OpRayleigh>() ; 
0066 
0067 
0068 **/
0069 
0070 template<typename T>
0071 inline T* U4Process::Get()
0072 {
0073     G4ProcessManager* mgr = GetManager(); 
0074     if(mgr == nullptr)
0075     {
0076         std::cerr << "U4Process::Get FAILED : MAYBE DO THIS LATER " << std::endl ; 
0077         return nullptr ; 
0078     }
0079     G4ProcessVector* procv = mgr->GetProcessList() ;
0080     G4int n = procv->entries() ;
0081 
0082     T* p = nullptr ; 
0083     int count = 0 ; 
0084 
0085     for(int i=0 ; i < n ; i++)
0086     {   
0087         G4VProcess* proc = (*procv)[i] ; 
0088         T* proc_ = dynamic_cast<T*>(proc) ; 
0089         if(proc_ != nullptr)
0090         {
0091             p = proc_ ; 
0092             count += 1 ; 
0093         }
0094     }
0095     assert( count == 0 || count == 1 ); 
0096     return p ; 
0097 }
0098 
0099 
0100 inline const char* U4Process::Name(const G4VProcess* proc)
0101 {
0102     const G4String& name = proc->GetProcessName() ;
0103     return name.c_str(); 
0104 }
0105 inline unsigned U4Process::ProcessType(const char* name)
0106 {
0107     unsigned type = U4Process_Unknown ; 
0108     if(strcmp(name, Transportation_) == 0)  type = U4Process_Transportation ; 
0109     if(strcmp(name, Scintillation_) == 0)   type = U4Process_Scintillation ; 
0110     if(strcmp(name, OpAbsorption_) == 0)    type = U4Process_OpAbsorption ; 
0111     if(strcmp(name, OpRayleigh_) == 0)      type = U4Process_OpRayleigh ; 
0112     if(strcmp(name, OpBoundary_) == 0)      type = U4Process_OpBoundary ; 
0113     return type ;  
0114 }
0115 
0116 inline bool U4Process::IsKnownProcess(unsigned type)
0117 { 
0118      return type != U4Process_Unknown ; 
0119 }
0120 inline bool U4Process::IsNormalProcess(unsigned type)
0121 {
0122     return type == U4Process_OpRayleigh || type == U4Process_OpAbsorption ; 
0123 } 
0124 inline bool U4Process::IsPeculiarProcess(unsigned type)
0125 {
0126     return type == U4Process_Transportation || type == U4Process_Scintillation || type == U4Process_OpBoundary ; 
0127 } 
0128 inline G4ProcessManager* U4Process::GetManager()
0129 {
0130     return G4OpticalPhoton::OpticalPhotonDefinition()->GetProcessManager() ; 
0131 }
0132 
0133 std::string U4Process::Desc()
0134 {
0135     G4ProcessManager* mgr = GetManager(); 
0136     assert(mgr); 
0137     std::stringstream ss ; 
0138     ss << "U4Process::Desc" ;
0139     G4ProcessVector* procv = mgr->GetProcessList() ;
0140     G4int n = procv->entries() ;
0141     ss << " n:[" << n << "]" << std::endl ;   
0142     for(int i=0 ; i < n ; i++)
0143     {   
0144         G4VProcess* proc = (*procv)[i] ; 
0145         const char* name = Name(proc);  
0146         unsigned type = ProcessType(name); 
0147         assert(IsKnownProcess(type)); 
0148         bool peculiar = IsPeculiarProcess(type); 
0149         G4double lenleft = proc->GetNumberOfInteractionLengthLeft() ; 
0150         bool unset = lenleft == -1. ; 
0151         bool expected = peculiar == unset ; 
0152         //assert( peculiar == unset );  
0153         // not obeyed when called from U4Recorder::PreUserTrackingAction_Optical but is from stepping
0154 
0155         ss << " (" << i << ")" 
0156            << " name " << name 
0157            << " type " << type
0158            << " lenleft " << lenleft
0159            << " unset " << unset
0160            << " peculiar " << peculiar 
0161            << " expected " << ( expected ? "YES" : "NO" )
0162            << std::endl
0163            ;   
0164     }   
0165     return ss.str();
0166 }
0167 
0168 
0169 /**
0170 U4Process::ClearNumberOfInteractionLengthLeft
0171 -----------------------------------------------
0172 
0173 This clears the interaction length left for "normal" optical processes : OpAbsorption and OpRayleigh 
0174 with no use of G4UniformRand.
0175 
0176 Formerly canonically called from tail of CManager::UserSteppingAction/CManager::prepareForNextStep
0177 when track is not terminating (ie not track status fStopAndKill so there is another step).
0178 
0179 The effect of this is to make the Geant4 random consumption more regular which 
0180 makes it easier to align with Opticks random consumption. 
0181 
0182 This provides a devious way to invoke the protected ClearNumberOfInteractionLengthLeft 
0183 via the public G4VDiscreteProcess::PostStepDoIt
0184 
0185 g4-;g4-cls G4VDiscreteProcess::
0186 
0187     112 G4VParticleChange* G4VDiscreteProcess::PostStepDoIt(
0188     113                             const G4Track& ,
0189     114                             const G4Step&
0190     115                             )
0191     116 { 
0192     117 //  clear NumberOfInteractionLengthLeft
0193     118     ClearNumberOfInteractionLengthLeft();
0194     119 
0195     120     return pParticleChange;
0196     121 }
0197 
0198 
0199 Noticed that returns from::
0200 
0201      InstrumentedG4OpBoundaryProcess::PostStepDoIt
0202      G4OpAbsorption::PostStepDoIt
0203      G4OpRayleigh::PostStepDoIt
0204 
0205 All do G4VDiscreteProcess::PostStepDoIt, eg::
0206 
0207      182 G4VParticleChange*
0208      183 InstrumentedG4OpBoundaryProcess::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
0209      184 {
0210      ...
0211      623         return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
0212      624 }
0213 
0214 
0215 Q: Why the artificial call when the PostStepDoIt returns do it anyhow ?
0216 A: Because only the winning process PostStepDoIt runs, and that changes 
0217    with the history : so you get a complicated consumption history 
0218    that would be difficult to align with.  
0219 
0220 **/
0221 
0222 void U4Process::ClearNumberOfInteractionLengthLeft(const G4Track& aTrack, const G4Step& aStep)
0223 {
0224     G4ProcessManager* mgr = GetManager(); 
0225     G4ProcessVector* procv = mgr->GetProcessList() ;
0226     for(int i=0 ; i < int(procv->entries()) ; i++)
0227     {
0228         G4VProcess* proc = (*procv)[i] ;
0229         unsigned type = ProcessType(Name(proc)) ;   
0230         if(IsNormalProcess(type))
0231         {
0232             G4VDiscreteProcess* dproc = dynamic_cast<G4VDiscreteProcess*>(proc) ;
0233             assert(dproc);
0234             dproc->G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
0235         }
0236     }
0237 }
0238 
0239 
0240