Geant4_10
G4ScoreSplittingProcess.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4ScoreSplittingProcess.hh 68733 2013-04-05 09:45:28Z gcosmo $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4ScoreSplittingProcess.hh
33 //
34 // Description:
35 // This process is used to split the length and energy
36 // of a step in a regular structure into sub-steps, and to
37 // call the scorers for each sub-volume.
38 // It invokes sensitive detectors assigned in the *mass*
39 // world.
40 //
41 // Design and first implementation: J. Apostolakis / M.Asai 2010
42 //---------------------------------------------------------------
43 
44 
45 #ifndef G4ScoreSplittingProcess_h
46 #define G4ScoreSplittingProcess_h 1
47 
48 #include "globals.hh"
49 class G4Step;
50 class G4Navigator;
52 class G4PathFinder;
53 class G4VTouchable;
54 class G4VPhysicalVolume;
55 class G4ParticleChange;
56 class G4EnergySplitter;
57 
58 #include "G4VProcess.hh"
59 #include "G4FieldTrack.hh"
60 #include "G4TouchableHandle.hh"
61 class G4TouchableHistory;
62 // #include "G4TouchableHistory.hh"
63 
64 //------------------------------------------
65 //
66 // G4ScoreSplittingProcess class
67 //
68 //------------------------------------------
69 
70 
71 // Class Description:
72 
74 {
75 public: // with description
76 
77  //------------------------
78  // Constructor/Destructor
79  //------------------------
80 
81  G4ScoreSplittingProcess(const G4String& processName = "ScoreSplittingProc",
83  virtual ~G4ScoreSplittingProcess();
84 
85  //--------------------------------------------------------------
86  // Process interface
87  //--------------------------------------------------------------
88 
89  void StartTracking(G4Track*);
90 
91  //------------------------------------------------------------------------
92  // GetPhysicalInteractionLength() and DoIt() methods for AtRest
93  //------------------------------------------------------------------------
94 
96  const G4Track& ,
98  );
99 
101  const G4Track& ,
102  const G4Step&
103  );
104 
105  //------------------------------------------------------------------------
106  // GetPhysicalInteractionLength() and DoIt() methods for AlongStep
107  //------------------------------------------------------------------------
108 
110  const G4Track&,
111  G4double ,
112  G4double ,
113  G4double&,
115  );
116 
118  const G4Track& ,
119  const G4Step&
120  );
121 
122  //-----------------------------------------------------------------------
123  // GetPhysicalInteractionLength() and DoIt() methods for PostStep
124  //-----------------------------------------------------------------------
125 
127  G4double previousStepSize,
129 
130  G4VParticleChange* PostStepDoIt(const G4Track& ,const G4Step& );
131 
132 private:
133  G4TouchableHistory* CreateTouchableForSubStep( G4int newVoxelNum, G4ThreeVector newPosition );
134 
135 private:
136  void CopyStepStart(const G4Step & step);
137 
138  G4Step * fSplitStep;
139  G4StepPoint *fSplitPreStepPoint;
140  G4StepPoint *fSplitPostStepPoint;
141 
142  G4VParticleChange dummyParticleChange;
143  G4ParticleChange xParticleChange;
144 
145  // G4TransportationManager* fTransportationManager;
146  // G4PathFinder* fPathFinder;
147 
148  // -------------------------------
149  // Touchables for the Split Step
150  // -------------------------------
151  G4TouchableHandle fOldTouchableH;
152  G4TouchableHandle fNewTouchableH;
153 
154  // Memory of Touchables of full step
155  G4TouchableHandle fInitialTouchableH;
156  G4TouchableHandle fFinalTouchableH;
157 
158  G4EnergySplitter *fpEnergySplitter;
159 
160  // ******************************************************
161  // For TESTS:
162  // ******************************************************
163 public:
164  void Verbose(const G4Step&) const;
165 };
166 
167 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4ScoreSplittingProcess(const G4String &processName="ScoreSplittingProc", G4ProcessType theType=fParameterisation)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
int G4int
Definition: G4Types.hh:78
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void Verbose(const G4Step &) const
Definition: G4Step.hh:76
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
double G4double
Definition: G4Types.hh:76
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4ForceCondition
G4GPILSelection
G4ProcessType