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 | // File name: G4PiMinusAbsorptionAtRest.hh |
---|
27 | // |
---|
28 | // Author: Maria Grazia Pia (pia@genova.infn.it) |
---|
29 | // |
---|
30 | // Creation date: 8 May 1998 |
---|
31 | // |
---|
32 | // Modifications: |
---|
33 | // MGP 4 Jul 1998 Changed excitation energy calculation |
---|
34 | // MGP 14 Sep 1998 Fixed excitation energy calculation |
---|
35 | // |
---|
36 | // ------------------------------------------------------------------- |
---|
37 | |
---|
38 | #include "G4ios.hh" |
---|
39 | |
---|
40 | #include "G4PiMinusAbsorptionAtRest.hh" |
---|
41 | |
---|
42 | #include "G4PiMinusStopLi.hh" |
---|
43 | #include "G4PiMinusStopC.hh" |
---|
44 | #include "G4PiMinusStopN.hh" |
---|
45 | #include "G4PiMinusStopO.hh" |
---|
46 | #include "G4PiMinusStopAl.hh" |
---|
47 | #include "G4PiMinusStopCu.hh" |
---|
48 | #include "G4PiMinusStopCo.hh" |
---|
49 | #include "G4PiMinusStopTa.hh" |
---|
50 | #include "G4PiMinusStopPb.hh" |
---|
51 | #include "G4StopTheoDeexcitation.hh" |
---|
52 | #include "G4StopDummyDeexcitation.hh" |
---|
53 | #include "G4DynamicParticle.hh" |
---|
54 | #include "G4DynamicParticleVector.hh" |
---|
55 | #include "G4NucleiPropertiesTable.hh" |
---|
56 | #include "Randomize.hh" |
---|
57 | #include "G4ThreeVector.hh" |
---|
58 | #include "G4LorentzVector.hh" |
---|
59 | |
---|
60 | // Constructor |
---|
61 | |
---|
62 | G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest(const G4String& processName, |
---|
63 | G4ProcessType aType) : |
---|
64 | G4VRestProcess (processName, aType) |
---|
65 | { |
---|
66 | SetProcessSubType(fHadronAtRest); |
---|
67 | |
---|
68 | _indexDeexcitation = 0; |
---|
69 | |
---|
70 | if (verboseLevel>0) |
---|
71 | { G4cout << GetProcessName() << " is created "<< G4endl; } |
---|
72 | } |
---|
73 | |
---|
74 | |
---|
75 | // Destructor |
---|
76 | |
---|
77 | G4PiMinusAbsorptionAtRest::~G4PiMinusAbsorptionAtRest() |
---|
78 | {} |
---|
79 | |
---|
80 | |
---|
81 | G4VParticleChange* G4PiMinusAbsorptionAtRest::AtRestDoIt(const G4Track& track, const G4Step& ) |
---|
82 | { |
---|
83 | const G4DynamicParticle* stoppedHadron = track.GetDynamicParticle(); |
---|
84 | |
---|
85 | // Check applicability |
---|
86 | if (! IsApplicable(*(stoppedHadron->GetDefinition()))) |
---|
87 | { |
---|
88 | G4cerr |
---|
89 | << "G4PiMinusAbsorptionAtRest: ERROR, particle must be a pion minus!" |
---|
90 | << G4endl; |
---|
91 | return NULL; |
---|
92 | } |
---|
93 | |
---|
94 | // Get the current material |
---|
95 | const G4Material* material = track.GetMaterial(); |
---|
96 | |
---|
97 | G4double A=-1; |
---|
98 | G4double Z=-1; |
---|
99 | G4double random = G4UniformRand(); |
---|
100 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
101 | unsigned int i; |
---|
102 | G4double sum = 0; |
---|
103 | G4double totalsum=0; |
---|
104 | for(i=0; i<material->GetNumberOfElements(); ++i) |
---|
105 | { |
---|
106 | if((*theElementVector)[i]->GetZ()!=1) totalsum+=material->GetFractionVector()[i]; |
---|
107 | } |
---|
108 | for (i = 0; i<material->GetNumberOfElements(); ++i) |
---|
109 | { |
---|
110 | if((*theElementVector)[i]->GetZ()!=1) sum += material->GetFractionVector()[i]; |
---|
111 | if ( sum/totalsum > random ) |
---|
112 | { |
---|
113 | A = (*theElementVector)[i]->GetA()*mole/g; |
---|
114 | Z = (*theElementVector)[i]->GetZ(); |
---|
115 | break; |
---|
116 | } |
---|
117 | } |
---|
118 | |
---|
119 | // Do the interaction with the nucleon cluster |
---|
120 | |
---|
121 | G4PiMinusStopMaterial* algorithm = LoadAlgorithm(static_cast<G4int>(Z)); |
---|
122 | G4PiMinusStopAbsorption stopAbsorption(algorithm,Z,A); |
---|
123 | stopAbsorption.SetVerboseLevel(verboseLevel); |
---|
124 | |
---|
125 | G4DynamicParticleVector* absorptionProducts = stopAbsorption.DoAbsorption(); |
---|
126 | |
---|
127 | // Deal with the leftover nucleus |
---|
128 | |
---|
129 | G4double pionEnergy = stoppedHadron->GetTotalEnergy(); |
---|
130 | G4double excitation = pionEnergy - stopAbsorption.Energy(); |
---|
131 | if (excitation < 0.) |
---|
132 | { |
---|
133 | G4Exception("G4PiMinusAbsorptionAtRest", "007", FatalException, |
---|
134 | "AtRestDoIt -- excitation energy < 0"); |
---|
135 | } |
---|
136 | if (verboseLevel>0) { G4cout << " excitation " << excitation << G4endl; } |
---|
137 | |
---|
138 | G4StopDeexcitationAlgorithm* nucleusAlgorithm = LoadNucleusAlgorithm(); |
---|
139 | G4StopDeexcitation stopDeexcitation(nucleusAlgorithm); |
---|
140 | |
---|
141 | G4double newZ = Z - stopAbsorption.NProtons(); |
---|
142 | G4double newN = A - Z - stopAbsorption.NNeutrons(); |
---|
143 | G4double newA = newZ + newN; |
---|
144 | G4double pNucleus = (stopAbsorption.RecoilMomentum()).mag(); |
---|
145 | G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,excitation,pNucleus); |
---|
146 | |
---|
147 | unsigned int nAbsorptionProducts = 0; |
---|
148 | if (absorptionProducts != 0) |
---|
149 | { nAbsorptionProducts = absorptionProducts->size(); } |
---|
150 | |
---|
151 | unsigned int nFragmentationProducts = 0; |
---|
152 | if (fragmentationProducts != 0) |
---|
153 | { nFragmentationProducts = fragmentationProducts->size(); } |
---|
154 | |
---|
155 | if (verboseLevel>0) |
---|
156 | { |
---|
157 | G4cout << "nAbsorptionProducts = " << nAbsorptionProducts |
---|
158 | << " nFragmentationProducts = " << nFragmentationProducts |
---|
159 | << G4endl; |
---|
160 | } |
---|
161 | |
---|
162 | // Deal with ParticleChange final state |
---|
163 | |
---|
164 | aParticleChange.Initialize(track); |
---|
165 | aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts + nFragmentationProducts)); |
---|
166 | |
---|
167 | for (i = 0; i<nAbsorptionProducts; i++) |
---|
168 | { aParticleChange.AddSecondary((*absorptionProducts)[i]); } |
---|
169 | |
---|
170 | // for (i = 0; i<nFragmentationProducts; i++) |
---|
171 | // { aParticleChange.AddSecondary(fragmentationProducts->at(i)); } |
---|
172 | for(i=0; i<nFragmentationProducts; i++) |
---|
173 | { |
---|
174 | G4DynamicParticle * aNew = |
---|
175 | new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(), |
---|
176 | (*fragmentationProducts)[i]->GetMomentum()); |
---|
177 | G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime()); |
---|
178 | aParticleChange.AddSecondary(aNew, newTime); |
---|
179 | delete (*fragmentationProducts)[i]; |
---|
180 | } |
---|
181 | |
---|
182 | if (fragmentationProducts != 0) delete fragmentationProducts; |
---|
183 | |
---|
184 | if (_indexDeexcitation == 1) aParticleChange.ProposeLocalEnergyDeposit(excitation); |
---|
185 | |
---|
186 | // Kill the absorbed pion |
---|
187 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
188 | |
---|
189 | return &aParticleChange; |
---|
190 | |
---|
191 | } |
---|
192 | |
---|
193 | G4PiMinusStopMaterial* G4PiMinusAbsorptionAtRest::LoadAlgorithm(int Z) |
---|
194 | { |
---|
195 | if (verboseLevel>0) |
---|
196 | { |
---|
197 | G4cout << "Load material algorithm " << Z << G4endl; |
---|
198 | } |
---|
199 | |
---|
200 | G4int index = 3; |
---|
201 | if (Z <= 3) { index = 3;} |
---|
202 | if (Z > 3 && Z<= 6) {index = 6;} |
---|
203 | if (Z == 7) {index = 7;} |
---|
204 | if (Z >= 8 && Z<= 11) {index = 8;} |
---|
205 | if (Z >= 12 && Z<= 18) {index = 13;} |
---|
206 | if (Z >=19 && Z<= 27) {index = 27;} |
---|
207 | if (Z >= 28 && Z<= 51) {index = 29;} |
---|
208 | if (Z >=52 ) {index = 73;} |
---|
209 | |
---|
210 | switch (index) |
---|
211 | { |
---|
212 | case 3: |
---|
213 | if (verboseLevel>0) |
---|
214 | { G4cout << " =================== Load Li algorithm " << G4endl; } |
---|
215 | return new G4PiMinusStopLi(); |
---|
216 | case 6: |
---|
217 | if (verboseLevel>0) |
---|
218 | { G4cout << " =================== Load C algorithm " << G4endl; } |
---|
219 | return new G4PiMinusStopC(); |
---|
220 | case 7: |
---|
221 | if (verboseLevel>0) |
---|
222 | { G4cout << " =================== Load N algorithm " << G4endl; } |
---|
223 | return new G4PiMinusStopN(); |
---|
224 | case 8: |
---|
225 | if (verboseLevel>0) |
---|
226 | { G4cout << " =================== Load O algorithm " << G4endl; } |
---|
227 | return new G4PiMinusStopO(); |
---|
228 | case 13: |
---|
229 | if (verboseLevel>0) |
---|
230 | { G4cout << " =================== Load Al algorithm " << G4endl; } |
---|
231 | return new G4PiMinusStopAl(); |
---|
232 | case 27: |
---|
233 | if (verboseLevel>0) |
---|
234 | { G4cout << " =================== Load Cu algorithm " << G4endl; } |
---|
235 | return new G4PiMinusStopCu(); |
---|
236 | case 29: |
---|
237 | if (verboseLevel>0) |
---|
238 | { G4cout << " =================== Load Co algorithm " << G4endl; } |
---|
239 | return new G4PiMinusStopCo(); |
---|
240 | case 73: |
---|
241 | if (verboseLevel>0) |
---|
242 | { G4cout << " =================== Load Ta algorithm " << G4endl; } |
---|
243 | return new G4PiMinusStopTa(); |
---|
244 | default: |
---|
245 | if (verboseLevel>0) |
---|
246 | { G4cout << " =================== Load default material algorithm " << G4endl; } |
---|
247 | return new G4PiMinusStopC(); |
---|
248 | } |
---|
249 | } |
---|
250 | |
---|
251 | G4StopDeexcitationAlgorithm* G4PiMinusAbsorptionAtRest::LoadNucleusAlgorithm() |
---|
252 | { |
---|
253 | |
---|
254 | switch (_indexDeexcitation) |
---|
255 | { |
---|
256 | case 0: |
---|
257 | if (verboseLevel>0) |
---|
258 | { G4cout << " =================== Load Theo deexcitation " << G4endl; } |
---|
259 | return new G4StopTheoDeexcitation(); |
---|
260 | case 1: |
---|
261 | if (verboseLevel>0) |
---|
262 | { G4cout << " =================== Load Dummy deexcitation " << G4endl; } |
---|
263 | return new G4StopDummyDeexcitation(); |
---|
264 | default: |
---|
265 | if (verboseLevel>0) |
---|
266 | { G4cout << " =================== Load default deexcitation " << G4endl; } |
---|
267 | return new G4StopTheoDeexcitation(); |
---|
268 | } |
---|
269 | } |
---|
270 | |
---|
271 | void G4PiMinusAbsorptionAtRest::SetDeexcitationAlgorithm(G4int index) |
---|
272 | { |
---|
273 | _indexDeexcitation = index; |
---|
274 | } |
---|