Skip to content
Snippets Groups Projects
SteppingAction.cpp 2.69 KiB
Newer Older
Markus Frank's avatar
Markus Frank committed
#include "DDG4/Geant4StepHandler.h"
#include "DDG4/Geant4Converter.h"

#include "SteppingAction.h"
#include "EventAction.h"
#include "G4Step.hh"
Markus Frank's avatar
Markus Frank committed
#include "G4PVPlacement.hh"
#include "G4VSensitiveDetector.hh"


using namespace DD4hep;
using namespace DD4hep::Simulation;

SteppingAction::SteppingAction(EventAction* evt) : eventaction(evt)					 
Markus Frank's avatar
Markus Frank committed

void SteppingAction::UserSteppingAction(const G4Step* aStep)   {  
  Geant4StepHandler step(aStep);
  Geant4Mapping&    mapping = Geant4Mapping::instance();
Markus Frank's avatar
Markus Frank committed
  SiMaterial     = G4Material::GetMaterial("Silicon");
  TPCGasMaterial = G4Material::GetMaterial("Argon");  

  // get volume of the current step
  G4VPhysicalVolume* volume = step.preVolume();
  G4Material* material = volume->GetLogicalVolume()->GetMaterial();
  
  // collect energy and track length step by step
  G4double edep = aStep->GetTotalEnergyDeposit();
  
  G4double stepl = 0.;
Markus Frank's avatar
Markus Frank committed
  if (step.trackDef()->GetPDGCharge() != 0.) stepl = aStep->GetStepLength();

  if (material == SiMaterial || material == TPCGasMaterial) { 
Markus Frank's avatar
Markus Frank committed
    eventaction->SumSensitive(edep, stepl, 0.);
Markus Frank's avatar
Markus Frank committed
    eventaction->SumSupport(edep, stepl, 0.);
  }
    
  //example of saving random number seed of this event, under condition
  //// if (condition) G4RunManager::GetRunManager()->rndmSaveThisEvent(); 
Markus Frank's avatar
Markus Frank committed

  Position pos1 = step.prePos();
  Position pos2 = step.postPos();
Markus Frank's avatar
Markus Frank committed
  
  if ( step.sd(step.pre) ) {
Markus Frank's avatar
Markus Frank committed
    ::printf("  Track:%08ld pos:%.0f Len:%.1f  SD:%s [%s] Deposit:%.0f keV\n",
	     long(step.track),pos2.R()/cm,len/cm,step.sdName(step.pre,"----"), step.preStepStatus(),
Markus Frank's avatar
Markus Frank committed
	     edep/keV);
  }
Markus Frank's avatar
Markus Frank committed
  ::printf("  Track:%08ld Pos:(%8f %8f %8f) -> (%f %f %f)  Mom:%7.0f %7.0f %7.0f \n",
	   long(step.track), pos1.X(), pos1.Y(), pos1.Z(), pos2.X(), pos2.Y(), pos2.Z(), mom.X(), mom.Y(), mom.Z());
Markus Frank's avatar
Markus Frank committed
  ::printf("                pre-Vol: %s  Status:%s  SD:%s\n",
	   step.volName(step.pre,"----"), step.preStepStatus(), step.sdName(step.pre,"----"));
  ::printf("                post-Vol:%s  Status:%s  SD:%s\n",
	   step.volName(step.post,"----"), step.postStepStatus(), step.sdName(step.post,"----"));
  const G4VPhysicalVolume* pv  = step.volume(step.post);
  Geometry::PlacedVolume place = mapping.placement(pv);
  if ( place.isValid() )   {
    if ( place.volume().isSensitive() )  {
      // Example code to access the physical vlume and the cell id
      Geometry::VolumeManager vm = mapping.lcdd().volumeManager();
      Geometry::VolumeManager::VolumeID cell_id = vm.lookupID(place);
      //const TGeoNode* tpv = pv.ptr();
      printf("           Found Sensitive TGeoNode:%s CellID: %lld!\n",place.name(),cell_id);
    }