Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-01 07:51:07

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 SplitProcess.cc
0027 /// \brief Implementation of the SplitProcess class
0028 ///
0029 /// Implementation of the variance reduction particle split class
0030 
0031 // This example is provided by the Geant4-DNA collaboration
0032 // Any report or published results obtained using the Geant4-DNA software
0033 // shall cite the following Geant4-DNA collaboration publication:
0034 // Med. Phys. 37 (2010) 4692-4708
0035 // The Geant4-DNA web site is available at http://geant4-dna.org
0036 //
0037 
0038 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0039 
0040 #include "SplitProcess.hh"
0041 
0042 #include "UserTrackInformation.hh"
0043 
0044 #include "G4LogicalVolume.hh"
0045 #include "G4Region.hh"
0046 #include "G4RegionStore.hh"
0047 #include "G4TouchableHandle.hh"
0048 #include "G4Track.hh"
0049 #include "G4VParticleChange.hh"
0050 
0051 #include <vector>
0052 
0053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0054 
0055 SplitProcess::SplitProcess(G4String regName, G4int nsplit) : fRegionName(regName), fNSplit(nsplit)
0056 {
0057   fRegion = G4RegionStore::GetInstance()->FindOrCreateRegion(fRegionName);
0058 }
0059 
0060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0061 
0062 SplitProcess::~SplitProcess() {}
0063 
0064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0065 
0066 G4VParticleChange* SplitProcess::PostStepDoIt(const G4Track& track, const G4Step& step)
0067 {
0068   G4VParticleChange* particleChange(0);
0069 
0070   if (step.GetPreStepPoint()->GetTouchableHandle()->GetVolume()->GetLogicalVolume()->GetRegion()
0071       != fRegion)
0072   {
0073     particleChange = pRegProcess->PostStepDoIt(track, step);
0074     assert(0 != particleChange);
0075     return particleChange;
0076   }
0077 
0078   if (fNSplit == 1) {
0079     particleChange = pRegProcess->PostStepDoIt(track, step);
0080     assert(0 != particleChange);
0081     return particleChange;
0082   }
0083 
0084   UserTrackInformation* parentInformation =
0085     (UserTrackInformation*)(step.GetTrack()->GetUserInformation());
0086   G4int initialSplitTrackID = parentInformation->GetSplitTrackID();
0087 
0088   if (initialSplitTrackID > 1) {
0089     particleChange = pRegProcess->PostStepDoIt(track, step);
0090     assert(0 != particleChange);
0091     return particleChange;
0092   }
0093 
0094   G4double weight = track.GetWeight() / fNSplit;
0095   G4int splitTrackID = 3;
0096 
0097   std::vector<G4Track*> secondaries;
0098   std::vector<G4int> vSplitTrack;
0099 
0100   for (int i = 0; i < fNSplit; i++) {
0101     particleChange = pRegProcess->PostStepDoIt(track, step);
0102     assert(0 != particleChange);
0103     particleChange->SetVerboseLevel(0);
0104     G4Track* newTrack = new G4Track(*(particleChange->GetSecondary(0)));
0105 
0106     secondaries.push_back(newTrack);
0107     vSplitTrack.push_back(splitTrackID);
0108 
0109     splitTrackID++;
0110   }
0111 
0112   parentInformation->SetSplitTrackID(2);
0113 
0114   particleChange->SetNumberOfSecondaries(secondaries.size());
0115   particleChange->SetSecondaryWeightByProcess(true);
0116 
0117   std::vector<G4Track*>::iterator iter = secondaries.begin();
0118   G4int i = 0;
0119   while (iter != secondaries.end()) {
0120     G4Track* newTrack = *iter;
0121     newTrack->SetWeight(weight);
0122 
0123     UserTrackInformation* secondaryInformation = new UserTrackInformation();
0124     secondaryInformation->SetSplitTrackID(vSplitTrack[i]);
0125     newTrack->SetUserInformation(secondaryInformation);
0126 
0127     particleChange->AddSecondary(newTrack);
0128 
0129     iter++;
0130     i++;
0131   }
0132 
0133   return particleChange;
0134 }
0135 
0136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......