#include <RadLengthAction.h>
Definition at line 25 of file RadLengthAction.h.
 
◆ RadLengthAction()
      
        
          | G4UA::RadLengthAction::RadLengthAction | ( | const Config & | config | ) |  | 
      
 
 
◆ BeginOfEventAction()
  
  | 
        
          | void G4UA::RadLengthAction::BeginOfEventAction | ( | const G4Event * | anEvent | ) |  |  | overridevirtual | 
 
Definition at line 169 of file RadLengthAction.cxx.
  173     G4PrimaryVertex *vert=anEvent->GetPrimaryVertex(0);
 
  175     G4PrimaryParticle *
part=vert->GetPrimary();
 
  176     G4ThreeVector 
mom=
part->GetMomentum();
 
  180     G4cout<<
"Begin Of Event"<<
" "<<(
double)vert->GetX0()<<
" "<<(
double)vert->GetY0()<<
" "<<(
double)vert->GetZ0()<< G4endl;
 
 
 
 
◆ BeginOfRunAction()
  
  | 
        
          | void G4UA::RadLengthAction::BeginOfRunAction | ( | const G4Run * |  | ) |  |  | overridevirtual | 
 
Definition at line 50 of file RadLengthAction.cxx.
   57     G4SDManager *SDM=G4SDManager::GetSDMpointer();
 
   58     m_SDMDT = (G4VSensitiveDetector*) SDM->FindSensitiveDetector(
"AMS");
 
   59     m_SDTGC = (G4VSensitiveDetector*) SDM->FindSensitiveDetector(
"TMS");
 
   60     m_SDCSC = (G4VSensitiveDetector*) SDM->FindSensitiveDetector(
"CMS");
 
   61     m_SDRPC = (G4VSensitiveDetector*) SDM->FindSensitiveDetector(
"RMS");
 
   70     G4LogicalVolumeStore *
store = G4LogicalVolumeStore::GetInstance();
 
   72     G4LogicalVolume* atlas = 
store->GetVolume(
"Atlas::Atlas");
 
   76     std::vector<G4LogicalVolume*> logvec;
 
   77     logvec.push_back(atlas);
 
   80     std::vector<G4VPhysicalVolume*> physvec;
 
   91       std::vector<G4LogicalVolume*> tmplogvec;
 
   94       for(G4LogicalVolume* logvol : logvec) {
 
   95         for(
unsigned int k=0; 
k<logvol->GetNoDaughters();
k++){
 
   97           tmplogvec.push_back(logvol->GetDaughter(
k)->GetLogicalVolume());
 
  100             physvec.push_back(logvol->GetDaughter(
k));
 
  113     for (G4VPhysicalVolume* physvol : physvec) {
 
  114       std::string fulldaughtername = physvol->GetName();
 
  115       std::string daughtername;
 
  116       std::string::size_type npos;
 
  117       npos=fulldaughtername.find(
"::");
 
  120       else daughtername = fulldaughtername;
 
  132       treeMap[
p.first]=
new TTree(TString(
p.first), TString(
p.first));
 
  137     if(
m_hSvc.retrieve().isFailure()) 
return;
 
  143       std::string 
filename= 
"/RadLengthAction/";
 
 
 
 
◆ EndOfEventAction()
  
  | 
        
          | void G4UA::RadLengthAction::EndOfEventAction | ( | const G4Event * |  | ) |  |  | overridevirtual | 
 
 
◆ EndOfRunAction()
  
  | 
        
          | void G4UA::RadLengthAction::EndOfRunAction | ( | const G4Run * |  | ) |  |  | overridevirtual | 
 
 
◆ fillVariables()
  
  | 
        
          | void G4UA::RadLengthAction::fillVariables | ( | const std::vector< double > & | varvec, |  
          |  |  | const std::string & | name |  
          |  | ) |  |  |  | private | 
 
 
◆ UserSteppingAction()
  
  | 
        
          | void G4UA::RadLengthAction::UserSteppingAction | ( | const G4Step * | aStep | ) |  |  | overridevirtual | 
 
Definition at line 218 of file RadLengthAction.cxx.
  228       const G4TouchableHistory* touchHist =
 
  229         static_cast<const G4TouchableHistory*
>(aStep->GetPreStepPoint()->GetTouchable());
 
  232       G4ThreeVector 
xyz = (G4ThreeVector) aStep->GetPreStepPoint()->GetPosition();
 
  237       double radl=aStep->GetPreStepPoint()->GetMaterial()->GetRadlen();
 
  238       double intl=aStep->GetPreStepPoint()->GetMaterial()->GetNuclearInterLength();
 
  241       double stepl=aStep->GetStepLength();
 
  242       G4cout<<aStep->GetPreStepPoint()->GetMaterial()->GetName()<<
" "<<radl<<
" "<<stepl<<
" "<<stepl/radl<<
" "<<(
double) 
xyz[0]<<
" "<<(
double) 
xyz[1]<<
" "<< (
double) 
xyz[2]<<G4endl;
 
  244       const unsigned int nVariablesToSave(12);
 
  245       std::vector<double> varvec(nVariablesToSave);
 
  246       varvec.at(0)  = (
double) aStep->GetDeltaEnergy();
 
  247       varvec.at(1)  = (
double) stepl/radl;
 
  248       varvec.at(2)  = (
double) stepl/intl;
 
  251       varvec.at(5)  = (
double) 
xyz.pseudoRapidity();
 
  268       G4VSensitiveDetector* SD = (G4VSensitiveDetector*) touchHist->GetVolume()->GetLogicalVolume()->GetSensitiveDetector();
 
 
 
 
◆ chargePrimary
  
  | 
        
          | double G4UA::RadLengthAction::chargePrimary = 0.0 |  | private | 
 
 
◆ etaPrimary
  
  | 
        
          | double G4UA::RadLengthAction::etaPrimary = 0.0 |  | private | 
 
 
◆ m_config
  
  | 
        
          | Config G4UA::RadLengthAction::m_config |  | private | 
 
 
◆ m_hSvc
◆ m_SDCSC
  
  | 
        
          | G4VSensitiveDetector* G4UA::RadLengthAction::m_SDCSC |  | private | 
 
 
◆ m_SDMDT
  
  | 
        
          | G4VSensitiveDetector* G4UA::RadLengthAction::m_SDMDT |  | private | 
 
 
◆ m_SDRPC
  
  | 
        
          | G4VSensitiveDetector* G4UA::RadLengthAction::m_SDRPC |  | private | 
 
 
◆ m_SDTGC
  
  | 
        
          | G4VSensitiveDetector* G4UA::RadLengthAction::m_SDTGC |  | private | 
 
 
◆ MuChamberPassed
  
  | 
        
          | bool G4UA::RadLengthAction::MuChamberPassed = false |  | private | 
 
 
◆ MuTriggerPassed
  
  | 
        
          | bool G4UA::RadLengthAction::MuTriggerPassed = false |  | private | 
 
 
◆ phiPrimary
  
  | 
        
          | double G4UA::RadLengthAction::phiPrimary = 0.0 |  | private | 
 
 
◆ topvolmap
  
  | 
        
          | std::map<std::string,G4VPhysicalVolume*> G4UA::RadLengthAction::topvolmap |  | private | 
 
 
◆ treeMap
  
  | 
        
          | std::map<std::string,TTree*> G4UA::RadLengthAction::treeMap |  | private | 
 
 
◆ variables
  
  | 
        
          | std::map<std::string,std::vector<double> > G4UA::RadLengthAction::variables |  | private | 
 
 
The documentation for this class was generated from the following files: