/[PAMELA software]/PamG4/src/SteppingAction.cc
ViewVC logotype

Contents of /PamG4/src/SteppingAction.cc

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (show annotations) (download) (vendor branch)
Tue Dec 13 16:20:39 2005 UTC (19 years ago) by cafagna
Branch: start, MAIN
CVS Tags: v1r0, bogo, HEAD
Changes since 1.1: +0 -0 lines
Geant4 code for ND and Calorimeter

1
2 #include "SteppingAction.hh"
3 #include "DetectorConstruction.hh"
4 #include "G4Polyline.hh"
5 #include "G4Colour.hh"
6 #include "G4SteppingManager.hh"
7 #include "G4VisAttributes.hh"
8 #include "G4VVisManager.hh"
9 #include "G4VTouchable.hh"
10 #include "G4TouchableHistory.hh"
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 #include "RunAction.hh"
14 #include "G4ios.hh"
15 #include "G4UnitsTable.hh"
16 #include "SteppingActionMessenger.hh"
17 #include "G4ParticleDefinition.hh"
18 #include "G4ParticleTypes.hh"
19 #include "G4UImanager.hh"
20
21 SteppingAction::SteppingAction(DetectorConstruction* myDC
22 , RunAction* myRU)
23 :myDetector(myDC) ,runAction(myRU), drawFlag(false)
24
25 {new SteppingActionMessenger(this);}
26
27 SteppingAction::~SteppingAction()
28 { }
29
30
31 void SteppingAction::UserSteppingAction(const G4Step* aStep)
32 { const G4SteppingManager* pSM = fpSteppingManager;
33 G4TouchableHistory* theTouchable
34 = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
35 G4VPhysicalVolume* currentV = theTouchable->GetVolume();
36 G4String name = currentV->GetName();
37 G4double time = aStep->GetTrack()->GetGlobalTime();
38
39 if(drawFlag)
40 {
41 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
42 if(pVVisManager) {
43
44 //----------Define a line segment
45 G4Polyline polyline;
46 G4double charge = pSM->GetTrack()->GetDefinition()->GetPDGCharge();
47 G4Colour colour;
48 if (charge < 0.) colour = G4Colour(1., 0., 0.);
49 else if (charge > 0.) colour = G4Colour(0., 0., 1.);
50 else colour = G4Colour(0., 1., 0.);
51 if ( pSM->GetTrack()->GetDefinition()==G4Neutron::NeutronDefinition() ) colour = G4Colour(1., 1., 0.);
52 G4VisAttributes attribs(colour);
53 polyline.SetVisAttributes(attribs );
54 polyline.push_back(pSM->GetStep()->GetPreStepPoint()->GetPosition());
55 polyline.push_back(pSM->GetStep()->GetPostStepPoint()->GetPosition());
56 pVVisManager->Draw(polyline);
57 }
58 };
59
60
61 if(name=="physHe3Ar"){
62 G4int copyNo = theTouchable ->GetReplicaNumber(1);
63 G4double EnergyD = aStep->GetTotalEnergyDeposit();
64
65 if(EnergyD>0.)
66 {
67 if(time>1.*ns && time<150000.*ns)
68 {
69 G4String pname = aStep->GetTrack()->GetDefinition()->GetParticleName();
70 G4cout<<copyNo<<pname<<" time = "<<time/ns<<" energy/keV "<<EnergyD/keV
71 <<G4endl;
72
73 if(runAction->tube[copyNo]==0)
74 {
75 runAction->Etube[copyNo] =+EnergyD;
76 runAction->TimeTube[copyNo]=time;
77 runAction->tube[copyNo]=1;
78 }
79
80 if(runAction->tube[copyNo]==1)
81 {
82 if(fabs(time - runAction->TimeTube[copyNo])>50000.*ns)
83 {
84 runAction->E1tube[copyNo] += EnergyD;
85 }
86 else
87 runAction->Etube[copyNo] += EnergyD;}
88 }
89 }
90 }
91 }
92
93 // 2005 by G.I.Vasilyev

  ViewVC Help
Powered by ViewVC 1.1.23