Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
SplitProcess.cc
이 파일의 문서화 페이지로 가기
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 // $ID$
35 
36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
37 
38 #include "SplitProcess.hh"
39 #include "UserTrackInformation.hh"
40 
41 #include "G4Track.hh"
42 #include "G4VParticleChange.hh"
43 #include "G4LogicalVolume.hh"
44 #include "G4TouchableHandle.hh"
45 #include "G4Region.hh"
46 #include "G4RegionStore.hh"
47 
48 #include <vector>
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
53 :fRegionName(regName), fNSplit(nsplit)
54 {
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61 }
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
66  G4VParticleChange* particleChange(0);
67 
69  ->GetRegion() != fRegion ) {
70  particleChange = pRegProcess->PostStepDoIt(track, step);
71  assert ( 0 != particleChange);
72  return particleChange;
73  }
74 
75  if ( fNSplit == 1 ) {
76  particleChange = pRegProcess->PostStepDoIt(track, step);
77  assert ( 0 != particleChange);
78  return particleChange;
79  }
80 
81  UserTrackInformation* parentInformation =
83  G4int initialSplitTrackID = parentInformation->GetSplitTrackID();
84 
85  if ( initialSplitTrackID > 1 ) {
86  particleChange = pRegProcess->PostStepDoIt(track, step);
87  assert ( 0 != particleChange);
88  return particleChange;
89  }
90 
91  G4double weight = track.GetWeight()/fNSplit;
92  G4int splitTrackID = 3;
93 
94  std::vector<G4Track*> secondaries;
95  std::vector<G4int> vSplitTrack;
96 
97  for ( int i = 0; i < fNSplit; i++ ) {
98  particleChange = pRegProcess->PostStepDoIt(track, step);
99  assert( 0 != particleChange);
100  particleChange->SetVerboseLevel(0);
101  G4Track* newTrack = new G4Track(*(particleChange->GetSecondary(0)));
102 
103  secondaries.push_back( newTrack );
104  vSplitTrack.push_back( splitTrackID );
105 
106  splitTrackID++;
107  }
108 
109  parentInformation->SetSplitTrackID(2);
110 
111  particleChange->SetNumberOfSecondaries(secondaries.size());
112  particleChange->SetSecondaryWeightByProcess(true);
113 
114  std::vector<G4Track*>::iterator iter = secondaries.begin();
115  G4int i = 0;
116  while( iter != secondaries.end() ) {
117  G4Track* newTrack = *iter;
118  newTrack->SetWeight(weight);
119 
120  UserTrackInformation* secondaryInformation = new UserTrackInformation();
121  secondaryInformation->SetSplitTrackID(vSplitTrack[i]);
122  newTrack->SetUserInformation(secondaryInformation);
123 
124  particleChange->AddSecondary(newTrack);
125 
126  iter++;
127  i++;
128  }
129 
130  return particleChange;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void AddSecondary(G4Track *aSecondary)
Definition of the SplitProcess class.
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
Definition: SplitProcess.cc:65
virtual ~SplitProcess()
Definition: SplitProcess.cc:60
G4VProcess * pRegProcess
G4LogicalVolume * GetLogicalVolume() const
void SetWeight(G4double aValue)
G4StepPoint * GetPreStepPoint() const
SplitProcess(G4String regName, G4int split)
Definition: SplitProcess.cc:52
void SetSecondaryWeightByProcess(G4bool)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void SetVerboseLevel(G4int vLevel)
const G4TouchableHandle & GetTouchableHandle() const
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
Definition: G4Step.hh:76
Definition of the UserTrackInformation class.
double weight
Definition: plottest35.C:25
static G4RegionStore * GetInstance()
int G4int
Definition: G4Types.hh:78
G4Track * GetSecondary(G4int anIndex) const
G4VUserTrackInformation * GetUserInformation() const
G4String fRegionName
Definition: SplitProcess.hh:48
G4Region * fRegion
Definition: SplitProcess.hh:47
void SetSplitTrackID(G4int splitTrackID)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4Region * GetRegion() const
void SetUserInformation(G4VUserTrackInformation *aValue) const
G4Region * FindOrCreateRegion(const G4String &name)
void SetNumberOfSecondaries(G4int totSecondaries)