File indexing completed on 2025-02-23 09:22:09
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037 #include "SplitProcess.hh"
0038
0039 #include "UserTrackInformation.hh"
0040
0041 #include "G4LogicalVolume.hh"
0042 #include "G4Region.hh"
0043 #include "G4RegionStore.hh"
0044 #include "G4TouchableHandle.hh"
0045 #include "G4Track.hh"
0046 #include "G4VParticleChange.hh"
0047
0048 #include <vector>
0049
0050
0051
0052 SplitProcess::SplitProcess(G4String regName, G4int nsplit) : fRegionName(regName), fNSplit(nsplit)
0053 {
0054 fRegion = G4RegionStore::GetInstance()->FindOrCreateRegion(fRegionName);
0055 }
0056
0057
0058
0059 SplitProcess::~SplitProcess() {}
0060
0061
0062
0063 G4VParticleChange* SplitProcess::PostStepDoIt(const G4Track& track, const G4Step& step)
0064 {
0065 G4VParticleChange* particleChange(0);
0066
0067 if (step.GetPreStepPoint()->GetTouchableHandle()->GetVolume()->GetLogicalVolume()->GetRegion()
0068 != fRegion)
0069 {
0070 particleChange = pRegProcess->PostStepDoIt(track, step);
0071 assert(0 != particleChange);
0072 return particleChange;
0073 }
0074
0075 if (fNSplit == 1) {
0076 particleChange = pRegProcess->PostStepDoIt(track, step);
0077 assert(0 != particleChange);
0078 return particleChange;
0079 }
0080
0081 UserTrackInformation* parentInformation =
0082 (UserTrackInformation*)(step.GetTrack()->GetUserInformation());
0083 G4int initialSplitTrackID = parentInformation->GetSplitTrackID();
0084
0085 if (initialSplitTrackID > 1) {
0086 particleChange = pRegProcess->PostStepDoIt(track, step);
0087 assert(0 != particleChange);
0088 return particleChange;
0089 }
0090
0091 G4double weight = track.GetWeight() / fNSplit;
0092 G4int splitTrackID = 3;
0093
0094 std::vector<G4Track*> secondaries;
0095 std::vector<G4int> vSplitTrack;
0096
0097 for (int i = 0; i < fNSplit; i++) {
0098 particleChange = pRegProcess->PostStepDoIt(track, step);
0099 assert(0 != particleChange);
0100 particleChange->SetVerboseLevel(0);
0101 G4Track* newTrack = new G4Track(*(particleChange->GetSecondary(0)));
0102
0103 secondaries.push_back(newTrack);
0104 vSplitTrack.push_back(splitTrackID);
0105
0106 splitTrackID++;
0107 }
0108
0109 parentInformation->SetSplitTrackID(2);
0110
0111 particleChange->SetNumberOfSecondaries(secondaries.size());
0112 particleChange->SetSecondaryWeightByProcess(true);
0113
0114 std::vector<G4Track*>::iterator iter = secondaries.begin();
0115 G4int i = 0;
0116 while (iter != secondaries.end()) {
0117 G4Track* newTrack = *iter;
0118 newTrack->SetWeight(weight);
0119
0120 UserTrackInformation* secondaryInformation = new UserTrackInformation();
0121 secondaryInformation->SetSplitTrackID(vSplitTrack[i]);
0122 newTrack->SetUserInformation(secondaryInformation);
0123
0124 particleChange->AddSecondary(newTrack);
0125
0126 iter++;
0127 i++;
0128 }
0129
0130 return particleChange;
0131 }
0132
0133