Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:22:22

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 //
0027 /// \file optical/LXe/src/LXeSteppingAction.cc
0028 /// \brief Implementation of the LXeSteppingAction class
0029 //
0030 //
0031 #include "LXeSteppingAction.hh"
0032 
0033 #include "LXeEventAction.hh"
0034 #include "LXePMTSD.hh"
0035 #include "LXeSteppingMessenger.hh"
0036 #include "LXeTrajectory.hh"
0037 #include "LXeUserTrackInformation.hh"
0038 
0039 #include "G4OpticalPhoton.hh"
0040 #include "G4ProcessManager.hh"
0041 #include "G4SDManager.hh"
0042 #include "G4Step.hh"
0043 #include "G4StepPoint.hh"
0044 #include "G4SteppingManager.hh"
0045 #include "G4Track.hh"
0046 #include "G4TrackStatus.hh"
0047 #include "G4VPhysicalVolume.hh"
0048 
0049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0050 
0051 LXeSteppingAction::LXeSteppingAction(LXeEventAction* ea) : fEventAction(ea)
0052 {
0053   fSteppingMessenger = new LXeSteppingMessenger(this);
0054 }
0055 
0056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0057 
0058 LXeSteppingAction::~LXeSteppingAction()
0059 {
0060   delete fSteppingMessenger;
0061 }
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 void LXeSteppingAction::UserSteppingAction(const G4Step* theStep)
0066 {
0067   G4Track* theTrack = theStep->GetTrack();
0068   const G4ParticleDefinition* part = theTrack->GetDefinition();
0069   G4int pdg = part->GetPDGEncoding();
0070 
0071   if (theTrack->GetCurrentStepNumber() == 1) fExpectedNextStatus = Undefined;
0072 
0073   auto trackInformation = static_cast<LXeUserTrackInformation*>(theTrack->GetUserInformation());
0074 
0075   G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
0076   G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
0077 
0078   G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
0079   G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
0080 
0081   G4OpBoundaryProcessStatus boundaryStatus = Undefined;
0082 
0083   // find the boundary process only once
0084   if (nullptr == fBoundary && pdg == -22) {
0085     G4ProcessManager* pm = part->GetProcessManager();
0086     G4int nprocesses = pm->GetProcessListLength();
0087     G4ProcessVector* pv = pm->GetProcessList();
0088     for (G4int i = 0; i < nprocesses; ++i) {
0089       if (nullptr != (*pv)[i] && (*pv)[i]->GetProcessName() == "OpBoundary") {
0090         fBoundary = dynamic_cast<G4OpBoundaryProcess*>((*pv)[i]);
0091         break;
0092       }
0093     }
0094   }
0095 
0096   if (theTrack->GetParentID() == 0) {
0097     // This is a primary track
0098     auto secondaries = theStep->GetSecondaryInCurrentStep();
0099     // If we haven't already found the conversion position and there were
0100     // secondaries generated, then search for it
0101     // since this is happening before the secondary is being tracked,
0102     // the vertex position has not been set yet (set in initial step)
0103     if (nullptr != secondaries && !fEventAction->IsConvPosSet()) {
0104       if (!secondaries->empty()) {
0105         for (auto& tr : *secondaries) {
0106           const G4VProcess* creator = tr->GetCreatorProcess();
0107           if (nullptr != creator) {
0108             G4int type = creator->GetProcessSubType();
0109             // 12 - photoeffect
0110             // 13 - Compton scattering
0111             // 14 - gamma conversion
0112             if (type >= 12 && type <= 14) {
0113               fEventAction->SetConvPos(tr->GetPosition());
0114             }
0115           }
0116         }
0117       }
0118     }
0119     if (fOneStepPrimaries && thePrePV->GetName() == "scintillator")
0120       theTrack->SetTrackStatus(fStopAndKill);
0121   }
0122 
0123   if (nullptr == thePostPV) {  // out of world
0124     fExpectedNextStatus = Undefined;
0125     return;
0126   }
0127 
0128   // Optical photon only
0129   if (pdg == -22) {
0130     if (thePrePV->GetName() == "Slab")
0131       // force drawing of photons in WLS slab
0132       trackInformation->SetForceDrawTrajectory(true);
0133     else if (thePostPV->GetName() == "expHall")
0134       // Kill photons entering expHall from something other than Slab
0135       theTrack->SetTrackStatus(fStopAndKill);
0136 
0137     // Was the photon absorbed by the absorption process
0138     auto proc = thePostPoint->GetProcessDefinedStep();
0139     if (nullptr != proc && proc->GetProcessName() == "OpAbsorption") {
0140       fEventAction->IncAbsorption();
0141       trackInformation->AddTrackStatusFlag(absorbed);
0142     }
0143     if (nullptr != fBoundary) boundaryStatus = fBoundary->GetStatus();
0144 
0145     if (thePostPoint->GetStepStatus() == fGeomBoundary) {
0146       // Check to see if the particle was actually at a boundary
0147       // Otherwise the boundary status may not be valid
0148       if (fExpectedNextStatus == StepTooSmall) {
0149         if (boundaryStatus != StepTooSmall) {
0150           G4cout << "LXeSteppingAction::UserSteppingAction(): "
0151                  << "trackID=" << theTrack->GetTrackID() << " parentID=" << theTrack->GetParentID()
0152                  << " " << part->GetParticleName() << " E(MeV)=" << theTrack->GetKineticEnergy()
0153                  << "n/ at " << theTrack->GetPosition() << " prePV: " << thePrePV->GetName()
0154                  << " postPV: " << thePostPV->GetName() << G4endl;
0155           G4ExceptionDescription ed;
0156           ed << "LXeSteppingAction: "
0157              << "No reallocation step after reflection!"
0158              << "Something is wrong with the surface normal or geometry";
0159           G4Exception("LXeSteppingAction:", "LXeExpl01", JustWarning, ed, "");
0160           return;
0161         }
0162       }
0163       fExpectedNextStatus = Undefined;
0164       switch (boundaryStatus) {
0165         case Absorption:
0166           trackInformation->AddTrackStatusFlag(boundaryAbsorbed);
0167           fEventAction->IncBoundaryAbsorption();
0168           break;
0169         case Detection:  // Note, this assumes that the volume causing detection
0170                          // is the photocathode because it is the only one with
0171                          // non-zero efficiency
0172         {
0173           // Trigger sensitive detector manually since photon is
0174           // absorbed but status was Detection
0175           G4SDManager* SDman = G4SDManager::GetSDMpointer();
0176           G4String sdName = "/LXeDet/pmtSD";
0177           LXePMTSD* pmtSD = (LXePMTSD*)SDman->FindSensitiveDetector(sdName);
0178           if (pmtSD) pmtSD->ProcessHits_boundary(theStep, nullptr);
0179           trackInformation->AddTrackStatusFlag(hitPMT);
0180           break;
0181         }
0182         case FresnelReflection:
0183         case TotalInternalReflection:
0184         case LambertianReflection:
0185         case LobeReflection:
0186         case SpikeReflection:
0187         case BackScattering:
0188           trackInformation->IncReflections();
0189           fExpectedNextStatus = StepTooSmall;
0190           break;
0191         default:
0192           break;
0193       }
0194       if (thePostPV->GetName() == "sphere") trackInformation->AddTrackStatusFlag(hitSphere);
0195     }
0196   }
0197 }