source: trunk/source/processes/scoring/include/G4ScoreSplittingProcess.hh@ 1353

Last change on this file since 1353 was 1350, checked in by garnier, 15 years ago

update to last version 4.9.4

File size: 5.6 KB
Line 
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,v 1.4 2010/11/22 14:30:49 japost Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
29//
30//
31//---------------------------------------------------------------
32//
33// G4ScoreSplittingProcess.hh
34//
35// Description:
36// This process is used to split the length and energy
37// of a step in a regular structure into sub-steps, and to
38// call the scorers for each sub-volume.
39// It invokes sensitive detectors assigned in the *mass*
40// world.
41//
42// Design and first implementation: J. Apostolakis / M.Asai 2010
43//---------------------------------------------------------------
44
45
46#ifndef G4ScoreSplittingProcess_h
47#define G4ScoreSplittingProcess_h 1
48
49#include "globals.hh"
50class G4Step;
51class G4Navigator;
52class G4TransportationManager;
53class G4PathFinder;
54class G4VTouchable;
55class G4VPhysicalVolume;
56class G4ParticleChange;
57class G4EnergySplitter;
58
59#include "G4VProcess.hh"
60#include "G4FieldTrack.hh"
61#include "G4TouchableHandle.hh"
62class G4TouchableHistory;
63// #include "G4TouchableHistory.hh"
64
65//------------------------------------------
66//
67// G4ScoreSplittingProcess class
68//
69//------------------------------------------
70
71
72// Class Description:
73
74class G4ScoreSplittingProcess : public G4VProcess
75{
76public: // with description
77
78 //------------------------
79 // Constructor/Destructor
80 //------------------------
81
82 G4ScoreSplittingProcess(const G4String& processName = "ScoreSplittingProc",
83 G4ProcessType theType = fParameterisation);
84 virtual ~G4ScoreSplittingProcess();
85
86 //--------------------------------------------------------------
87 // Process interface
88 //--------------------------------------------------------------
89
90 void StartTracking(G4Track*);
91
92 //------------------------------------------------------------------------
93 // GetPhysicalInteractionLength() and DoIt() methods for AtRest
94 //------------------------------------------------------------------------
95
96 G4double AtRestGetPhysicalInteractionLength(
97 const G4Track& ,
98 G4ForceCondition*
99 );
100
101 G4VParticleChange* AtRestDoIt(
102 const G4Track& ,
103 const G4Step&
104 );
105
106 //------------------------------------------------------------------------
107 // GetPhysicalInteractionLength() and DoIt() methods for AlongStep
108 //------------------------------------------------------------------------
109
110 G4double AlongStepGetPhysicalInteractionLength(
111 const G4Track&,
112 G4double ,
113 G4double ,
114 G4double&,
115 G4GPILSelection*
116 );
117
118 G4VParticleChange* AlongStepDoIt(
119 const G4Track& ,
120 const G4Step&
121 );
122
123 //-----------------------------------------------------------------------
124 // GetPhysicalInteractionLength() and DoIt() methods for PostStep
125 //-----------------------------------------------------------------------
126
127 G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
128 G4double previousStepSize,
129 G4ForceCondition* condition);
130
131 G4VParticleChange* PostStepDoIt(const G4Track& ,const G4Step& );
132
133private:
134 G4TouchableHistory* CreateTouchableForSubStep( G4int newVoxelNum, G4ThreeVector newPosition );
135
136private:
137 void CopyStepStart(const G4Step & step);
138
139 G4Step * fSplitStep;
140 G4StepPoint *fSplitPreStepPoint;
141 G4StepPoint *fSplitPostStepPoint;
142
143 G4VParticleChange dummyParticleChange;
144 G4ParticleChange xParticleChange;
145
146 // G4TransportationManager* fTransportationManager;
147 // G4PathFinder* fPathFinder;
148
149 // -------------------------------
150 // Touchables for the Split Step
151 // -------------------------------
152 G4TouchableHandle fOldTouchableH;
153 G4TouchableHandle fNewTouchableH;
154
155 // Memory of Touchables of full step
156 G4TouchableHandle fInitialTouchableH;
157 G4TouchableHandle fFinalTouchableH;
158
159 G4EnergySplitter *fpEnergySplitter;
160
161 // ******************************************************
162 // For TESTS:
163 // ******************************************************
164public:
165 void Verbose(const G4Step&) const;
166};
167
168#endif
Note: See TracBrowser for help on using the repository browser.