source: trunk/source/processes/hadronic/models/incl/src/G4InclLightIonInterface.cc

Last change on this file was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

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