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 | // $Id: G4InclCascadeInterface.cc,v 1.10 2007/12/10 16:32:02 gunter Exp $ |
---|
27 | // Translation of INCL4.2/ABLA V3 |
---|
28 | // Pekka Kaitaniemi, HIP (translation) |
---|
29 | // Christelle Schmidt, IPNL (fission code) |
---|
30 | // Alain Boudard, CEA (contact person INCL/ABLA) |
---|
31 | // Aatos Heikkinen, HIP (project coordination) |
---|
32 | |
---|
33 | //#define DEBUGINCL 1 |
---|
34 | |
---|
35 | #include "G4InclCascadeInterface.hh" |
---|
36 | #include "math.h" |
---|
37 | #include "G4GenericIon.hh" |
---|
38 | #include "CLHEP/Random/Random.h" |
---|
39 | |
---|
40 | G4InclCascadeInterface::G4InclCascadeInterface() |
---|
41 | { |
---|
42 | hazard = new G4Hazard(); |
---|
43 | const G4long* table_entry = CLHEP::HepRandom::getTheSeeds(); // Get random seed from CLHEP. |
---|
44 | hazard->ial = (*table_entry); |
---|
45 | |
---|
46 | varntp = new G4VarNtp(); |
---|
47 | calincl = new G4Calincl(); |
---|
48 | ws = new G4Ws(); |
---|
49 | mat = new G4Mat(); |
---|
50 | incl = new G4Incl(hazard, calincl, ws, mat, varntp); |
---|
51 | |
---|
52 | verboseLevel = 0; |
---|
53 | } |
---|
54 | |
---|
55 | G4InclCascadeInterface::~G4InclCascadeInterface() |
---|
56 | { |
---|
57 | delete hazard; |
---|
58 | delete varntp; |
---|
59 | delete calincl; |
---|
60 | delete ws; |
---|
61 | delete mat; |
---|
62 | delete incl; |
---|
63 | } |
---|
64 | |
---|
65 | G4HadFinalState* G4InclCascadeInterface::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus) |
---|
66 | { |
---|
67 | G4int maxTries = 200; |
---|
68 | |
---|
69 | G4int particleI, n = 0; |
---|
70 | |
---|
71 | G4int bulletType = 0; |
---|
72 | |
---|
73 | // Print diagnostic messages: 0 = silent, 1 and 2 = verbose |
---|
74 | verboseLevel = 0; |
---|
75 | |
---|
76 | // Increase the event number: |
---|
77 | eventNumber++; |
---|
78 | |
---|
79 | if (verboseLevel > 1) { |
---|
80 | G4cout << " >>> G4InclCascadeInterface::ApplyYourself called" << G4endl; |
---|
81 | } |
---|
82 | |
---|
83 | if(verboseLevel > 1) { |
---|
84 | G4cout <<"G4InclCascadeInterface: Now processing INCL4 event number:" << eventNumber << G4endl; |
---|
85 | } |
---|
86 | |
---|
87 | // INCL4 needs the energy in units MeV |
---|
88 | G4double bulletE = aTrack.GetKineticEnergy() * MeV; |
---|
89 | |
---|
90 | #ifdef DEBUGINCL |
---|
91 | G4cout <<"Bullet energy = " << bulletE / MeV << G4endl; |
---|
92 | #endif |
---|
93 | |
---|
94 | G4double targetA = theNucleus.GetN(); |
---|
95 | G4double targetZ = theNucleus.GetZ(); |
---|
96 | |
---|
97 | G4double eKin; |
---|
98 | G4double momx = 0.0, momy = 0.0, momz = 0.0; |
---|
99 | G4DynamicParticle *cascadeParticle = 0; |
---|
100 | G4ParticleDefinition *aParticleDefinition = 0; |
---|
101 | |
---|
102 | // INCL assumes the projectile particle is going in the direction of |
---|
103 | // the Z-axis. Here we construct proper rotation to convert the |
---|
104 | // momentum vectors of the outcoming particles to the original |
---|
105 | // coordinate system. |
---|
106 | G4LorentzVector projectileMomentum = aTrack.Get4Momentum(); |
---|
107 | G4LorentzRotation toZ; |
---|
108 | toZ.rotateZ(-projectileMomentum.phi()); |
---|
109 | toZ.rotateY(-projectileMomentum.theta()); |
---|
110 | G4LorentzRotation toLabFrame = toZ.inverse(); |
---|
111 | |
---|
112 | theResult.Clear(); // Make sure the output data structure is clean. |
---|
113 | |
---|
114 | // Map Geant4 particle types to corresponding INCL4 types. |
---|
115 | enum bulletParticleType {nucleus = 0, proton = 1, neutron = 2, pionPlus = 3, pionZero = 4, |
---|
116 | pionMinus = 5, deuteron = 6, triton = 7, he3 = 8, he4 = 9}; |
---|
117 | |
---|
118 | // Coding particles for use with INCL4 and ABLA |
---|
119 | if (aTrack.GetDefinition() == G4Proton::Proton() ) bulletType = proton; |
---|
120 | if (aTrack.GetDefinition() == G4Neutron::Neutron() ) bulletType = neutron; |
---|
121 | if (aTrack.GetDefinition() == G4PionPlus::PionPlus() ) bulletType = pionPlus; |
---|
122 | if (aTrack.GetDefinition() == G4PionMinus::PionMinus() ) bulletType = pionMinus; |
---|
123 | if (aTrack.GetDefinition() == G4PionZero::PionZero() ) bulletType = pionZero; |
---|
124 | |
---|
125 | #ifdef DEBUGINCL |
---|
126 | G4int baryonBullet = 0, chargeBullet = 0; |
---|
127 | if(bulletType == proton || bulletType == neutron) baryonBullet = 1; |
---|
128 | if(bulletType == proton || bulletType == pionPlus) chargeBullet = 1; |
---|
129 | if(bulletType == pionMinus) chargeBullet = -1; |
---|
130 | G4int baryonNumber = int(std::floor(targetA)) + baryonBullet; |
---|
131 | G4int chargeNumber = int(std::floor(targetZ)) + chargeBullet; |
---|
132 | G4double mass = aTrack.GetDefinition()->GetPDGMass(); |
---|
133 | G4double amass = theNucleus.AtomicMass(targetA, targetZ); |
---|
134 | G4double eKinSum = bulletE; |
---|
135 | G4LorentzVector labv = G4LorentzVector(0.0, 0.0, std::sqrt(bulletE*(bulletE + 2.*mass)), bulletE + mass + amass); |
---|
136 | G4cout <<"Energy in the beginning = " << labv.e() / MeV << G4endl; |
---|
137 | #endif |
---|
138 | |
---|
139 | for(int i = 0; i < 15; i++) { |
---|
140 | calincl->f[i] = 0.0; // Initialize INCL input data |
---|
141 | } |
---|
142 | |
---|
143 | // Check wheter the input is acceptable. |
---|
144 | if((bulletType != 0) && ((targetA != 1) && (targetZ != 1))) { |
---|
145 | calincl->f[0] = targetA; // Target mass number |
---|
146 | calincl->f[1] = targetZ; // Charge number |
---|
147 | calincl->f[6] = bulletType; // Type |
---|
148 | calincl->f[2] = bulletE; // Energy [MeV] |
---|
149 | calincl->f[5] = 1.0; // Time scaling |
---|
150 | calincl->f[4] = 45.0; // Nuclear potential |
---|
151 | |
---|
152 | ws->nosurf = -2; // Nucleus surface, -2 = Woods-Saxon |
---|
153 | ws->xfoisa = 8; |
---|
154 | ws->npaulstr = 0; |
---|
155 | |
---|
156 | int nTries = 0; |
---|
157 | varntp->ntrack = 0; |
---|
158 | |
---|
159 | mat->nbmat = 1; |
---|
160 | mat->amat[0] = int(calincl->f[0]); |
---|
161 | mat->zmat[0] = int(calincl->f[1]); |
---|
162 | |
---|
163 | incl->initIncl(true); |
---|
164 | |
---|
165 | while((varntp->ntrack <= 0) && (nTries < maxTries)) { // Loop until we produce real cascade |
---|
166 | nTries++; |
---|
167 | if(verboseLevel > 1) { |
---|
168 | G4cout <<"G4InclCascadeInterface: Try number = " << nTries << G4endl; |
---|
169 | } |
---|
170 | incl->processEventIncl(); |
---|
171 | |
---|
172 | if(verboseLevel > 1) { |
---|
173 | G4cout <<"G4InclCascadeInterface: number of tracks = " << varntp->ntrack <<G4endl; |
---|
174 | } |
---|
175 | } |
---|
176 | |
---|
177 | if(verboseLevel > 1) { |
---|
178 | /** |
---|
179 | * Diagnostic output |
---|
180 | */ |
---|
181 | G4cout <<"G4InclCascadeInterface: Bullet type: " << bulletType << G4endl; |
---|
182 | G4cout <<"G4Incl4AblaCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl; |
---|
183 | |
---|
184 | G4cout <<"G4InclCascadeInterface: Target A: " << targetA << G4endl; |
---|
185 | G4cout <<"G4InclCascadeInterface: Target Z: " << targetZ << G4endl; |
---|
186 | |
---|
187 | if(verboseLevel > 3) { |
---|
188 | diagdata <<"G4InclCascadeInterface: Bullet type: " << bulletType << G4endl; |
---|
189 | diagdata <<"G4InclCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl; |
---|
190 | |
---|
191 | diagdata <<"G4InclCascadeInterface: Target A: " << targetA << G4endl; |
---|
192 | diagdata <<"G4InclCascadeInterface: Target Z: " << targetZ << G4endl; |
---|
193 | } |
---|
194 | |
---|
195 | for(particleI = 0; particleI < varntp->ntrack; particleI++) { |
---|
196 | G4cout << n << " " << calincl->f[6] << " " << calincl->f[2] << " "; |
---|
197 | G4cout << varntp->massini << " " << varntp->mzini << " "; |
---|
198 | G4cout << varntp->exini << " " << varntp->mulncasc << " " << varntp->mulnevap << " " << varntp->mulntot << " "; |
---|
199 | G4cout << varntp->bimpact << " " << varntp->jremn << " " << varntp->kfis << " " << varntp->estfis << " "; |
---|
200 | G4cout << varntp->izfis << " " << varntp->iafis << " " << varntp->ntrack << " " << varntp->itypcasc[particleI] << " "; |
---|
201 | G4cout << varntp->avv[particleI] << " " << varntp->zvv[particleI] << " " << varntp->enerj[particleI] << " "; |
---|
202 | G4cout << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl; |
---|
203 | // For diagnostic output |
---|
204 | if(verboseLevel > 3) { |
---|
205 | diagdata << n << " " << calincl->f[6] << " " << calincl->f[2] << " "; |
---|
206 | diagdata << varntp->massini << " " << varntp->mzini << " "; |
---|
207 | diagdata << varntp->exini << " " << varntp->mulncasc << " " << varntp->mulnevap << " " << varntp->mulntot << " "; |
---|
208 | diagdata << varntp->bimpact << " " << varntp->jremn << " " << varntp->kfis << " " << varntp->estfis << " "; |
---|
209 | diagdata << varntp->izfis << " " << varntp->iafis << " " << varntp->ntrack << " "; |
---|
210 | diagdata << varntp->itypcasc[particleI] << " "; |
---|
211 | diagdata << varntp->avv[particleI] << " " << varntp->zvv[particleI] << " " << varntp->enerj[particleI] << " "; |
---|
212 | diagdata << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl; |
---|
213 | } |
---|
214 | } |
---|
215 | } |
---|
216 | |
---|
217 | // Check whether a valid cascade was produced. |
---|
218 | // If not return the original bullet particle with the same momentum. |
---|
219 | if(varntp->ntrack <= 0) { |
---|
220 | if(verboseLevel > 1) { |
---|
221 | G4cout <<"WARNING G4InclCascadeInterface: No cascade. Returning original particle with original momentum." << G4endl; |
---|
222 | G4cout <<"\t Reached maximum trials of 200 to produce inelastic scattering." << G4endl; |
---|
223 | } |
---|
224 | |
---|
225 | theResult.SetStatusChange(stopAndKill); |
---|
226 | |
---|
227 | if(bulletType == proton) { |
---|
228 | aParticleDefinition = G4Proton::ProtonDefinition(); |
---|
229 | } |
---|
230 | if(bulletType == neutron) { |
---|
231 | aParticleDefinition = G4Neutron::NeutronDefinition(); |
---|
232 | } |
---|
233 | if(bulletType == pionPlus) { |
---|
234 | aParticleDefinition = G4PionPlus::PionPlusDefinition(); |
---|
235 | } |
---|
236 | if(bulletType == pionZero) { |
---|
237 | aParticleDefinition = G4PionZero::PionZeroDefinition(); |
---|
238 | } |
---|
239 | if(bulletType == pionMinus) { |
---|
240 | aParticleDefinition = G4PionMinus::PionMinusDefinition(); |
---|
241 | } |
---|
242 | |
---|
243 | cascadeParticle = new G4DynamicParticle(); |
---|
244 | cascadeParticle->SetDefinition(aParticleDefinition); |
---|
245 | cascadeParticle->Set4Momentum(aTrack.Get4Momentum()); |
---|
246 | theResult.AddSecondary(cascadeParticle); |
---|
247 | } |
---|
248 | |
---|
249 | // Convert INCL4 output to Geant4 compatible data structures. |
---|
250 | // Elementary particles are converted to G4DynamicParticle. |
---|
251 | theResult.SetStatusChange(stopAndKill); |
---|
252 | |
---|
253 | #ifdef DEBUGINCL |
---|
254 | G4cout << "E [MeV]" << std::setw(12) << " Ekin [MeV]" << std::setw(12) << " E* [MeV]" << std::setw(12) << "Px [MeV]" << std::setw(12) << " Py [MeV]" << std::setw(12) << "Pz [MeV]" << std::setw(12) << "Pt [MeV]" << std::setw(12) << "A" << std::setw(12) << "Z" << G4endl; |
---|
255 | #endif |
---|
256 | |
---|
257 | for(particleI = 0; particleI < varntp->ntrack; particleI++) { // Loop through the INCL4+ABLA output. |
---|
258 | // Get energy/momentum and construct momentum vector in INCL4 coordinates. |
---|
259 | momx = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::cos(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV; |
---|
260 | momy = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::sin(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV; |
---|
261 | momz = varntp->plab[particleI]*std::cos(varntp->tetlab[particleI]*CLHEP::pi/180.0)*MeV; |
---|
262 | |
---|
263 | eKin = varntp->enerj[particleI] * MeV; |
---|
264 | |
---|
265 | G4ThreeVector momDirection(momx, momy, momz); // Direction of the particle. |
---|
266 | momDirection = momDirection.unit(); |
---|
267 | if(verboseLevel > 2) { |
---|
268 | G4cout <<"G4InclCascadeInterface: " << G4endl; |
---|
269 | G4cout <<"A = " << varntp->avv[particleI] << " Z = " << varntp->zvv[particleI] << G4endl; |
---|
270 | G4cout <<"eKin = " << eKin << " MeV" << G4endl; |
---|
271 | G4cout <<"px = " << momDirection.x() << " py = " << momDirection.y() <<" pz = " << momDirection.z() << G4endl; |
---|
272 | } |
---|
273 | |
---|
274 | G4int particleIdentified = 0; // Check particle ID. |
---|
275 | |
---|
276 | if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 1)) { // Proton |
---|
277 | cascadeParticle = |
---|
278 | new G4DynamicParticle(G4Proton::ProtonDefinition(), momDirection, eKin); |
---|
279 | particleIdentified++; |
---|
280 | #ifdef DEBUGINCL |
---|
281 | baryonNumber--; |
---|
282 | chargeNumber--; |
---|
283 | #endif |
---|
284 | } |
---|
285 | |
---|
286 | if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 0)) { // Neutron |
---|
287 | cascadeParticle = |
---|
288 | new G4DynamicParticle(G4Neutron::NeutronDefinition(), momDirection, eKin); |
---|
289 | particleIdentified++; |
---|
290 | #ifdef DEBUGINCL |
---|
291 | baryonNumber--; |
---|
292 | #endif |
---|
293 | } |
---|
294 | |
---|
295 | if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 1)) { // PionPlus |
---|
296 | cascadeParticle = |
---|
297 | new G4DynamicParticle(G4PionPlus::PionPlusDefinition(), momDirection, eKin); |
---|
298 | particleIdentified++; |
---|
299 | #ifdef DEBUGINCL |
---|
300 | chargeNumber--; |
---|
301 | #endif |
---|
302 | } |
---|
303 | |
---|
304 | if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 0)) { // PionZero |
---|
305 | cascadeParticle = |
---|
306 | new G4DynamicParticle(G4PionZero::PionZeroDefinition(), momDirection, eKin); |
---|
307 | particleIdentified++; |
---|
308 | } |
---|
309 | |
---|
310 | if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == -1)) { // PionMinus |
---|
311 | cascadeParticle = |
---|
312 | new G4DynamicParticle(G4PionMinus::PionMinusDefinition(), momDirection, eKin); |
---|
313 | particleIdentified++; |
---|
314 | #ifdef DEBUGINCL |
---|
315 | chargeNumber++; |
---|
316 | #endif |
---|
317 | } |
---|
318 | |
---|
319 | if((varntp->avv[particleI] > 1) && (varntp->zvv[particleI] >= 1)) { // Nucleus fragment |
---|
320 | G4ParticleDefinition * aIonDef = 0; |
---|
321 | G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable(); |
---|
322 | |
---|
323 | G4int A = G4int(varntp->avv[particleI]); |
---|
324 | G4int Z = G4int(varntp->zvv[particleI]); |
---|
325 | G4double excitationE = G4double(varntp->exini) * MeV; |
---|
326 | |
---|
327 | if(verboseLevel > 1) { |
---|
328 | G4cout <<"Finding ion: A = " << A << " Z = " << Z << " E* = " << excitationE/MeV << G4endl; |
---|
329 | } |
---|
330 | aIonDef = theTableOfParticles->GetIon(Z, A, excitationE); |
---|
331 | |
---|
332 | if(aIonDef == 0) { |
---|
333 | if(verboseLevel > 1) { |
---|
334 | G4cout <<"G4InclCascadeInterface: " << G4endl; |
---|
335 | G4cout <<"FATAL ERROR: aIonDef = 0" << G4endl; |
---|
336 | G4cout <<"A = " << A << " Z = " << Z << " E* = " << excitationE << G4endl; |
---|
337 | } |
---|
338 | } |
---|
339 | |
---|
340 | if(aIonDef != 0) { // If the ion was identified add it to output. |
---|
341 | cascadeParticle = |
---|
342 | new G4DynamicParticle(aIonDef, momDirection, eKin); |
---|
343 | particleIdentified++; |
---|
344 | #ifdef DEBUGINCL |
---|
345 | baryonNumber = baryonNumber - A; |
---|
346 | chargeNumber = chargeNumber - Z; |
---|
347 | #endif |
---|
348 | } |
---|
349 | } |
---|
350 | |
---|
351 | if(particleIdentified == 1) { // Particle identified properly. |
---|
352 | cascadeParticle->Set4Momentum(cascadeParticle->Get4Momentum()*=toLabFrame); |
---|
353 | #ifdef DEBUGINCL |
---|
354 | G4ParticleDefinition *pd = cascadeParticle->GetDefinition(); |
---|
355 | G4LorentzVector fm = cascadeParticle->Get4Momentum(); |
---|
356 | G4ThreeVector mom = cascadeParticle->GetMomentum(); |
---|
357 | G4double m = pd->GetPDGMass(); |
---|
358 | G4double p = mom.mag(); |
---|
359 | labv -= fm; |
---|
360 | G4double px = mom.x() * MeV; |
---|
361 | G4double py = mom.y() * MeV; |
---|
362 | G4double pz = mom.z() * MeV; |
---|
363 | G4double pt = std::sqrt(px*px+py*py); |
---|
364 | G4double e = fm.e(); |
---|
365 | eKinSum -= cascadeParticle->GetKineticEnergy() * MeV; |
---|
366 | G4double exE; |
---|
367 | if(varntp->avv[particleI] > 1) { |
---|
368 | exE = varntp->exini; |
---|
369 | } |
---|
370 | else { |
---|
371 | exE = 0.0; |
---|
372 | } |
---|
373 | G4cout << fm.e() / MeV |
---|
374 | << std::setw(12) << cascadeParticle->GetKineticEnergy() / MeV |
---|
375 | << std::setw(12) << exE / MeV |
---|
376 | << std::setw(12) << mom.x() / MeV |
---|
377 | << std::setw(12) << mom.y() / MeV |
---|
378 | << std::setw(12) << mom.z() / MeV |
---|
379 | << std::setw(12) << pt / MeV |
---|
380 | << std::setw(12) << varntp->avv[particleI] |
---|
381 | << std::setw(12) << varntp->zvv[particleI] << G4endl; |
---|
382 | #endif |
---|
383 | theResult.AddSecondary(cascadeParticle); // Put data into G4HadFinalState. |
---|
384 | } |
---|
385 | else { // Particle identification failed. |
---|
386 | if(particleIdentified > 1) { // Particle was identified as more than one particle type. |
---|
387 | if(verboseLevel > 1) { |
---|
388 | G4cout <<"G4InclCascadeInterface: One outcoming particle was identified as"; |
---|
389 | G4cout <<"more than one particle type. This is probably due to a bug in the interface." << G4endl; |
---|
390 | G4cout <<"Particle A:" << varntp->avv[particleI] << "Z: " << varntp->zvv[particleI] << G4endl; |
---|
391 | G4cout << "(particleIdentified =" << particleIdentified << ")" << G4endl; |
---|
392 | } |
---|
393 | } |
---|
394 | } |
---|
395 | } |
---|
396 | #ifdef DEBUGINCL |
---|
397 | G4cout <<"--------------------------------------------------------------------------------" << G4endl; |
---|
398 | G4double pt = std::sqrt(std::pow(labv.x(), 2) + std::pow(labv.y(), 2)); |
---|
399 | G4cout << labv.e() / MeV << std::setw(12) << eKinSum / MeV << std::setw(12) << labv.x() << std::setw(12) << labv.y() << std::setw(12) << labv.z() << std::setw(12) << pt / MeV << std::setw(12) << baryonNumber << std::setw(12) << chargeNumber << " totals" << G4endl; |
---|
400 | G4cout << G4endl; |
---|
401 | |
---|
402 | if(verboseLevel > 3) { |
---|
403 | if(baryonNumber != 0) { |
---|
404 | G4cout <<"WARNING G4InclCascadeInterface: Baryon number conservation violated." << G4endl; |
---|
405 | G4cout <<"Baryon number balance after the event: " << baryonNumber << G4endl; |
---|
406 | if(baryonNumber < 0) { |
---|
407 | G4cout <<"Too many baryons produced." << G4endl; |
---|
408 | } else { |
---|
409 | G4cout <<"Too few baryons produced." << G4endl; |
---|
410 | } |
---|
411 | } |
---|
412 | } |
---|
413 | #endif |
---|
414 | |
---|
415 | varntp->ntrack = 0; // Clean up the number of generated particles in the event. |
---|
416 | } |
---|
417 | /** |
---|
418 | * Report unsupported features. |
---|
419 | * (Check bullet, target, energy range) |
---|
420 | */ |
---|
421 | else { // If the bullet type was not recognized by the interface, it will be returned back without any interaction. |
---|
422 | theResult.SetStatusChange(stopAndKill); |
---|
423 | |
---|
424 | G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable(); |
---|
425 | cascadeParticle = new G4DynamicParticle(theTableOfParticles->FindParticle(aTrack.GetDefinition()), aTrack.Get4Momentum()); |
---|
426 | |
---|
427 | theResult.AddSecondary(cascadeParticle); |
---|
428 | |
---|
429 | if(verboseLevel > 1) { |
---|
430 | G4cout <<"ERROR G4InclCascadeInterface: Processing event number (internal) failed " << eventNumber << G4endl; |
---|
431 | } |
---|
432 | if(verboseLevel > 3) { |
---|
433 | diagdata <<"ERROR G4InclCascadeInterface: Processing event number (internal) failed " << eventNumber << G4endl; |
---|
434 | } |
---|
435 | |
---|
436 | if(bulletType == 0) { |
---|
437 | if(verboseLevel > 1) { |
---|
438 | G4cout <<"G4InclCascadeInterface: Unknown bullet type" << G4endl; |
---|
439 | G4cout <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl; |
---|
440 | } |
---|
441 | if(verboseLevel > 3) { |
---|
442 | diagdata <<"G4InclCascadeInterface: Unknown bullet type" << G4endl; |
---|
443 | diagdata <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl; |
---|
444 | } |
---|
445 | } |
---|
446 | |
---|
447 | if((targetA == 1) && (targetZ == 1)) { // Unsupported target |
---|
448 | if(verboseLevel > 1) { |
---|
449 | G4cout <<"Unsupported target: " << G4endl; |
---|
450 | G4cout <<"Target A: " << targetA << G4endl; |
---|
451 | G4cout <<"TargetZ: " << targetZ << G4endl; |
---|
452 | } |
---|
453 | if(verboseLevel > 3) { |
---|
454 | diagdata <<"Unsupported target: " << G4endl; |
---|
455 | diagdata <<"Target A: " << targetA << G4endl; |
---|
456 | diagdata <<"TargetZ: " << targetZ << G4endl; |
---|
457 | } |
---|
458 | } |
---|
459 | |
---|
460 | if(bulletE < 100) { // INCL does not support E < 100 MeV. |
---|
461 | if(verboseLevel > 1) { |
---|
462 | G4cout <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl; |
---|
463 | G4cout <<"WARNING: Returning the original bullet with original energy back to Geant4." << G4endl; |
---|
464 | } |
---|
465 | if(verboseLevel > 3) { |
---|
466 | diagdata <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl; |
---|
467 | } |
---|
468 | } |
---|
469 | |
---|
470 | if(verboseLevel > 3) { |
---|
471 | diagdata <<"WARNING: returning the original bullet with original energy back to Geant4." << G4endl; |
---|
472 | } |
---|
473 | } |
---|
474 | |
---|
475 | return &theResult; |
---|
476 | } |
---|
477 | |
---|
478 | G4ReactionProductVector* G4InclCascadeInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) { |
---|
479 | return 0; |
---|
480 | } |
---|
481 | |
---|
482 | |
---|