source: trunk/source/processes/hadronic/models/incl/src/G4InclAblaLightIonInterface.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 30.4 KB
Line 
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: G4InclAblaLightIonInterface.cc,v 1.14 2010/10/29 06:48:43 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#include <vector>
34
35#include "G4InclAblaLightIonInterface.hh"
36#include "G4FermiBreakUp.hh"
37#include "math.h"
38#include "G4GenericIon.hh"
39#include "CLHEP/Random/Random.h"
40
41G4InclAblaLightIonInterface::G4InclAblaLightIonInterface()
42{
43 hazard = new G4Hazard();
44
45 const G4long* table_entry = CLHEP::HepRandom::getTheSeeds(); // Get random seed from CLHEP.
46 hazard->ial = (*table_entry);
47
48 varntp = new G4VarNtp();
49 calincl = 0;
50 ws = new G4Ws();
51 mat = new G4Mat();
52 incl = new G4Incl(hazard, calincl, ws, mat, varntp);
53 useProjectileSpectator = true;
54 useFermiBreakup = true;
55 incl->setUseProjectileSpectators(useProjectileSpectator);
56 if(!getenv("G4INCLABLANOFERMIBREAKUP")) { // Use Fermi Break-up by default if it is NOT explicitly disabled
57 incl->setUseFermiBreakUp(true);
58 useFermiBreakup = true;
59 }
60 verboseLevel = 0;
61 if(getenv("G4INCLVERBOSE")) {
62 verboseLevel = 1;
63 }
64}
65
66G4InclAblaLightIonInterface::~G4InclAblaLightIonInterface()
67{
68 delete hazard;
69 delete varntp;
70 delete calincl;
71 delete ws;
72 delete mat;
73 delete incl;
74}
75
76G4HadFinalState* G4InclAblaLightIonInterface::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
77{
78 // const G4bool useFermiBreakup = false;
79 G4int maxTries = 200;
80
81 G4int particleI;
82
83 G4int baryonNumberBalanceInINCL = 0;
84 G4int chargeNumberBalanceInINCL = 0;
85
86 G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
87
88 // Increase the event number:
89 eventNumber++;
90
91 // Clean up the INCL input
92 if(calincl != 0) {
93 delete calincl;
94 calincl = 0;
95 }
96
97 if (verboseLevel > 1) {
98 G4cout << " >>> G4InclAblaLightIonInterface::ApplyYourself called" << G4endl;
99 }
100
101 if(verboseLevel > 1) {
102 G4cout <<"G4InclAblaLightIonInterface: Now processing INCL4 event number:" << eventNumber << G4endl;
103 }
104
105 if(verboseLevel > 1) {
106 G4Calincl::printProjectileTargetInfo(aTrack, theNucleus);
107 }
108
109 // Inverse kinematics for targets with Z = 1 and A = 1
110 // if(false) {
111 G4LorentzRotation toBreit = aTrack.Get4Momentum().boostVector();
112
113 if(theNucleus.GetZ_asInt() == 1 && theNucleus.GetA_asInt() == 1 && G4Calincl::canUseInverseKinematics(aTrack, theNucleus)) {
114 G4ParticleDefinition *oldTargetDef = theTableOfParticles->GetIon(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt(), 0.0);
115 const G4ParticleDefinition *oldProjectileDef = aTrack.GetDefinition();
116
117 G4int oldTargetA = oldTargetDef->GetAtomicMass();
118 G4int newTargetA = oldProjectileDef->GetAtomicMass();
119 G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
120
121 if(newTargetA > 0 && newTargetZ > 0 && oldTargetDef != 0 && oldProjectileDef != 0) {
122 G4Nucleus swappedTarget(oldProjectileDef->GetAtomicMass(), oldProjectileDef->GetAtomicNumber());
123
124 // G4cout <<"Original projectile kinE = " << aTrack.GetKineticEnergy() / MeV << G4endl;
125
126 // We need the same energy/nucleon.
127 G4double projectileE = ((aTrack.GetKineticEnergy() / MeV) / newTargetA) * oldTargetA * MeV;
128
129 // G4cout <<"projectileE = " << projectileE << G4endl;
130 G4DynamicParticle swappedProjectileParticle(oldTargetDef, G4ThreeVector(0.0, 0.0, 1.0), projectileE);
131 const G4LorentzVector swapped4Momentum = (swappedProjectileParticle.Get4Momentum()*=toBreit);
132 swappedProjectileParticle.Set4Momentum(swapped4Momentum);
133 const G4HadProjectile swappedProjectile(swappedProjectileParticle);
134 // G4cout <<"New projectile kinE = " << swappedProjectile.GetKineticEnergy() / MeV << G4endl;
135 calincl = new G4InclInput(swappedProjectile, swappedTarget, true);
136 } else {
137 G4cout <<"Badly defined target after swapping. Falling back to normal (non-swapped) mode." << G4endl;
138 calincl = new G4InclInput(aTrack, theNucleus, false);
139 }
140 } else {
141 calincl = new G4InclInput(aTrack, theNucleus, false);
142 }
143
144 G4double eKin;
145 G4double momx = 0.0, momy = 0.0, momz = 0.0;
146 G4DynamicParticle *cascadeParticle = 0;
147 G4ParticleDefinition *aParticleDefinition = 0;
148
149 // INCL assumes the projectile particle is going in the direction of
150 // the Z-axis. Here we construct proper rotation to convert the
151 // momentum vectors of the outcoming particles to the original
152 // coordinate system.
153 G4LorentzVector projectileMomentum = aTrack.Get4Momentum();
154 G4LorentzRotation toZ;
155 toZ.rotateZ(-projectileMomentum.phi());
156 toZ.rotateY(-projectileMomentum.theta());
157 G4LorentzRotation toLabFrame = toZ.inverse();
158
159 /*
160 G4cout <<"Projectile theta = " << projectileMomentum.theta() << " phi = " << projectileMomentum.phi() << G4endl;
161 G4cout <<"Projectile momentum "
162 << "(px = " << projectileMomentum.px()
163 << ", py = " << projectileMomentum.py()
164 << ", pz = " << projectileMomentum.pz() << ")" << G4endl;
165 G4cout << "Projectile energy = " << bulletE << " MeV" << G4endl;
166 */
167
168 G4FermiBreakUp *fermiBreakUp = new G4FermiBreakUp();
169 G4FragmentVector *theSpectatorFermiBreakupResult = 0;
170 G4FragmentVector *theFermiBreakupResult = 0;
171
172 theResult.Clear(); // Make sure the output data structure is clean.
173
174 std::vector<G4DynamicParticle*> result; // Temporary list for the results
175
176 // Map Geant4 particle types to corresponding INCL4 types.
177 enum bulletParticleType {nucleus = 0, proton = 1, neutron = 2, pionPlus = 3, pionZero = 4,
178 pionMinus = 5, deuteron = 6, triton = 7, he3 = 8, he4 = 9,
179 c12 = -12}; // Carbon beam support.
180
181 G4int bulletType = calincl->bulletType();
182 chargeNumberBalanceInINCL = calincl->targetZ();
183 baryonNumberBalanceInINCL = calincl->targetA();
184
185 // G4cout <<"Type of the projectile (INCL projectile code): " << bulletType << G4endl;
186
187 if(bulletType == proton) {
188 chargeNumberBalanceInINCL += 1;
189 baryonNumberBalanceInINCL += 1;
190 } else if(bulletType == neutron) {
191 baryonNumberBalanceInINCL += 1;
192 } else if(bulletType == pionPlus) { //Note: positive pion doesn't contribute to the baryon and charge number counters
193 chargeNumberBalanceInINCL += 1;
194 } else if(bulletType == pionMinus) {
195 chargeNumberBalanceInINCL -= 1;
196 } else if(bulletType == deuteron) {
197 chargeNumberBalanceInINCL += 1;
198 baryonNumberBalanceInINCL += 2;
199 } else if(bulletType == triton) {
200 chargeNumberBalanceInINCL += 1;
201 baryonNumberBalanceInINCL += 3;
202 } else if(bulletType == he3) {
203 chargeNumberBalanceInINCL += 2;
204 baryonNumberBalanceInINCL += 3;
205 } else if(bulletType == he4) {
206 chargeNumberBalanceInINCL += 2;
207 baryonNumberBalanceInINCL += 4;
208 } if(bulletType == c12) {
209 chargeNumberBalanceInINCL += 6;
210 baryonNumberBalanceInINCL += 12;
211 } if(bulletType == -666) {
212 chargeNumberBalanceInINCL += calincl->extendedProjectileZ();
213 baryonNumberBalanceInINCL += calincl->extendedProjectileA();
214 }
215
216 // Check wheter the input is acceptable.
217 if((bulletType != 0) && ((calincl->targetA() != 1) && (calincl->targetZ() != 1))) {
218 ws->nosurf = -2; // Nucleus surface, -2 = Woods-Saxon
219 ws->xfoisa = 8;
220 ws->npaulstr = 0;
221
222 int nTries = 0;
223 varntp->ntrack = 0;
224
225 mat->nbmat = 1;
226 mat->amat[0] = int(calincl->targetA());
227 mat->zmat[0] = int(calincl->targetA());
228
229 incl->setInput(calincl);
230 incl->initIncl(true);
231
232 while((varntp->ntrack <= 0) && (nTries < maxTries)) { // Loop until we produce real cascade
233 nTries++;
234 if(verboseLevel > 1) {
235 G4cout <<"G4InclAblaLightIonInterface: Try number = " << nTries << G4endl;
236 }
237 incl->processEventInclAbla(calincl, eventNumber);
238
239 if(verboseLevel > 1) {
240 G4cout <<"G4InclAblaLightIonInterface: number of tracks = " << varntp->ntrack <<G4endl;
241 }
242 }
243
244 if(verboseLevel > 1) {
245 /**
246 * Diagnostic output
247 */
248 G4cout <<"G4InclAblaLightIonInterface: Bullet type: " << calincl->bulletType() << G4endl;
249 G4cout <<"G4Incl4AblaCascadeInterface: Bullet energy: " << calincl->bulletE() << " MeV" << G4endl;
250 if(bulletType == -666) {
251 G4cout <<" Extended projectile: A = " << calincl->extendedProjectileA()
252 <<" Z = " << calincl->extendedProjectileZ() << G4endl;
253 }
254
255 G4cout <<"G4InclAblaLightIonInterface: Target A: " << calincl->targetA() << G4endl;
256 G4cout <<"G4InclAblaLightIonInterface: Target Z: " << calincl->targetZ() << G4endl;
257
258 if(verboseLevel > 3) {
259 diagdata <<"G4InclAblaLightIonInterface: Bullet type: " << calincl->bulletType() << G4endl;
260 diagdata <<"G4InclAblaLightIonInterface: Bullet energy: " << calincl->bulletE() << " MeV" << G4endl;
261
262 diagdata <<"G4InclAblaLightIonInterface: Target A: " << calincl->targetA() << G4endl;
263 diagdata <<"G4InclAblaLightIonInterface: Target Z: " << calincl->targetZ() << G4endl;
264 }
265 }
266
267 // Check whether a valid cascade was produced.
268 // If not return the original bullet particle with the same momentum.
269 if(varntp->ntrack <= 0) {
270 if(verboseLevel > 1) {
271 G4cout <<"WARNING G4InclAblaLightIonInterface: No cascade. Returning original particle with original momentum." << G4endl;
272 G4cout <<"\t Reached maximum trials of 200 to produce inelastic scattering." << G4endl;
273 }
274
275 theResult.SetStatusChange(stopAndKill);
276
277 if(bulletType == proton) {
278 aParticleDefinition = G4Proton::ProtonDefinition();
279 } else if(bulletType == neutron) {
280 aParticleDefinition = G4Neutron::NeutronDefinition();
281 } else if(bulletType == pionPlus) {
282 aParticleDefinition = G4PionPlus::PionPlusDefinition();
283 } else if(bulletType == pionZero) {
284 aParticleDefinition = G4PionZero::PionZeroDefinition();
285 } else if(bulletType == pionMinus) {
286 aParticleDefinition = G4PionMinus::PionMinusDefinition();
287 } else if(bulletType == deuteron) {
288 aParticleDefinition = G4Deuteron::DeuteronDefinition();
289 } else if(bulletType == triton) {
290 aParticleDefinition = G4Triton::TritonDefinition();
291 } else if(bulletType == he3) {
292 aParticleDefinition = G4He3::He3Definition();
293 } else if(bulletType == he4) {
294 aParticleDefinition = G4Alpha::AlphaDefinition();
295 } else { // Particle was not recognized. Probably an unsupported particle was given as input
296 aParticleDefinition = 0;
297 }
298
299 if(aParticleDefinition != 0) {
300 cascadeParticle = new G4DynamicParticle();
301 cascadeParticle->SetDefinition(aParticleDefinition);
302 cascadeParticle->Set4Momentum(aTrack.Get4Momentum());
303 result.push_back(cascadeParticle);
304 }
305 }
306
307 // Convert INCL4 output to Geant4 compatible data structures.
308 // Elementary particles are converted to G4DynamicParticle.
309 theResult.SetStatusChange(stopAndKill);
310
311 for(particleI = 0; particleI <= varntp->ntrack; particleI++) { // Loop through the INCL4+ABLA output.
312 // Get energy/momentum and construct momentum vector in INCL4 coordinates.
313 // if(varntp->itypcasc[particleI] == -1) continue; // Avoid nucleons that are part of the spectator
314 if(varntp->avv[particleI] == 0 && varntp->zvv[particleI] == 0) continue;
315 momx = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::cos(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
316 momy = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::sin(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
317 momz = varntp->plab[particleI]*std::cos(varntp->tetlab[particleI]*CLHEP::pi/180.0)*MeV;
318
319 eKin = varntp->enerj[particleI] * MeV;
320
321 G4ThreeVector momDirection(momx, momy, momz); // Direction of the particle.
322 momDirection = momDirection.unit();
323 if(verboseLevel > 2) {
324 G4cout <<"G4InclAblaLightIonInterface: " << G4endl;
325 G4cout <<"A = " << varntp->avv[particleI] << " Z = " << varntp->zvv[particleI] << G4endl;
326 G4cout <<"eKin = " << eKin << " MeV" << G4endl;
327 G4cout <<"px = " << momDirection.x() << " py = " << momDirection.y() <<" pz = " << momDirection.z() << G4endl;
328 }
329
330 G4int particleIdentified = 0; // Check particle ID.
331
332 if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 1)) { // Proton
333 cascadeParticle =
334 new G4DynamicParticle(G4Proton::ProtonDefinition(), momDirection, eKin);
335 particleIdentified++;
336 baryonNumberBalanceInINCL -= 1;
337 chargeNumberBalanceInINCL -= 1;
338 }
339
340 if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 0)) { // Neutron
341 cascadeParticle =
342 new G4DynamicParticle(G4Neutron::NeutronDefinition(), momDirection, eKin);
343 particleIdentified++;
344 baryonNumberBalanceInINCL -= 1;
345 }
346
347 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 1)) { // PionPlus
348 cascadeParticle =
349 new G4DynamicParticle(G4PionPlus::PionPlusDefinition(), momDirection, eKin);
350 particleIdentified++;
351 chargeNumberBalanceInINCL -= 1;
352 }
353
354 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 0)) { // PionZero
355 cascadeParticle =
356 new G4DynamicParticle(G4PionZero::PionZeroDefinition(), momDirection, eKin);
357 particleIdentified++;
358 chargeNumberBalanceInINCL -= 0;
359 }
360
361 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == -1)) { // PionMinus
362 cascadeParticle =
363 new G4DynamicParticle(G4PionMinus::PionMinusDefinition(), momDirection, eKin);
364 particleIdentified++;
365 chargeNumberBalanceInINCL -= -1;
366 }
367
368 if((varntp->avv[particleI] > 1) && (varntp->zvv[particleI] >= 1)) { // Nucleus fragment
369 G4ParticleDefinition * aIonDef = 0;
370
371 G4int A = G4int(varntp->avv[particleI]);
372 G4int Z = G4int(varntp->zvv[particleI]);
373 G4double excitationE = G4double(varntp->exini) * MeV;
374
375 if(verboseLevel > 1) {
376 G4cout <<"Finding ion: A = " << A << " Z = " << Z << " E* = " << excitationE/MeV << G4endl;
377 }
378 aIonDef = theTableOfParticles->GetIon(Z, A, excitationE);
379
380 if(aIonDef == 0) {
381 if(verboseLevel > 1) {
382 G4cout <<"G4InclAblaLightIonInterface: " << G4endl;
383 G4cout <<"FATAL ERROR: aIonDef = 0" << G4endl;
384 G4cout <<"A = " << A << " Z = " << Z << " E* = " << excitationE << G4endl;
385 }
386 }
387
388 if(aIonDef != 0) { // If the ion was identified add it to output.
389 cascadeParticle =
390 new G4DynamicParticle(aIonDef, momDirection, eKin);
391 particleIdentified++;
392 baryonNumberBalanceInINCL -= A;
393 chargeNumberBalanceInINCL -= Z;
394 }
395 }
396
397 if(particleIdentified == 1) { // Particle identified properly.
398 cascadeParticle->Set4Momentum(cascadeParticle->Get4Momentum()*=toLabFrame);
399 result.push_back(cascadeParticle);
400 }
401 else { // Particle identification failed.
402 if(particleIdentified > 1) { // Particle was identified as more than one particle type.
403 if(verboseLevel > 1) {
404 G4cout <<"G4InclAblaLightIonInterface: One outcoming particle was identified as";
405 G4cout <<"more than one particle type. This is probably due to a bug in the interface." << G4endl;
406 G4cout <<"Particle A:" << varntp->avv[particleI] << "Z: " << varntp->zvv[particleI] << G4endl;
407 G4cout << "(particleIdentified =" << particleIdentified << ")" << G4endl;
408 }
409 }
410 }
411 }
412
413 // Spectator nucleus Fermi break-up
414 if(useFermiBreakup && useProjectileSpectator && varntp->masp > 1) {
415 baryonNumberBalanceInINCL -= G4int(varntp->masp);
416 G4double nuclearMass = G4NucleiProperties::GetNuclearMass(G4int(varntp->masp), G4int(varntp->mzsp)) + varntp->exsp * MeV;
417 // Use momentum scaling to compensate for different masses in G4 and INCL:
418 G4double momentumScaling = G4InclUtils::calculate4MomentumScaling(G4int(varntp->masp),
419 G4int(varntp->mzsp),
420 varntp->exsp,
421 varntp->spectatorT,
422 varntp->spectatorP1,
423 varntp->spectatorP2,
424 varntp->spectatorP3);
425 G4LorentzVector p4(momentumScaling * varntp->spectatorP1 * MeV, momentumScaling * varntp->spectatorP2 * MeV,
426 momentumScaling * varntp->spectatorP3 * MeV,
427 varntp->spectatorT * MeV + nuclearMass);
428 // Four-momentum, baryon number and charge balance:
429 G4LorentzVector fourMomentumBalance = p4;
430 G4int baryonNumberBalance = G4int(varntp->masp);
431 chargeNumberBalanceInINCL -= G4int(varntp->mzsp);
432 G4int chargeBalance = G4int(varntp->mzsp);
433
434 G4LorentzRotation toFragmentZ;
435 // Assume that Fermi breakup uses Z as the direction of the projectile
436 toFragmentZ.rotateZ(-p4.theta());
437 toFragmentZ.rotateY(-p4.phi());
438 G4LorentzRotation toFragmentLab = toFragmentZ.inverse();
439 // p4 *= toFragmentZ;
440
441 G4LorentzVector p4rest = p4;
442 // p4rest.boost(-p4.boostVector());
443 if(verboseLevel > 0) {
444 G4cout <<"Spectator nucleus:" << G4endl;
445 G4cout <<"p4: " << G4endl;
446 G4cout <<" px: " << p4.px() <<" py: " << p4.py() <<" pz: " << p4.pz() << G4endl;
447 G4cout <<" E = " << p4.e() << G4endl;
448 G4cout <<"p4rest: " << G4endl;
449 G4cout <<" px: " << p4rest.px() <<" py: " << p4rest.py() <<" pz: " << p4rest.pz() << G4endl;
450 G4cout <<" E = " << p4rest.e() << G4endl;
451 }
452 G4Fragment theSpectatorNucleus(G4int(varntp->masp), G4int(varntp->mzsp), p4rest);
453 theSpectatorFermiBreakupResult = fermiBreakUp->BreakItUp(theSpectatorNucleus);
454 if(theSpectatorFermiBreakupResult != 0) {
455 G4FragmentVector::iterator fragment;
456 for(fragment = theSpectatorFermiBreakupResult->begin(); fragment != theSpectatorFermiBreakupResult->end(); fragment++) {
457 G4ParticleDefinition *theFragmentDefinition = 0;
458 if((*fragment)->GetA_asInt() == 1 && (*fragment)->GetZ_asInt() == 0) { // Neutron
459 theFragmentDefinition = G4Neutron::NeutronDefinition();
460 } else if ((*fragment)->GetA_asInt() == 1 && (*fragment)->GetZ_asInt() == 1) {
461 theFragmentDefinition = G4Proton::ProtonDefinition();
462 } else {
463 theFragmentDefinition = theTableOfParticles->GetIon((*fragment)->GetZ_asInt(), (*fragment)->GetA_asInt(), (*fragment)->GetExcitationEnergy());
464 }
465 if(theFragmentDefinition != 0) {
466 G4DynamicParticle *theFragment = new G4DynamicParticle(theFragmentDefinition, (*fragment)->GetMomentum());
467 G4LorentzVector labMomentum = theFragment->Get4Momentum();
468 // labMomentum.boost(p4.boostVector());
469 // labMomentum *= toFragmentLab;
470 // labMomentum *= toLabFrame;
471 theFragment->Set4Momentum(labMomentum);
472 fourMomentumBalance -= theFragment->Get4Momentum();
473 baryonNumberBalance -= theFragmentDefinition->GetAtomicMass();
474 chargeBalance -= theFragmentDefinition->GetAtomicNumber();
475 if(verboseLevel > 0) {
476 G4cout <<"Resulting fragment: " << G4endl;
477 G4cout <<" kinetic energy = " << theFragment->GetKineticEnergy() / MeV << " MeV" << G4endl;
478 G4cout <<" momentum = " << theFragment->GetMomentum().mag() / MeV << " MeV" << G4endl;
479 }
480 result.push_back(theFragment);
481 } else {
482 G4cout <<"G4InclAblaCascadeInterface: Error. Fragment produced by Fermi break-up does not exist."
483 << G4endl;
484 G4cout <<"Resulting fragment: " << G4endl;
485 G4cout <<" Z = " << (*fragment)->GetZ_asInt() << G4endl;
486 G4cout <<" A = " << (*fragment)->GetA_asInt() << G4endl;
487 G4cout <<" Excitation : " << (*fragment)->GetExcitationEnergy() / MeV << " MeV" << G4endl;
488 G4cout <<" momentum = " << (*fragment)->GetMomentum().mag() / MeV << " MeV" << G4endl;
489 }
490 }
491 if(std::abs(fourMomentumBalance.mag() / MeV) > 0.1 * MeV) {
492 G4cout <<"Four-momentum balance after spectator nucleus Fermi break-up:" << G4endl;
493 G4cout <<"Magnitude: " << fourMomentumBalance.mag() / MeV << " MeV" << G4endl;
494 G4cout <<"Vector components (px, py, pz, E) = ("
495 << fourMomentumBalance.px() << ", "
496 << fourMomentumBalance.py() << ", "
497 << fourMomentumBalance.pz() << ", "
498 << fourMomentumBalance.e() << ")" << G4endl;
499 }
500 if(baryonNumberBalance != 0) {
501 G4cout <<"Event " << eventNumber << ": Baryon number balance after spectator nucleus Fermi break-up: " << baryonNumberBalance << G4endl;
502 }
503 if(chargeBalance != 0) {
504 G4cout <<"Event " << eventNumber <<": Charge balance after spectator nucleus Fermi break-up: " << chargeBalance << G4endl;
505 }
506 }
507 }
508
509 // Finally do Fermi break-up if needed
510 if(varntp->needsFermiBreakup && varntp->massini > 0) {
511 baryonNumberBalanceInINCL -= G4int(varntp->massini);
512 chargeNumberBalanceInINCL -= G4int(varntp->mzini);
513 // Call Fermi Break-up
514 G4double nuclearMass = G4NucleiProperties::GetNuclearMass(G4int(varntp->massini), G4int(varntp->mzini)) + varntp->exini * MeV;
515 G4LorentzVector fragmentMomentum(varntp->pxrem * MeV, varntp->pyrem * MeV, varntp->pzrem * MeV,
516 varntp->erecrem * MeV + nuclearMass);
517 G4double momentumScaling = G4InclUtils::calculate4MomentumScaling(G4int(varntp->massini), G4int(varntp->mzini),
518 varntp->exini,
519 varntp->erecrem,
520 varntp->pxrem,
521 varntp->pyrem,
522 varntp->pzrem);
523 G4LorentzVector p4(momentumScaling * varntp->pxrem * MeV, momentumScaling * varntp->pyrem * MeV,
524 momentumScaling * varntp->pzrem * MeV,
525 varntp->erecrem + nuclearMass);
526
527 // For four-momentum, baryon number and charge conservation check:
528 G4LorentzVector fourMomentumBalance = p4;
529 G4int baryonNumberBalance = G4int(varntp->massini);
530 G4int chargeBalance = G4int(varntp->mzini);
531
532 G4LorentzRotation toFragmentZ;
533 toFragmentZ.rotateZ(-p4.theta());
534 toFragmentZ.rotateY(-p4.phi());
535 G4LorentzRotation toFragmentLab = toFragmentZ.inverse();
536 // p4 *= toFragmentZ;
537
538 G4LorentzVector p4rest = p4;
539 // p4rest.boost(-p4.boostVector());
540 if(verboseLevel > 0) {
541 G4cout <<"Cascade remnant nucleus:" << G4endl;
542 G4cout <<"p4: " << G4endl;
543 G4cout <<" px: " << p4.px() <<" py: " << p4.py() <<" pz: " << p4.pz() << G4endl;
544 G4cout <<" E = " << p4.e() << G4endl;
545
546 G4cout <<"p4rest: " << G4endl;
547 G4cout <<" px: " << p4rest.px() <<" py: " << p4rest.py() <<" pz: " << p4rest.pz() << G4endl;
548 G4cout <<" E = " << p4rest.e() << G4endl;
549 }
550
551 G4Fragment theCascadeRemnant(G4int(varntp->massini), G4int(varntp->mzini), p4rest);
552 theFermiBreakupResult = fermiBreakUp->BreakItUp(theCascadeRemnant);
553 if(theFermiBreakupResult != 0) {
554 G4FragmentVector::iterator fragment;
555 for(fragment = theFermiBreakupResult->begin(); fragment != theFermiBreakupResult->end(); fragment++) {
556 G4ParticleDefinition *theFragmentDefinition = 0;
557 if((*fragment)->GetA_asInt() == 1 && (*fragment)->GetZ_asInt() == 0) { // Neutron
558 theFragmentDefinition = G4Neutron::NeutronDefinition();
559 } else if ((*fragment)->GetA_asInt() == 1 && (*fragment)->GetZ_asInt() == 1) {
560 theFragmentDefinition = G4Proton::ProtonDefinition();
561 } else {
562 theFragmentDefinition = theTableOfParticles->GetIon((*fragment)->GetZ_asInt(), (*fragment)->GetA_asInt(), (*fragment)->GetExcitationEnergy());
563 }
564
565 if(theFragmentDefinition != 0) {
566 G4DynamicParticle *theFragment = new G4DynamicParticle(theFragmentDefinition, (*fragment)->GetMomentum());
567 G4LorentzVector labMomentum = theFragment->Get4Momentum();
568 // labMomentum.boost(p4.boostVector());
569 // labMomentum *= toFragmentLab;
570 // labMomentum *= toLabFrame;
571 theFragment->Set4Momentum(labMomentum);
572 fourMomentumBalance -= theFragment->Get4Momentum();
573 baryonNumberBalance -= theFragmentDefinition->GetAtomicMass();
574 chargeBalance -= theFragmentDefinition->GetAtomicNumber();
575 if(verboseLevel > 0) {
576 G4cout <<"Resulting fragment: " << G4endl;
577 G4cout <<" kinetic energy = " << theFragment->GetKineticEnergy() / MeV << " MeV" << G4endl;
578 G4cout <<" momentum = " << theFragment->GetMomentum().mag() / MeV << " MeV" << G4endl;
579 }
580 result.push_back(theFragment);
581 } else {
582 G4cout <<"G4InclAblaCascadeInterface: Error. Fragment produced by Fermi break-up does not exist." << G4endl;
583 G4cout <<"Resulting fragment: " << G4endl;
584 G4cout <<" Z = " << (*fragment)->GetZ_asInt() << G4endl;
585 G4cout <<" A = " << (*fragment)->GetA_asInt() << G4endl;
586 G4cout <<" Excitation : " << (*fragment)->GetExcitationEnergy() / MeV << " MeV" << G4endl;
587 G4cout <<" momentum = " << (*fragment)->GetMomentum().mag() / MeV << " MeV" << G4endl;
588 }
589 }
590 if(std::abs(fourMomentumBalance.mag() / MeV) > 0.1 * MeV) {
591 G4cout <<"Four-momentum balance after remnant nucleus Fermi break-up:" << G4endl;
592 G4cout <<"Magnitude: " << fourMomentumBalance.mag() / MeV << " MeV" << G4endl;
593 G4cout <<"Vector components (px, py, pz, E) = ("
594 << fourMomentumBalance.px() << ", "
595 << fourMomentumBalance.py() << ", "
596 << fourMomentumBalance.pz() << ", "
597 << fourMomentumBalance.e() << ")" << G4endl;
598 }
599 if(baryonNumberBalance != 0) {
600 G4cout <<"Baryon number balance after remnant nucleus Fermi break-up: " << baryonNumberBalance << G4endl;
601 }
602 if(chargeBalance != 0) {
603 G4cout <<"Charge balance after remnant nucleus Fermi break-up: " << chargeBalance << G4endl;
604 }
605 }
606 }
607
608 varntp->ntrack = 0; // Clean up the number of generated particles in the event.
609
610 if(baryonNumberBalanceInINCL != 0 && verboseLevel > 1) {
611 G4cout <<"Event " << eventNumber <<": G4InclAblaLightIonInterface: Baryon number conservation problem in INCL detected!" << G4endl;
612 G4cout <<"Baryon number balance: " << baryonNumberBalanceInINCL << G4endl;
613 if(baryonNumberBalanceInINCL < 0) {
614 G4cout <<"Event " << eventNumber <<": Too many outcoming baryons!" << G4endl;
615 } else if(baryonNumberBalanceInINCL > 0) {
616 G4cout <<"Event " << eventNumber <<": Too few outcoming baryons!" << G4endl;
617 }
618 }
619
620 if(chargeNumberBalanceInINCL != 0 && verboseLevel > 1) {
621 G4cout <<"Event " << eventNumber <<": G4InclAblaLightIonInterface: Charge number conservation problem in INCL detected!" << G4endl;
622 G4cout <<"Event " << eventNumber <<": Charge number balance: " << chargeNumberBalanceInINCL << G4endl;
623 }
624 }
625 /**
626 * Report unsupported features.
627 * (Check bullet, target, energy range)
628 */
629 else { // If the bullet type was not recognized by the interface, it will be returned back without any interaction.
630 theResult.SetStatusChange(stopAndKill);
631
632 G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
633 cascadeParticle = new G4DynamicParticle(theTableOfParticles->FindParticle(aTrack.GetDefinition()), aTrack.Get4Momentum());
634
635 result.push_back(cascadeParticle);
636
637 if(verboseLevel > 1) {
638 G4cout <<"G4InclAblaLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
639 }
640 if(verboseLevel > 3) {
641 diagdata <<"G4InclAblaLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
642 }
643
644 if(bulletType == 0) {
645 if(verboseLevel > 1) {
646 G4cout <<"G4InclAblaLightIonInterface: Unknown bullet type" << G4endl;
647 G4cout <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
648 }
649 if(verboseLevel > 3) {
650 diagdata <<"G4InclAblaLightIonInterface: Unknown bullet type" << G4endl;
651 diagdata <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
652 }
653 }
654
655 if((calincl->targetA() == 1) && (calincl->targetZ() == 1)) { // Unsupported target
656 if(verboseLevel > 1) {
657 G4cout <<"Unsupported target: " << G4endl;
658 G4cout <<"Target A: " << calincl->targetA() << G4endl;
659 G4cout <<"TargetZ: " << calincl->targetZ() << G4endl;
660 }
661 if(verboseLevel > 3) {
662 diagdata <<"Unsupported target: " << G4endl;
663 diagdata <<"Target A: " << calincl->targetA() << G4endl;
664 diagdata <<"TargetZ: " << calincl->targetZ() << G4endl;
665 }
666 }
667
668 if(calincl->bulletE() < 100) { // INCL does not support E < 100 MeV.
669 if(verboseLevel > 1) {
670 G4cout <<"Unsupported bullet energy: " << calincl->bulletE() << " MeV. (Lower limit is 100 MeV)." << G4endl;
671 G4cout <<"WARNING: Returning the original bullet with original energy back to Geant4." << G4endl;
672 }
673 if(verboseLevel > 3) {
674 diagdata <<"Unsupported bullet energy: " << calincl->bulletE() << " MeV. (Lower limit is 100 MeV)." << G4endl;
675 }
676 }
677
678 if(verboseLevel > 3) {
679 diagdata <<"WARNING: returning the original bullet with original energy back to Geant4." << G4endl;
680 }
681 }
682
683 // Finally copy the accumulated secondaries into the result collection:
684 G4ThreeVector boostVector = aTrack.Get4Momentum().boostVector();
685 G4LorentzRotation boostBack = toBreit.inverse();
686
687 for(std::vector<G4DynamicParticle*>::iterator i = result.begin(); i != result.end(); ++i) {
688 // If the calculation was performed in inverse kinematics we have to
689 // convert the result back...
690 if(calincl->isInverseKinematics()) {
691 G4LorentzVector mom = (*i)->Get4Momentum();
692 mom.setPz(-1.0 * mom.pz()); // Reverse the z-component of the momentum vector
693 mom *= boostBack;
694 (*i)->Set4Momentum(mom);
695 }
696 theResult.AddSecondary((*i));
697 }
698
699 delete calincl;
700 calincl = 0;
701 return &theResult;
702}
703
704G4ReactionProductVector* G4InclAblaLightIonInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) {
705 return 0;
706}
707
708
Note: See TracBrowser for help on using the repository browser.