Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:50

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // This example is provided by the Geant4-DNA collaboration
0027 // chem6 example is derived from chem4 and chem5 examples
0028 //
0029 // Any report or published results obtained using the Geant4-DNA software
0030 // shall cite the following Geant4-DNA collaboration publication:
0031 // J. Appl. Phys. 125 (2019) 104301
0032 // Med. Phys. 45 (2018) e722-e739
0033 // J. Comput. Phys. 274 (2014) 841-882
0034 // Med. Phys. 37 (2010) 4692-4708
0035 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
0036 // The Geant4-DNA web site is available at http://geant4-dna.org
0037 //
0038 // Authors: W. G. Shin and S. Incerti (CENBG, France)
0039 //
0040 // $Id$
0041 //
0042 /// \file PrimaryKiller.cc
0043 /// \brief Implementation of the PrimaryKiller class
0044 
0045 #include "PrimaryKiller.hh"
0046 
0047 #include <G4Event.hh>
0048 #include <G4RunManager.hh>
0049 #include <G4SystemOfUnits.hh>
0050 #include <G4UIcmdWith3VectorAndUnit.hh>
0051 #include <G4UIcmdWithADoubleAndUnit.hh>
0052 #include <G4UnitsTable.hh>
0053 
0054 /** \file PrimaryKiller.cc
0055     \class PrimaryKiller
0056 
0057     Kill the primary particle:
0058     - either after a given energy loss
0059     - or after the primary particle has reached a given energy
0060  */
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0063 
0064 PrimaryKiller::PrimaryKiller(G4String name, G4int depth)
0065   : G4VPrimitiveScorer(name, depth), G4UImessenger()
0066 {
0067   fELoss = 0.;  // cumulated energy for current event
0068 
0069   fELossRange_Min = DBL_MAX;  // fELoss from which the primary is killed
0070   fELossRange_Max = DBL_MAX;  // fELoss from which the event is aborted
0071   fKineticE_Min = 0;  // kinetic energy below which the primary is killed
0072   fPhantomSize = G4ThreeVector(1 * km, 1 * km, 1 * km);
0073 
0074   fpELossUI = new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin", this);
0075   fpAbortEventIfELossUpperThan = new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMax", this);
0076   fpMinKineticE = new G4UIcmdWithADoubleAndUnit("/primaryKiller/minKineticE", this);
0077   fpSizeUI = new G4UIcmdWith3VectorAndUnit("/primaryKiller/setSize", this);
0078   fpSizeUI->SetDefaultUnit("um");
0079 }
0080 
0081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0082 
0083 PrimaryKiller::~PrimaryKiller()
0084 {
0085   delete fpELossUI;
0086   delete fpAbortEventIfELossUpperThan;
0087   delete fpSizeUI;
0088 }
0089 
0090 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0091 
0092 void PrimaryKiller::SetNewValue(G4UIcommand* command, G4String newValue)
0093 {
0094   if (command == fpELossUI) {
0095     fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue);
0096   }
0097   else if (command == fpAbortEventIfELossUpperThan) {
0098     fELossRange_Max = fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue);
0099   }
0100   else if (command == fpSizeUI) {
0101     fPhantomSize = fpSizeUI->GetNew3VectorValue(newValue);
0102   }
0103 }
0104 
0105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0106 
0107 G4bool PrimaryKiller::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0108 {
0109   const G4Track* track = aStep->GetTrack();
0110   G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
0111 
0112   if (std::abs(pos.x()) > fPhantomSize.getX() / 2 || std::abs(pos.y()) > fPhantomSize.getY() / 2
0113       || std::abs(pos.z()) > fPhantomSize.getZ() / 2)
0114   {
0115     ((G4Track*)track)->SetTrackStatus(fStopAndKill);
0116     return false;
0117   }
0118 
0119   if (track->GetTrackID() != 1 || track->GetParticleDefinition()->GetPDGEncoding() != 11)
0120     return FALSE;
0121 
0122   //-------------------
0123 
0124   double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy();
0125 
0126   G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy() - kineticE;
0127 
0128   if (eLoss == 0.) return FALSE;
0129 
0130   //-------------------
0131 
0132   fELoss += eLoss;
0133 
0134   if (fELoss > fELossRange_Max) {
0135     G4RunManager::GetRunManager()->AbortEvent();
0136     /*    int eventID =
0137          G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID();
0138 
0139         G4cout << " * PrimaryKiller: aborts event " << eventID <<" energy loss "
0140                   "is too large. \n"
0141                << " * Energy loss by primary is: "
0142                << G4BestUnit(fELoss, "Energy")
0143                << ". Event is aborted if the Eloss is > "
0144                << G4BestUnit(fELossRange_Max, "Energy")
0145                << G4endl;
0146      */
0147   }
0148 
0149   if (fELoss >= fELossRange_Min || kineticE <= fKineticE_Min) {
0150     ((G4Track*)track)->SetTrackStatus(fStopAndKill);
0151     //     G4cout << "kill track at : "<<'\n';
0152     //           << G4BestUnit(kineticE, "Energy")
0153     //           << ", E loss is: "
0154     //           << G4BestUnit(fELoss, "Energy")
0155     //           << " /fELossMax: "
0156     //           << G4BestUnit(fELossMax, "Energy")
0157     //           << ", EThreshold is: "
0158     //           << G4BestUnit(fEThreshold, "Energy")
0159     //           << G4endl;
0160   }
0161 
0162   return TRUE;
0163 }
0164 
0165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0166 
0167 void PrimaryKiller::Initialize(G4HCofThisEvent* /*HCE*/)
0168 {
0169   fELoss = 0.;
0170 }
0171 
0172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0173 
0174 void PrimaryKiller::EndOfEvent(G4HCofThisEvent*) {}
0175 
0176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0177 
0178 void PrimaryKiller::DrawAll()
0179 {
0180   ;
0181 }
0182 
0183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0184 
0185 void PrimaryKiller::PrintAll() {}
0186 
0187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....