Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PrimaryKiller Class Reference

#include <PrimaryKiller.hh>

Inheritance diagram for PrimaryKiller:
Collaboration diagram for PrimaryKiller:

Public Member Functions

 PrimaryKiller (G4String name, G4int depth=0)
 
virtual ~PrimaryKiller ()
 
void SetEnergyThreshold (double energy)
 
void SetMinLossEnergyLimit (double energy)
 
void SetMaxLossEnergyLimit (double energy)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()
 
G4int GetCollectionID (G4int)
 
void SetUnit (const G4String &unit)
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
void SetMultiFunctionalDetector (G4MultiFunctionalDetector *d)
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
void SetFilter (G4VSDFilter *f)
 
G4VSDFilterGetFilter () const
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel () const
 
void SetNijk (G4int i, G4int j, G4int k)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Protected Member Functions

virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 
- Protected Member Functions inherited from G4VPrimitiveScorer
virtual G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int verboseLevel
 
G4int indexDepth
 
G4String unitName
 
G4double unitValue
 
G4int fNi
 
G4int fNj
 
G4int fNk
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir
 
G4String baseDirName
 
G4bool commandsShouldBeInMaster
 

Detailed Description

Kill the primary particle:

  • either after a given energy loss
  • or after the primary particle has reached a given energy

Definition at line 52 of file PrimaryKiller.hh.

Constructor & Destructor Documentation

PrimaryKiller::PrimaryKiller ( G4String  name,
G4int  depth = 0 
)

Definition at line 56 of file PrimaryKiller.cc.

57 : G4VPrimitiveScorer(name,depth),
59 {
60  fELoss = 0.; // cumulated energy for current event
61 
62  fELossRange_Min = DBL_MAX; // fELoss from which the primary is killed
63  fELossRange_Max = DBL_MAX; // fELoss from which the event is aborted
64  fKineticE_Min = 0; // kinetic energy below which the primary is killed
65 
66  fpELossUI = new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin",this);
67  fpAbortEventIfELossUpperThan =
68  new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMax", this);
69  fpMinKineticE =
70  new G4UIcmdWithADoubleAndUnit("/primaryKiller/minKineticE", this);
71 }
#define DBL_MAX
Definition: templates.hh:83
G4VPrimitiveScorer(G4String name, G4int depth=0)
PrimaryKiller::~PrimaryKiller ( )
virtual

Definition at line 75 of file PrimaryKiller.cc.

76 {;}

Member Function Documentation

void PrimaryKiller::clear ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 158 of file PrimaryKiller.cc.

159 {
160  fELoss = 0.;
161 }

Here is the caller graph for this function:

void PrimaryKiller::DrawAll ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 165 of file PrimaryKiller.cc.

166 {;}
void PrimaryKiller::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 152 of file PrimaryKiller.cc.

153 {
154 }
void PrimaryKiller::Initialize ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 145 of file PrimaryKiller.cc.

146 {
147  clear();
148 }
virtual void clear()

Here is the call graph for this function:

void PrimaryKiller::PrintAll ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 170 of file PrimaryKiller.cc.

171 {}
G4bool PrimaryKiller::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 94 of file PrimaryKiller.cc.

95 {
96  const G4Track* track = aStep->GetTrack();
97 
98  if(track->GetTrackID() != 1) return FALSE;
99 
100  //-------------------
101 
102  double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy();
103 
104  G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy()
105  - kineticE;
106 
107  if ( eLoss == 0. ) return FALSE;
108 
109  //-------------------
110 
111  fELoss+=eLoss;
112 
113  if(fELoss > fELossRange_Max){
115  int eventID =
117  G4cout << " * PrimaryKiller: aborts event " << eventID <<", energy loss "
118  "is too large. \n"
119  << " * Energy loss by primary is: "
120  << G4BestUnit(fELoss, "Energy")
121  << ". Event is aborted if the Eloss is > "
122  << G4BestUnit(fELossRange_Max, "Energy")
123  << G4endl;
124 
125  }
126 
127  if(fELoss >= fELossRange_Min || kineticE <= fKineticE_Min){
128  ((G4Track*)track)->SetTrackStatus(fStopAndKill);
129 // G4cout << "kill track at : "
130 // << G4BestUnit(kineticE, "Energy")
131 // << ", E loss is: "
132 // << G4BestUnit(fELoss, "Energy")
133 // << " /fELossMax: "
134 // << G4BestUnit(fELossMax, "Energy")
135 // << ", EThreshold is: "
136 // << G4BestUnit(fEThreshold, "Energy")
137 // << G4endl;
138  }
139 
140  return TRUE;
141 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4int GetEventID() const
Definition: G4Event.hh:151
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
#define FALSE
Definition: globals.hh:52
#define TRUE
Definition: globals.hh:55
G4int GetTrackID() const
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
G4StepPoint * GetPostStepPoint() const
static G4EventManager * GetEventManager()
#define G4endl
Definition: G4ios.hh:61
virtual void AbortEvent()
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
const G4Event * GetConstCurrentEvent()
G4Track * GetTrack() const

Here is the call graph for this function:

void PrimaryKiller::SetEnergyThreshold ( double  energy)
inline

Set energy under which the particle should be killed

Definition at line 73 of file PrimaryKiller.hh.

73  {
74  fKineticE_Min = energy;
75  }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void PrimaryKiller::SetMaxLossEnergyLimit ( double  energy)
inline

Set the energy loss from which the event is aborted

Definition at line 85 of file PrimaryKiller.hh.

85  {
86  fELossRange_Max = energy;
87  }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void PrimaryKiller::SetMinLossEnergyLimit ( double  energy)
inline

Set the energy loss from which the primary is killed

Definition at line 79 of file PrimaryKiller.hh.

79  {
80  fELossRange_Min = energy;
81  }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void PrimaryKiller::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Method related to G4UImessenger used to control energy cuts through macro file

Reimplemented from G4UImessenger.

Definition at line 80 of file PrimaryKiller.cc.

82 {
83  if(command == fpELossUI){
84  fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue);
85  }
86  else if(command == fpAbortEventIfELossUpperThan){
87  fELossRange_Max =
88  fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue);
89  }
90 }
static G4double GetNewDoubleValue(const char *paramString)

Here is the call graph for this function:


The documentation for this class was generated from the following files: