/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSSteppingAction.cc

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 25.12.2003
00004    Copyright (c) 2003 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 */
00006 
00007 //====================================================
00008 //  Class description here ...
00009 //
00010 //====================================================
00011 
00012 #include "BDSSteppingAction.hh"
00013 #include "BDSExecOptions.hh"
00014 #include "BDSGlobalConstants.hh"
00015 
00016 #include "G4AffineTransform.hh"
00017 #include "G4NavigationHistory.hh"
00018 #include "G4Track.hh"
00019 #include "G4VProcess.hh"
00020 
00021 //====================================================
00022 
00023 BDSSteppingAction::BDSSteppingAction():_step(NULL)
00024 { 
00025 }
00026 
00027 //====================================================
00028 
00029 BDSSteppingAction::~BDSSteppingAction()
00030 {}
00031 
00032 //====================================================
00033 
00034 
00035 void BDSSteppingAction::UserSteppingAction(const G4Step* ThisStep){
00036   _step = ThisStep;
00037   if(BDSExecOptions::Instance()->GetVerboseStep()) {
00038     VerboseSteppingAction();
00039   }
00040   if (BDSGlobalConstants::Instance()->GetThresholdCutPhotons() > 0 || BDSGlobalConstants::Instance()->GetThresholdCutCharged() > 0) {
00041     ThresholdCutSteppingAction();
00042   }
00043 }
00044 
00045 void BDSSteppingAction::ThresholdCutSteppingAction(){
00046   // -------------  kill tracks according to cuts -------------------
00047   G4int pdgNr = _step->GetTrack()->GetParticleDefinition()->GetPDGEncoding();
00048   // this cuts apply to default region
00049   if (pdgNr == 22) {
00050     if(_step->GetTrack()->GetKineticEnergy()<BDSGlobalConstants::Instance()->GetThresholdCutPhotons())
00051       {
00052         _step->GetTrack()->SetTrackStatus(fStopAndKill);
00053       }
00054   } else if (abs(pdgNr) == 11) {
00055     //note this is 'thresholdcutcarged' but only works on electrons...
00056     if(_step->GetTrack()->GetKineticEnergy()<BDSGlobalConstants::Instance()->GetThresholdCutCharged())
00057       {
00058         _step->GetTrack()->SetTrackStatus(fStopAndKill);
00059       }
00060   }
00061 }
00062   
00063 void BDSSteppingAction::VerboseSteppingAction()
00064 { 
00065   // ------------  output in case of verbose step ---------------------
00066 
00067   int ID=_step->GetTrack()->GetTrackID();
00068   int G4precision = G4cout.precision();
00069   G4cout.precision(10);
00070   G4cout<<"This volume="<< _step->GetTrack()->GetVolume()->GetName()<<G4endl;
00071         
00072   G4LogicalVolume* LogVol=_step->GetTrack()->GetVolume()->GetLogicalVolume();
00073   G4cout<<"This log volume="<<LogVol->GetName() <<G4endl;
00074         
00075   G4cout<<"ID="<<ID<<" part="<<
00076     _step->GetTrack()->GetDefinition()->GetParticleName()<<
00077     "Energy="<<_step->GetTrack()->GetTotalEnergy()/CLHEP::GeV<<
00078     " mom Px="
00079         <<_step->GetTrack()->GetMomentum()[0]/CLHEP::GeV<<
00080     " Py="<<_step->GetTrack()->GetMomentum()[1]/CLHEP::GeV<<
00081     " Pz="<<_step->GetTrack()->GetMomentum()[2]/CLHEP::GeV<<" vol="<<
00082     _step->GetTrack()->GetVolume()->GetName()<<G4endl;
00083         
00084   G4cout<<" Global Position="<<_step->GetTrack()->GetPosition()<<G4endl;
00085   G4AffineTransform tf(_step->GetPreStepPoint()->GetTouchableHandle()->GetHistory()->GetTopTransform());
00086   G4cout<<" Local Position="<< tf.TransformPoint(_step->GetTrack()->GetPosition()) <<G4endl;
00087 
00088   if(_step->GetTrack()->GetMaterial()->GetName() !="LCVacuum")
00089     G4cout<<"material="<<_step->GetTrack()->GetMaterial()->GetName()<<G4endl;
00090 
00091   G4VProcess* proc=(G4VProcess*)(_step->GetPostStepPoint()->
00092                                  GetProcessDefinedStep() );
00093 
00094   if(proc)G4cout<<" post-step process="<<proc->GetProcessName()<<G4endl<<G4endl;
00095 
00096 
00097   proc=(G4VProcess*)(_step->GetPreStepPoint()->
00098                      GetProcessDefinedStep() );
00099 
00100   if(proc)G4cout<<" pre-step process="<<proc->GetProcessName()<<G4endl<<G4endl;
00101 
00102   // set precision back
00103   G4cout.precision(G4precision);
00104 }
00105   
00106 //====================================================

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7