source: trunk/source/processes/hadronic/models/incl/src/G4InclAblaCascadeInterface.cc @ 1150

Last change on this file since 1150 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 20.2 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: G4InclAblaCascadeInterface.cc,v 1.10 2007/12/10 16:31:55 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 "G4InclAblaCascadeInterface.hh"
36#include "math.h"
37#include "G4GenericIon.hh"
38#include "CLHEP/Random/Random.h"
39
40G4InclAblaCascadeInterface::G4InclAblaCascadeInterface()
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
55G4InclAblaCascadeInterface::~G4InclAblaCascadeInterface()
56{
57  delete hazard;
58  delete varntp;
59  delete calincl;
60  delete ws;
61  delete mat;
62  delete incl;
63}
64
65G4HadFinalState* G4InclAblaCascadeInterface::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 << " >>> G4InclAblaCascadeInterface::ApplyYourself called" << G4endl;
81  }
82
83  if(verboseLevel > 1) {
84    G4cout <<"G4InclAblaCascadeInterface: 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  G4LorentzVector labvA = G4LorentzVector(0.0, 0.0, 0.0, 0.0);
137  G4cout <<"Energy in the beginning = " << labv.e() / MeV << G4endl;
138#endif
139
140 for(int i = 0; i < 15; i++) {
141   calincl->f[i] = 0.0; // Initialize INCL input data
142 }
143
144  // Check wheter the input is acceptable.
145  if((bulletType != 0) && ((targetA != 1) && (targetZ != 1))) {
146    calincl->f[0] = targetA;    // Target mass number
147    calincl->f[1] = targetZ;    // Charge number
148    calincl->f[6] = bulletType; // Type
149    calincl->f[2] = bulletE;    // Energy [MeV]
150    calincl->f[5] = 1.0;        // Time scaling
151    calincl->f[4] = 45.0;       // Nuclear potential
152
153    ws->nosurf = -2;  // Nucleus surface, -2 = Woods-Saxon
154    ws->xfoisa = 8;
155    ws->npaulstr = 0;
156
157    int nTries = 0;
158    varntp->ntrack = 0;
159
160    mat->nbmat = 1;
161    mat->amat[0] = int(calincl->f[0]);
162    mat->zmat[0] = int(calincl->f[1]);
163
164    incl->initIncl(true);
165
166    while((varntp->ntrack <= 0) && (nTries < maxTries)) { // Loop until we produce real cascade
167      nTries++;
168      if(verboseLevel > 1) {
169        G4cout <<"G4InclAblaCascadeInterface: Try number = " << nTries << G4endl; 
170      }
171      incl->processEventInclAbla(eventNumber);
172
173      if(verboseLevel > 1) {
174        G4cout <<"G4InclAblaCascadeInterface: number of tracks = " <<  varntp->ntrack <<G4endl;
175      }
176    }
177
178    if(verboseLevel > 1) {
179      /**
180       * Diagnostic output
181       */
182      G4cout <<"G4InclAblaCascadeInterface: Bullet type: " << bulletType << G4endl;
183      G4cout <<"G4Incl4AblaCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
184
185      G4cout <<"G4InclAblaCascadeInterface: Target A:  " << targetA << G4endl;
186      G4cout <<"G4InclAblaCascadeInterface: Target Z:  " << targetZ << G4endl;
187
188      if(verboseLevel > 3) {
189        diagdata <<"G4InclAblaCascadeInterface: Bullet type: " << bulletType << G4endl;
190        diagdata <<"G4InclAblaCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
191       
192        diagdata <<"G4InclAblaCascadeInterface: Target A:  " << targetA << G4endl;
193        diagdata <<"G4InclAblaCascadeInterface: Target Z:  " << targetZ << G4endl;
194      }
195
196      for(particleI = 0; particleI < varntp->ntrack; particleI++) {
197        G4cout << n                       << " " << calincl->f[6]             << " " << calincl->f[2] << " ";
198        G4cout << varntp->massini         << " " << varntp->mzini             << " ";
199        G4cout << varntp->exini           << " " << varntp->mulncasc          << " " << varntp->mulnevap          << " " << varntp->mulntot << " ";
200        G4cout << varntp->bimpact         << " " << varntp->jremn             << " " << varntp->kfis              << " " << varntp->estfis << " ";
201        G4cout << varntp->izfis           << " " << varntp->iafis             << " " << varntp->ntrack            << " " << varntp->itypcasc[particleI] << " ";
202        G4cout << varntp->avv[particleI]  << " " << varntp->zvv[particleI]    << " " << varntp->enerj[particleI]  << " ";
203        G4cout << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
204        // For diagnostic output
205        if(verboseLevel > 3) {
206          diagdata << n                       << " " << calincl->f[6]             << " " << calincl->f[2] << " ";
207          diagdata << varntp->massini         << " " << varntp->mzini             << " ";
208          diagdata << varntp->exini           << " " << varntp->mulncasc          << " " << varntp->mulnevap          << " " << varntp->mulntot << " ";
209          diagdata << varntp->bimpact         << " " << varntp->jremn             << " " << varntp->kfis              << " " << varntp->estfis << " ";
210          diagdata << varntp->izfis           << " " << varntp->iafis             << " " << varntp->ntrack            << " ";
211          diagdata                                                                       << varntp->itypcasc[particleI] << " ";
212          diagdata << varntp->avv[particleI]  << " " << varntp->zvv[particleI]    << " " << varntp->enerj[particleI]  << " ";
213          diagdata << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
214        }
215      }
216    }
217
218    // Check whether a valid cascade was produced.
219    // If not return the original bullet particle with the same momentum.
220    if(varntp->ntrack <= 0) {
221      if(verboseLevel > 1) {
222        G4cout <<"WARNING G4InclAblaCascadeInterface: No cascade. Returning original particle with original momentum." << G4endl;
223        G4cout <<"\t Reached maximum trials of 200 to produce inelastic scattering." << G4endl;
224      }
225
226      theResult.SetStatusChange(stopAndKill);
227     
228      if(bulletType == proton) {
229        aParticleDefinition = G4Proton::ProtonDefinition();
230      }
231      if(bulletType == neutron) {
232        aParticleDefinition = G4Neutron::NeutronDefinition();
233      }
234      if(bulletType == pionPlus) {
235        aParticleDefinition = G4PionPlus::PionPlusDefinition();
236      }
237      if(bulletType == pionZero) {
238        aParticleDefinition = G4PionZero::PionZeroDefinition();
239      }
240      if(bulletType == pionMinus) {
241        aParticleDefinition = G4PionMinus::PionMinusDefinition();
242      }
243
244      cascadeParticle = new G4DynamicParticle();
245      cascadeParticle->SetDefinition(aParticleDefinition);
246      cascadeParticle->Set4Momentum(aTrack.Get4Momentum());
247      theResult.AddSecondary(cascadeParticle);
248    }
249
250    // Convert INCL4 output to Geant4 compatible data structures.
251    // Elementary particles are converted to G4DynamicParticle.
252    theResult.SetStatusChange(stopAndKill);
253
254#ifdef DEBUGINCL
255    G4cout << "E [MeV]" << std::setw(12)
256           << " Ekin [MeV]" << std::setw(12)
257           << "Px [MeV]" << std::setw(12)
258           << " Py [MeV]" << std::setw(12)
259           << "Pz [MeV]" << std::setw(12)
260           << "Pt [MeV]" << std::setw(12)
261           << "A" << std::setw(12)
262           << "Z" << G4endl;
263#endif
264
265    for(particleI = 0; particleI < varntp->ntrack; particleI++) { // Loop through the INCL4+ABLA output.
266      // Get energy/momentum and construct momentum vector in INCL4 coordinates.
267      momx = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::cos(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
268      momy = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::sin(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
269      momz = varntp->plab[particleI]*std::cos(varntp->tetlab[particleI]*CLHEP::pi/180.0)*MeV;
270
271      eKin = varntp->enerj[particleI] * MeV;
272
273      G4ThreeVector momDirection(momx, momy, momz); // Direction of the particle.
274      momDirection = momDirection.unit();
275      if(verboseLevel > 2) {
276        G4cout <<"G4InclAblaCascadeInterface: " << G4endl;
277        G4cout <<"A    = " << varntp->avv[particleI] << " Z = "  << varntp->zvv[particleI] << G4endl;
278        G4cout <<"eKin = " << eKin                   << " MeV"   << G4endl;
279        G4cout <<"px   = " << momDirection.x()       << " py = " << momDirection.y()       <<" pz = " << momDirection.z() << G4endl;
280      }
281
282      G4int particleIdentified = 0; // Check particle ID.
283
284      if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 1)) { // Proton
285        cascadeParticle = 
286          new G4DynamicParticle(G4Proton::ProtonDefinition(), momDirection, eKin);
287        particleIdentified++;
288#ifdef DEBUGINCL
289        baryonNumber--;
290        chargeNumber--;
291#endif
292      }
293
294      if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 0)) { // Neutron
295        cascadeParticle = 
296          new G4DynamicParticle(G4Neutron::NeutronDefinition(), momDirection, eKin);
297        particleIdentified++;
298#ifdef DEBUGINCL
299        baryonNumber--;
300#endif
301      }
302
303      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 1)) { // PionPlus
304        cascadeParticle = 
305          new G4DynamicParticle(G4PionPlus::PionPlusDefinition(), momDirection, eKin);
306        particleIdentified++;
307#ifdef DEBUGINCL
308        chargeNumber--;
309#endif
310      }
311
312      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 0)) { // PionZero
313        cascadeParticle = 
314          new G4DynamicParticle(G4PionZero::PionZeroDefinition(), momDirection, eKin);
315        particleIdentified++;
316      }
317
318      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == -1)) { // PionMinus
319        cascadeParticle = 
320          new G4DynamicParticle(G4PionMinus::PionMinusDefinition(), momDirection, eKin);
321        particleIdentified++;
322#ifdef DEBUGINCL
323        chargeNumber++;
324#endif
325      }
326
327      if((varntp->avv[particleI] > 1) && (varntp->zvv[particleI] >= 1)) { // Nucleus fragment
328        G4ParticleDefinition * aIonDef = 0;
329        G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
330
331        G4int A = G4int(varntp->avv[particleI]);
332        G4int Z = G4int(varntp->zvv[particleI]);
333        G4double excitationE = G4double(varntp->exini) * MeV;
334
335        if(verboseLevel > 1) {
336          G4cout <<"Finding ion: A = " << A << " Z = " << Z << " E* = " << excitationE/MeV << G4endl;
337        }
338        aIonDef = theTableOfParticles->GetIon(Z, A, excitationE);
339       
340        if(aIonDef == 0) {
341          if(verboseLevel > 1) {
342            G4cout <<"G4InclAblaCascadeInterface: " << G4endl;
343            G4cout <<"FATAL ERROR: aIonDef = 0" << G4endl;
344            G4cout <<"A = " << A << " Z = " << Z << " E* = " << excitationE << G4endl;
345          }
346        }
347
348        if(aIonDef != 0) { // If the ion was identified add it to output.
349          cascadeParticle =
350            new G4DynamicParticle(aIonDef, momDirection, eKin);
351          particleIdentified++;
352#ifdef DEBUGINCL
353          baryonNumber = baryonNumber - A;
354          chargeNumber = chargeNumber - Z;
355#endif
356        }
357      }
358       
359      if(particleIdentified == 1) { // Particle identified properly.
360        cascadeParticle->Set4Momentum(cascadeParticle->Get4Momentum()*=toLabFrame);
361#ifdef DEBUGINCL
362        G4ParticleDefinition *pd  = cascadeParticle->GetDefinition();
363        G4LorentzVector fm  = cascadeParticle->Get4Momentum();
364        G4ThreeVector mom = cascadeParticle->GetMomentum();
365        G4double m = pd->GetPDGMass();
366        G4double p = mom.mag();
367        labv -= fm;
368        if(varntp->avv[particleI] > 1) {
369          labvA += fm;
370        }
371        G4double px = mom.x() * MeV;
372        G4double py = mom.y() * MeV;
373        G4double pz = mom.z() * MeV;
374        G4double pt = std::sqrt(px*px+py*py);
375        G4double e  = fm.e();
376        eKinSum -= cascadeParticle->GetKineticEnergy() * MeV;
377        G4double exE;
378        if(varntp->avv[particleI] > 1) {
379          exE = varntp->exini;
380        }
381        else {
382          exE = 0.0;
383        }
384        G4cout << fm.e() / MeV
385               << std::setw(12) << cascadeParticle->GetKineticEnergy() / MeV
386               << std::setw(12) << mom.x() / MeV
387               << std::setw(12) << mom.y() / MeV
388               << std::setw(12) << mom.z() / MeV
389               << std::setw(12) << pt / MeV
390               << std::setw(12) << varntp->avv[particleI]
391               << std::setw(12) << varntp->zvv[particleI] << G4endl;
392#endif
393        theResult.AddSecondary(cascadeParticle); // Put data into G4HadFinalState.
394      }
395      else { // Particle identification failed.
396        if(particleIdentified > 1) { // Particle was identified as more than one particle type.
397          if(verboseLevel > 1) {
398            G4cout <<"G4InclAblaCascadeInterface: One outcoming particle was identified as";
399            G4cout <<"more than one particle type. This is probably due to a bug in the interface." << G4endl;
400            G4cout <<"Particle A:" << varntp->avv[particleI] << "Z: " << varntp->zvv[particleI] << G4endl;
401            G4cout << "(particleIdentified =" << particleIdentified << ")"  << G4endl;
402          }
403        }
404      }
405    }
406
407#ifdef DEBUGINCL
408    G4cout <<"--------------------------------------------------------------------------------" << G4endl;
409    G4double pt = std::sqrt(std::pow(labv.x(), 2) + std::pow(labv.y(), 2));
410    G4double ptA = std::sqrt(std::pow(labvA.x(), 2) + std::pow(labvA.y(), 2));
411    G4cout << labv.e() / MeV << std::setw(12)
412           << eKinSum  / MeV << std::setw(12)
413           << labv.x() / MeV << std::setw(12)
414           << labv.y() / MeV << std::setw(12)
415           << labv.z() / MeV << std::setw(12)
416           << pt       / MeV << std::setw(12)
417           << baryonNumber << std::setw(12)
418           << chargeNumber << " totals" << G4endl;
419    G4cout << " - " << std::setw(12)
420           << " - " << std::setw(12)
421           << labvA.x() / MeV << std::setw(12)
422           << labvA.y() / MeV << std::setw(12)
423           << labvA.z() / MeV << std::setw(12)
424           << ptA       / MeV << std::setw(12)
425           << " - " << std::setw(12) << " - " << " totals ABLA" << G4endl;
426    G4cout << G4endl;
427 
428    if(verboseLevel > 3) {
429      if(baryonNumber != 0) {
430        G4cout <<"WARNING G4InclCascadeInterface: Baryon number conservation violated." << G4endl;
431        G4cout <<"Baryon number balance after the event: " << baryonNumber << G4endl;
432        if(baryonNumber < 0) {
433          G4cout <<"Too many baryons produced." << G4endl;
434        } else {
435          G4cout <<"Too few baryons produced." << G4endl;
436        }
437      }
438    }
439#endif
440   
441    varntp->ntrack = 0; // Clean up the number of generated particles in the event.
442  }
443  /**
444   * Report unsupported features.
445   * (Check bullet, target, energy range)
446   */
447  else { // If the bullet type was not recognized by the interface, it will be returned back without any interaction.
448    theResult.SetStatusChange(stopAndKill);
449
450    G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
451    cascadeParticle = new G4DynamicParticle(theTableOfParticles->FindParticle(aTrack.GetDefinition()), aTrack.Get4Momentum());
452
453    theResult.AddSecondary(cascadeParticle);
454
455    if(verboseLevel > 1) {
456      G4cout <<"ERROR G4InclAblaCascadeInterface: Processing event number (internal) failed " << eventNumber << G4endl;
457    }
458    if(verboseLevel > 3) {
459      diagdata <<"ERROR G4InclAblaCascadeInterface: Error processing event number (internal) failed " << eventNumber << G4endl;
460    }
461
462    if(bulletType == 0) {
463      if(verboseLevel > 1) {
464        G4cout <<"G4InclAblaCascadeInterface: Unknown bullet type" << G4endl;     
465        G4cout <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
466      }
467      if(verboseLevel > 3) {
468        diagdata <<"G4InclAblaCascadeInterface: Unknown bullet type" << G4endl;     
469        diagdata <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
470      }
471    }
472
473    if((targetA == 1) && (targetZ == 1)) { // Unsupported target
474      if(verboseLevel > 1) {
475        G4cout <<"Unsupported target: " << G4endl;
476        G4cout <<"Target A: " << targetA << G4endl;
477        G4cout <<"TargetZ: " << targetZ << G4endl;
478      }
479      if(verboseLevel > 3) {
480        diagdata <<"Unsupported target: " << G4endl;
481        diagdata <<"Target A: " << targetA << G4endl;
482       diagdata <<"TargetZ: " << targetZ << G4endl;
483      }
484    }
485
486    if(bulletE < 100) { // INCL does not support E < 100 MeV.
487      if(verboseLevel > 1) {
488        G4cout <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
489        G4cout <<"WARNING: Returning the original bullet with original energy back to Geant4." << G4endl;
490      }
491      if(verboseLevel > 3) {
492        diagdata <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
493      }
494    }
495
496    if(verboseLevel > 3) {
497      diagdata <<"WARNING: returning the original bullet with original energy back to Geant4." << G4endl;
498    }
499  }
500
501  return &theResult;
502} 
503
504G4ReactionProductVector* G4InclAblaCascadeInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) {
505  return 0;
506}
507
508
Note: See TracBrowser for help on using the repository browser.