Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:17:43

0001 #include "PHG4CEmcTestBeamSteppingAction.h"
0002 
0003 #include "PHG4CEmcTestBeamDetector.h"
0004 
0005 #include <g4main/PHG4Hit.h>
0006 #include <g4main/PHG4HitContainer.h>
0007 #include <g4main/PHG4Hitv1.h>
0008 #include <g4main/PHG4Shower.h>
0009 #include <g4main/PHG4SteppingAction.h>  // for PHG4SteppingAction
0010 #include <g4main/PHG4TrackUserInfoV1.h>
0011 
0012 #include <phool/getClass.h>
0013 
0014 #include <Geant4/G4ParticleDefinition.hh>      // for G4ParticleDefinition
0015 #include <Geant4/G4ReferenceCountedHandle.hh>  // for G4ReferenceCountedHandle
0016 #include <Geant4/G4Step.hh>
0017 #include <Geant4/G4StepPoint.hh>   // for G4StepPoint
0018 #include <Geant4/G4StepStatus.hh>  // for fGeomBoundary, fUndefined
0019 #include <Geant4/G4String.hh>      // for G4String
0020 #include <Geant4/G4SystemOfUnits.hh>
0021 #include <Geant4/G4ThreeVector.hh>            // for G4ThreeVector
0022 #include <Geant4/G4TouchableHandle.hh>        // for G4TouchableHandle
0023 #include <Geant4/G4Track.hh>                  // for G4Track
0024 #include <Geant4/G4TrackStatus.hh>            // for fStopAndKill
0025 #include <Geant4/G4Types.hh>                  // for G4double
0026 #include <Geant4/G4VTouchable.hh>             // for G4VTouchable
0027 #include <Geant4/G4VUserTrackInformation.hh>  // for G4VUserTrackInformation
0028 
0029 #include <iostream>
0030 #include <string>  // for operator+, string, ope...
0031 
0032 class G4VPhysicalVolume;
0033 class PHCompositeNode;
0034 
0035 using namespace std;
0036 //____________________________________________________________________________..
0037 PHG4CEmcTestBeamSteppingAction::PHG4CEmcTestBeamSteppingAction(PHG4CEmcTestBeamDetector* detector)
0038   : PHG4SteppingAction(detector->GetName())
0039   , detector_(detector)
0040   , hits_(nullptr)
0041   , absorberhits_(nullptr)
0042   , hit(nullptr)
0043 {
0044 }
0045 
0046 //____________________________________________________________________________..
0047 bool PHG4CEmcTestBeamSteppingAction::UserSteppingAction(const G4Step* aStep, bool /*was_used*/)
0048 {
0049   G4TouchableHandle touch = aStep->GetPreStepPoint()->GetTouchableHandle();
0050   // get volume of the current step
0051   G4VPhysicalVolume* volume = touch->GetVolume();
0052 
0053   int whichactive = detector_->IsInCEmcTestBeam(volume);
0054 
0055   if (!whichactive)
0056   {
0057     return false;
0058   }
0059 
0060   // collect energy and track length step by step
0061   G4double edep = aStep->GetTotalEnergyDeposit() / GeV;
0062   G4double eion = (aStep->GetTotalEnergyDeposit() - aStep->GetNonIonizingEnergyDeposit()) / GeV;
0063 
0064   const G4Track* aTrack = aStep->GetTrack();
0065 
0066   // if this block stops everything, just put all kinetic energy into edep
0067   if (detector_->IsBlackHole())
0068   {
0069     edep = aTrack->GetKineticEnergy() / GeV;
0070     G4Track* killtrack = const_cast<G4Track*>(aTrack);
0071     killtrack->SetTrackStatus(fStopAndKill);
0072   }
0073 
0074   int tower_id = touch->GetCopyNumber(2);  // tower number
0075   // make sure we are in a volume
0076   if (detector_->IsActive())
0077   {
0078     bool geantino = false;
0079     // the check for the pdg code speeds things up, I do not want to make
0080     // an expensive string compare for every track when we know
0081     // geantino or chargedgeantino has pid=0
0082     if (aTrack->GetParticleDefinition()->GetPDGEncoding() == 0 &&
0083         aTrack->GetParticleDefinition()->GetParticleName().find("geantino") != string::npos)
0084     {
0085       geantino = true;
0086     }
0087     G4StepPoint* prePoint = aStep->GetPreStepPoint();
0088     G4StepPoint* postPoint = aStep->GetPostStepPoint();
0089     //       cout << "track id " << aTrack->GetTrackID() << endl;
0090     //       cout << "time prepoint: " << prePoint->GetGlobalTime() << endl;
0091     //       cout << "time postpoint: " << postPoint->GetGlobalTime() << endl;
0092     switch (prePoint->GetStepStatus())
0093     {
0094     case fGeomBoundary:
0095     case fUndefined:
0096       hit = new PHG4Hitv1();
0097       hit->set_layer((unsigned int) tower_id);
0098       hit->set_scint_id(touch->GetCopyNumber(1));  // the copy number of the sandwich
0099       // here we set the entrance values in cm
0100       hit->set_x(0, prePoint->GetPosition().x() / cm);
0101       hit->set_y(0, prePoint->GetPosition().y() / cm);
0102       hit->set_z(0, prePoint->GetPosition().z() / cm);
0103       // time in ns
0104       hit->set_t(0, prePoint->GetGlobalTime() / nanosecond);
0105       // set the track ID
0106       {
0107         hit->set_trkid(aTrack->GetTrackID());
0108         if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
0109         {
0110           if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
0111           {
0112             hit->set_trkid(pp->GetUserTrackId());
0113             hit->set_shower_id(pp->GetShower()->get_id());
0114           }
0115         }
0116       }
0117 
0118       // set the initial energy deposit
0119       hit->set_edep(0);
0120       hit->set_eion(0);  // only implemented for v5 otherwise empty
0121       PHG4HitContainer* hitcontainer;
0122       // here we do things which are different between scintillator and absorber hits
0123       if (whichactive > 0)  // return of isinCEmcTestDetector, > 0 hit in scintillator, < 0 hit in absorber
0124       {
0125         hitcontainer = hits_;
0126       }
0127       else
0128       {
0129         hitcontainer = absorberhits_;
0130       }
0131       // here we set what is common for scintillator and absorber hits
0132       hitcontainer->AddHit(tower_id, hit);
0133       if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
0134       {
0135         if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
0136         {
0137           pp->GetShower()->add_g4hit_id(hitcontainer->GetID(), hit->get_hit_id());
0138         }
0139       }
0140       break;
0141     default:
0142       break;
0143     }
0144     // here we just update the exit values, it will be overwritten
0145     // for every step until we leave the volume or the particle
0146     // ceases to exist
0147     hit->set_x(1, postPoint->GetPosition().x() / cm);
0148     hit->set_y(1, postPoint->GetPosition().y() / cm);
0149     hit->set_z(1, postPoint->GetPosition().z() / cm);
0150 
0151     hit->set_t(1, postPoint->GetGlobalTime() / nanosecond);
0152     // sum up the energy to get total deposited
0153     hit->set_edep(hit->get_edep() + edep);
0154     hit->set_eion(hit->get_eion() + eion);
0155     if (geantino)
0156     {
0157       hit->set_edep(-1);  // only energy=0 g4hits get dropped, this way geantinos survive the g4hit compression
0158       hit->set_eion(-1);
0159     }
0160     if (edep > 0)
0161     {
0162       if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
0163       {
0164         if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
0165         {
0166           pp->SetKeep(1);  // we want to keep the track
0167         }
0168       }
0169     }
0170 
0171     //       hit->identify();
0172     // return true to indicate the hit was used
0173     return true;
0174   }
0175   else
0176   {
0177     return false;
0178   }
0179 }
0180 
0181 //____________________________________________________________________________..
0182 void PHG4CEmcTestBeamSteppingAction::SetInterfacePointers(PHCompositeNode* topNode)
0183 {
0184   string hitnodename;
0185   string absorbernodename;
0186   if (detector_->SuperDetector() != "NONE")
0187   {
0188     hitnodename = "G4HIT_" + detector_->SuperDetector();
0189     absorbernodename = "G4HIT_ABSORBER_" + detector_->SuperDetector();
0190   }
0191   else
0192   {
0193     hitnodename = "G4HIT_" + detector_->GetName();
0194     absorbernodename = "G4HIT_ABSORBER_" + detector_->GetName();
0195   }
0196 
0197   // now look for the map and grab a pointer to it.
0198   hits_ = findNode::getClass<PHG4HitContainer>(topNode, hitnodename);
0199   absorberhits_ = findNode::getClass<PHG4HitContainer>(topNode, absorbernodename);
0200 
0201   // if we do not find the node it's messed up.
0202   if (!hits_)
0203   {
0204     std::cout << "PHG4CEmcTestBeamSteppingAction::SetTopNode - unable to find " << hitnodename << std::endl;
0205   }
0206   if (!absorberhits_)
0207   {
0208     if (Verbosity() > 0)
0209     {
0210       cout << "PHG4HcalSteppingAction::SetTopNode - unable to find " << absorbernodename << endl;
0211     }
0212   }
0213 }