Geant4_10
G4UserSpecialCuts.cc
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: G4UserSpecialCuts.cc 68048 2013-03-13 14:34:07Z gcosmo $
28 //
29 // --------------------------------------------------------------
30 // History
31 //
32 // 15-04-98 first implementation, mma
33 // 07-04-03 migrade to cut per region (V.Ivanchenko)
34 // 18-09-03 substitute manager for the loss tables (V.Ivanchenko)
35 // 23-01-04 add protection for charged geantino in range cut (H.Kurashige)
36 // 09-09-04 tracking cut applied only if Rmin or Emin > DBL_MIN
37 // 21-01-11 changed order of checks: 1st energy tracking cut, 2nd track
38 // length, 3d time cut, 4th range for charged particles
39 // with non-zero mass; removed string comparisons (V.Ivanchenko)
40 // --------------------------------------------------------------
41 
42 #include "G4UserSpecialCuts.hh"
44 
45 #include "G4PhysicalConstants.hh"
46 #include "G4Step.hh"
47 #include "G4UserLimits.hh"
48 #include "G4VParticleChange.hh"
49 #include "G4LossTableManager.hh"
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 
54  : G4VProcess(aName, fGeneral )
55 {
56  // set Process Sub Type
57  SetProcessSubType(static_cast<int>(USER_SPECIAL_CUTS));
58 
59  if (verboseLevel>0)
60  {
61  G4cout << GetProcessName() << " is created "<< G4endl;
62  }
63  theLossTableManager = G4LossTableManager::Instance();
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69 {
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75  : G4VProcess(right)
76 {
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
83  G4double, // previousStepSize
85 {
86  // condition is set to "Not Forced"
87  *condition = NotForced;
88 
89  G4double ProposedStep = DBL_MAX;
90  G4UserLimits* pUserLimits =
92  if (pUserLimits)
93  {
94  // check max kinetic energy first
95  //
96  G4double Ekine = aTrack.GetKineticEnergy();
97  if(Ekine <= pUserLimits->GetUserMinEkine(aTrack)) { return 0.; }
98 
99  // max track length
100  //
101  ProposedStep = (pUserLimits->GetUserMaxTrackLength(aTrack)
102  - aTrack.GetTrackLength());
103  if (ProposedStep < 0.) { return 0.; }
104 
105  // max time limit
106  //
107  G4double tlimit = pUserLimits->GetUserMaxTime(aTrack);
108  if(tlimit < DBL_MAX) {
109  G4double beta = (aTrack.GetDynamicParticle()->GetTotalMomentum())
110  /(aTrack.GetTotalEnergy());
111  G4double dTime = (tlimit - aTrack.GetGlobalTime());
112  G4double temp = beta*c_light*dTime;
113  if (temp < 0.) { return 0.; }
114  if (ProposedStep > temp) { ProposedStep = temp; }
115  }
116 
117  // min remaining range
118  // (only for charged particle except for chargedGeantino)
119  //
120  G4double Rmin = pUserLimits->GetUserMinRange(aTrack);
121  if (Rmin > DBL_MIN) {
122  G4ParticleDefinition* Particle = aTrack.GetDefinition();
123  if ( (Particle->GetPDGCharge() != 0.) && (Particle->GetPDGMass() > 0.0))
124  {
125  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
126  G4double RangeNow = theLossTableManager->GetRange(Particle,Ekine,couple);
127  G4double temp = RangeNow - Rmin;
128  if (temp < 0.) { return 0.; }
129  if (ProposedStep > temp) { ProposedStep = temp; }
130  }
131  }
132  }
133  return ProposedStep;
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
140  const G4Step& )
141 //
142 // Kill the current particle, if requested by G4UserLimits
143 //
144 {
145  aParticleChange.Initialize(aTrack);
149  return &aParticleChange;
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
154 
155 
156 
157 
158 
159 
160 
161 
162 
163 
164 
G4double condition(const G4ErrorSymMatrix &m)
G4ParticleDefinition * GetDefinition() const
G4UserSpecialCuts(const G4String &processName="UserSpecialCut")
virtual ~G4UserSpecialCuts()
static G4LossTableManager * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetUserMaxTrackLength(const G4Track &)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetTotalMomentum() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4double GetUserMaxTime(const G4Track &)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTrackLength() const
G4UserLimits * GetUserLimits() const
virtual void Initialize(const G4Track &)
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
void ProposeEnergy(G4double finalEnergy)
#define DBL_MIN
Definition: templates.hh:75
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
virtual G4double GetUserMinRange(const G4Track &)
G4VPhysicalVolume * GetVolume() const
G4double GetTotalEnergy() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
float c_light
Definition: hepunit.py:257