Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DMXMaxTimeCuts.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 // --------------------------------------------------------------
28 // GEANT 4 - Underground Dark Matter Detector Advanced Example
29 //
30 // For information related to this code contact: Alex Howard
31 // e-mail: alexander.howard@cern.ch
32 // --------------------------------------------------------------
33 // Comments
34 //
35 // Underground Advanced
36 // by A. Howard and H. Araujo
37 // (27th November 2001)
38 //
39 // MaxTimeCuts program
40 // --------------------------------------------------------------
41 
42 #include "DMXMaxTimeCuts.hh"
43 
44 #include "G4PhysicalConstants.hh"
45 #include "G4Step.hh"
46 #include "G4UserLimits.hh"
47 #include "G4VParticleChange.hh"
48 
50  : DMXSpecialCuts(aName)
51 {
52  if (verboseLevel>1) {
53  G4cout << GetProcessName() << " is created "<< G4endl;
54  }
56 }
57 
59 {}
60 
62  : DMXSpecialCuts()
63 {}
64 
65 
67  const G4Track& aTrack,
68  G4double ,
70  )
71 {
72  // condition is set to "Not Forced"
73  *condition = NotForced;
74 
75  G4double proposedStep = DBL_MAX;
76  // get the pointer to UserLimits
77  G4UserLimits* pUserLimits = aTrack.GetVolume()->GetLogicalVolume()->GetUserLimits();
78  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
79 
80  // can apply cuts for specific particles - use if(particleDef):
81  // G4ParticleDefinition* aParticleDef = aTrack.GetDefinition();
82 
83  // G4cout << " Time: " << pUserLimits->GetUserMaxTime(aTrack) << G4endl;
84 
85  if (pUserLimits) {
86  G4double temp = DBL_MAX;
87  //max time limit
88  G4double dTime= (pUserLimits->GetUserMaxTime(aTrack) - aTrack.GetGlobalTime());
89  if (dTime < 0. ) {
90  proposedStep = 0.;
91  } else {
92  G4double beta = (aParticle->GetTotalMomentum())/(aParticle->GetTotalEnergy());
93  temp = beta*c_light*dTime;
94  if (proposedStep > temp) proposedStep = temp;
95  }
96 
97  }
98  return proposedStep;
99 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetTotalMomentum() const
G4GLOB_DLL std::ostream G4cout
virtual G4double GetUserMaxTime(const G4Track &)
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4UserLimits * GetUserLimits() const
virtual ~DMXMaxTimeCuts()
G4LogicalVolume * GetLogicalVolume() const
DMXMaxTimeCuts(const G4String &processName="DMXMaxTimeCuts")
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetProcessType(G4ProcessType)
Definition: G4VProcess.hh:420
double G4double
Definition: G4Types.hh:76
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
float c_light
Definition: hepunit.py:257