Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-03-29 07:51:37

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 /// \file PrimaryKiller.cc
0027 /// \brief Implementation of the PrimaryKiller class
0028 
0029 // This example is provided by the Geant4-DNA collaboration
0030 // chem6 example is derived from chem4 and chem5 examples
0031 //
0032 // Any report or published results obtained using the Geant4-DNA software
0033 // shall cite the following Geant4-DNA collaboration publication:
0034 // J. Appl. Phys. 125 (2019) 104301
0035 // Med. Phys. 45 (2018) e722-e739
0036 // J. Comput. Phys. 274 (2014) 841-882
0037 // Med. Phys. 37 (2010) 4692-4708
0038 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
0039 // The Geant4-DNA web site is available at http://geant4-dna.org
0040 //
0041 // Authors: W. G. Shin and S. Incerti (CENBG, France)
0042 
0043 #include "PrimaryKiller.hh"
0044 
0045 #include <G4Event.hh>
0046 #include <G4RunManager.hh>
0047 #include <G4UIcmdWith3VectorAndUnit.hh>
0048 #include <G4UIcmdWithADoubleAndUnit.hh>
0049 #include <G4UnitsTable.hh>
0050 
0051 /** \file PrimaryKiller.cc
0052     \class PrimaryKiller
0053 
0054     Kill the primary particle:
0055     - either after a given energy loss
0056     - or after the primary particle has reached a given energy
0057  */
0058 
0059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0060 
0061 PrimaryKiller::PrimaryKiller(const G4String &name, const G4int depth)
0062   : G4VPrimitiveScorer(name, depth), G4UImessenger() {
0063   fpELossUI = std::make_unique<G4UIcmdWithADoubleAndUnit>
0064       ("/primaryKiller/eLossMin", this);
0065   fpAbortEventIfELossUpperThan = std::make_unique<G4UIcmdWithADoubleAndUnit>
0066       ("/primaryKiller/eLossMax", this);
0067   fpMinKineticE = std::make_unique<G4UIcmdWithADoubleAndUnit>
0068       ("/primaryKiller/minKineticE", this);
0069   fpSizeUI = std::make_unique<G4UIcmdWith3VectorAndUnit>
0070       ("/primaryKiller/setSize", this);
0071   fpSizeUI->SetDefaultUnit("um");
0072 }
0073 
0074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0075 
0076 void PrimaryKiller::SetNewValue(G4UIcommand *command, G4String newValue) {
0077   if (command == fpELossUI.get()) {
0078     fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue);
0079   } else if (command == fpAbortEventIfELossUpperThan.get()) {
0080     fELossRange_Max = fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue);
0081   } else if (command == fpSizeUI.get()) {
0082     fPhantomSize = fpSizeUI->GetNew3VectorValue(newValue);
0083   }
0084 }
0085 
0086 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0087 
0088 G4bool PrimaryKiller::ProcessHits(G4Step *aStep, G4TouchableHistory *) {
0089   const G4Track *track = aStep->GetTrack();
0090   G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
0091 
0092   if (std::abs(pos.x()) > fPhantomSize.getX() / 2 || std::abs(pos.y()) > fPhantomSize.getY() / 2
0093       || std::abs(pos.z()) > fPhantomSize.getZ() / 2) {
0094     const_cast<G4Track *>(track)->SetTrackStatus(fStopAndKill);
0095     return false;
0096   }
0097 
0098   if (track->GetTrackID() != 1 || track->GetParticleDefinition()->GetPDGEncoding() != 11) {
0099     return FALSE;
0100   }
0101 
0102   //-------------------
0103 
0104   const G4double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy();
0105 
0106   const G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy() - kineticE;
0107 
0108   if (eLoss == 0.) { return FALSE; }
0109 
0110   //-------------------
0111 
0112   fELoss += eLoss;
0113 
0114   if (fELoss > fELossRange_Max) {
0115     G4RunManager::GetRunManager()->AbortEvent();
0116     /*    int eventID =
0117          G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID();
0118 
0119         G4cout << " * PrimaryKiller: aborts event " << eventID <<" energy loss "
0120                   "is too large. \n"
0121                << " * Energy loss by primary is: "
0122                << G4BestUnit(fELoss, "Energy")
0123                << ". Event is aborted if the Eloss is > "
0124                << G4BestUnit(fELossRange_Max, "Energy")
0125                << G4endl;
0126      */
0127   }
0128 
0129   if (fELoss >= fELossRange_Min || kineticE <= fKineticE_Min) {
0130     const_cast<G4Track *>(track)->SetTrackStatus(fStopAndKill);
0131     //     G4cout << "kill track at : "<<'\n';
0132     //           << G4BestUnit(kineticE, "Energy")
0133     //           << ", E loss is: "
0134     //           << G4BestUnit(fELoss, "Energy")
0135     //           << " /fELossMax: "
0136     //           << G4BestUnit(fELossMax, "Energy")
0137     //           << ", EThreshold is: "
0138     //           << G4BestUnit(fEThreshold, "Energy")
0139     //           << G4endl;
0140   }
0141 
0142   return TRUE;
0143 }
0144 
0145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....