Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-05-02 07:43:10

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file LXeSteppingAction.cc
0027 /// \brief Implementation of the LXeSteppingAction class
0028 
0029 #include "LXeSteppingAction.hh"
0030 
0031 #include "LXeEventAction.hh"
0032 #include "LXePMTSD.hh"
0033 #include "LXeSteppingMessenger.hh"
0034 #include "LXeTrajectory.hh"
0035 #include "LXeUserTrackInformation.hh"
0036 
0037 #include "G4OpticalPhoton.hh"
0038 #include "G4ProcessManager.hh"
0039 #include "G4SDManager.hh"
0040 #include "G4Step.hh"
0041 #include "G4StepPoint.hh"
0042 #include "G4SteppingManager.hh"
0043 #include "G4Track.hh"
0044 #include "G4TrackStatus.hh"
0045 #include "G4VPhysicalVolume.hh"
0046 
0047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0048 
0049 LXeSteppingAction::LXeSteppingAction(LXeEventAction* ea) : fEventAction(ea)
0050 {
0051   fSteppingMessenger = new LXeSteppingMessenger(this);
0052 }
0053 
0054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0055 
0056 LXeSteppingAction::~LXeSteppingAction()
0057 {
0058   delete fSteppingMessenger;
0059 }
0060 
0061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0062 
0063 void LXeSteppingAction::UserSteppingAction(const G4Step* theStep)
0064 {
0065   G4Track* theTrack = theStep->GetTrack();
0066   const G4ParticleDefinition* part = theTrack->GetDefinition();
0067   G4int pdg = part->GetPDGEncoding();
0068 
0069   if (theTrack->GetCurrentStepNumber() == 1) fExpectedNextStatus = Undefined;
0070 
0071   auto trackInformation = static_cast<LXeUserTrackInformation*>(theTrack->GetUserInformation());
0072 
0073   G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
0074   G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
0075 
0076   G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
0077   G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
0078 
0079   G4OpBoundaryProcessStatus boundaryStatus = Undefined;
0080 
0081   // find the boundary process only once
0082   if (nullptr == fBoundary && pdg == -22) {
0083     G4ProcessManager* pm = part->GetProcessManager();
0084     G4int nprocesses = pm->GetProcessListLength();
0085     G4ProcessVector* pv = pm->GetProcessList();
0086     for (G4int i = 0; i < nprocesses; ++i) {
0087       if (nullptr != (*pv)[i] && (*pv)[i]->GetProcessName() == "OpBoundary") {
0088         fBoundary = dynamic_cast<G4OpBoundaryProcess*>((*pv)[i]);
0089         break;
0090       }
0091     }
0092   }
0093 
0094   if (theTrack->GetParentID() == 0) {
0095     // This is a primary track
0096     auto secondaries = theStep->GetSecondaryInCurrentStep();
0097     // If we haven't already found the conversion position and there were
0098     // secondaries generated, then search for it
0099     // since this is happening before the secondary is being tracked,
0100     // the vertex position has not been set yet (set in initial step)
0101     if (nullptr != secondaries && !fEventAction->IsConvPosSet()) {
0102       if (!secondaries->empty()) {
0103         for (auto& tr : *secondaries) {
0104           const G4VProcess* creator = tr->GetCreatorProcess();
0105           if (nullptr != creator) {
0106             G4int type = creator->GetProcessSubType();
0107             // 12 - photoeffect
0108             // 13 - Compton scattering
0109             // 14 - gamma conversion
0110             if (type >= 12 && type <= 14) {
0111               fEventAction->SetConvPos(tr->GetPosition());
0112             }
0113           }
0114         }
0115       }
0116     }
0117     if (fOneStepPrimaries && thePrePV->GetName() == "scintillator")
0118       theTrack->SetTrackStatus(fStopAndKill);
0119   }
0120 
0121   if (nullptr == thePostPV) {  // out of world
0122     fExpectedNextStatus = Undefined;
0123     return;
0124   }
0125 
0126   // Optical photon only
0127   if (pdg == -22) {
0128     if (thePrePV->GetName() == "Slab")
0129       // force drawing of photons in WLS slab
0130       trackInformation->SetForceDrawTrajectory(true);
0131     else if (thePostPV->GetName() == "expHall")
0132       // Kill photons entering expHall from something other than Slab
0133       theTrack->SetTrackStatus(fStopAndKill);
0134 
0135     // Was the photon absorbed by the absorption process
0136     auto proc = thePostPoint->GetProcessDefinedStep();
0137     if (nullptr != proc && proc->GetProcessName() == "OpAbsorption") {
0138       fEventAction->IncAbsorption();
0139       trackInformation->AddTrackStatusFlag(absorbed);
0140     }
0141     if (nullptr != fBoundary) boundaryStatus = fBoundary->GetStatus();
0142 
0143     if (thePostPoint->GetStepStatus() == fGeomBoundary) {
0144       // Check to see if the particle was actually at a boundary
0145       // Otherwise the boundary status may not be valid
0146       if (fExpectedNextStatus == StepTooSmall) {
0147         if (boundaryStatus != StepTooSmall) {
0148           G4cout << "LXeSteppingAction::UserSteppingAction(): "
0149                  << "trackID=" << theTrack->GetTrackID() << " parentID=" << theTrack->GetParentID()
0150                  << " " << part->GetParticleName() << " E(MeV)=" << theTrack->GetKineticEnergy()
0151                  << "n/ at " << theTrack->GetPosition() << " prePV: " << thePrePV->GetName()
0152                  << " postPV: " << thePostPV->GetName() << G4endl;
0153           G4ExceptionDescription ed;
0154           ed << "LXeSteppingAction: "
0155              << "No reallocation step after reflection!"
0156              << "Something is wrong with the surface normal or geometry";
0157           G4Exception("LXeSteppingAction:", "LXeExpl01", JustWarning, ed, "");
0158           return;
0159         }
0160       }
0161       fExpectedNextStatus = Undefined;
0162       switch (boundaryStatus) {
0163         case Absorption:
0164           trackInformation->AddTrackStatusFlag(boundaryAbsorbed);
0165           fEventAction->IncBoundaryAbsorption();
0166           break;
0167         case Detection:  // Note, this assumes that the volume causing detection
0168                          // is the photocathode because it is the only one with
0169                          // non-zero efficiency
0170         {
0171           // Trigger sensitive detector manually since photon is
0172           // absorbed but status was Detection
0173           G4SDManager* SDman = G4SDManager::GetSDMpointer();
0174           G4String sdName = "/LXeDet/pmtSD";
0175           LXePMTSD* pmtSD = (LXePMTSD*)SDman->FindSensitiveDetector(sdName);
0176           if (pmtSD) pmtSD->ProcessHits_boundary(theStep, nullptr);
0177           trackInformation->AddTrackStatusFlag(hitPMT);
0178           break;
0179         }
0180         case FresnelReflection:
0181         case TotalInternalReflection:
0182         case LambertianReflection:
0183         case LobeReflection:
0184         case SpikeReflection:
0185         case BackScattering:
0186           trackInformation->IncReflections();
0187           fExpectedNextStatus = StepTooSmall;
0188           break;
0189         default:
0190           break;
0191       }
0192       if (thePostPV->GetName() == "sphere") trackInformation->AddTrackStatusFlag(hitSphere);
0193     }
0194   }
0195 }