aboutsummaryrefslogtreecommitdiffstats
path: root/src/G4BeamTestSC4SD.cxx
blob: f9ef766702244ebee59e1465a6b59bf59aa2c4f0 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90

#include <sstream>

#include "G4UnitsTable.hh"

#include "G4BeamTestSC4SD.h"
#include "G4HCofThisEvent.hh"
#include "G4Step.hh"
#include "G4ThreeVector.hh"
#include "G4SDManager.hh"
#include "G4ios.hh"


//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

G4BeamTestSC4SD::G4BeamTestSC4SD(const G4String& name, const G4String& hitsCollectionName)
  : G4VSensitiveDetector(name)
      , fHitsCollection(NULL)
{
   collectionName.insert(hitsCollectionName);
}

//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

G4BeamTestSC4SD::~G4BeamTestSC4SD()
{
}

//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

void G4BeamTestSC4SD::Initialize(G4HCofThisEvent* hce)
{
  // Create hits collection

  fHitsCollection
    = new G4BeamTestSiHitsCollection(SensitiveDetectorName, collectionName[0]);

  // Add this collection in hce

  G4int hcID
    = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
  hce->AddHitsCollection( hcID, fHitsCollection );
}

//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

G4bool G4BeamTestSC4SD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
{
  //  G4double stepl = 10;
  // aStep->SetStepLength(stepl);

  G4String name=aStep->GetTrack()->GetDefinition()->GetParticleName();


  if (name == "pi-" || name == "e-") {
    G4cout << " Particle_name hit SC4 = " <<  name  << G4endl;
    // total energy
    G4double etot = aStep->GetTrack()->GetTotalEnergy();
    // energy deposit
    G4double edep = aStep->GetTotalEnergyDeposit();

    G4BeamTestSiHit* newHit = new G4BeamTestSiHit();
    // G4StepPoint* preStepPoint = aStep->GetPreStepPoint();

    newHit->SetTrackID  (aStep->GetTrack()->GetTrackID());
    newHit->SetEdep(edep);
    newHit->SetPos (aStep->GetPostStepPoint()->GetPosition());
    newHit-> SetTime (aStep->GetPreStepPoint()->GetProperTime());
    fHitsCollection->insert( newHit );

    return true;
  }
}

//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

void G4BeamTestSC4SD::EndOfEvent(G4HCofThisEvent*)
{
  if ( verboseLevel>1 ) {
    G4int nofHits = fHitsCollection->entries();
    G4cout << G4endl
	   << "-------->Hits Collection: in this event they are " << nofHits
	   << " hits in the tracker chambers: " << G4endl;
      for ( G4int i=0; i<nofHits; i++ ) (*fHitsCollection)[i]->Print();

  }
}

//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......