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

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

import all except CVS

File size: 16.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: G4InclLightIonInterface.cc,v 1.10 2007/12/10 16:32:07 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 "G4InclLightIonInterface.hh"
34#include "math.h"
35#include "G4GenericIon.hh"
36#include "CLHEP/Random/Random.h"
37
38G4InclLightIonInterface::G4InclLightIonInterface()
39{
40  hazard = new G4Hazard();
41  const G4long* table_entry = CLHEP::HepRandom::getTheSeeds(); // Get random seed from CLHEP.
42  hazard->ial = (*table_entry);
43
44  varntp = new G4VarNtp();
45  calincl = new G4Calincl();
46  ws = new G4Ws();
47  mat = new G4Mat();
48  incl = new G4Incl(hazard, calincl, ws, mat, varntp);
49
50  verboseLevel = 0;
51}
52
53G4InclLightIonInterface::~G4InclLightIonInterface()
54{
55  delete hazard;
56  delete varntp;
57  delete calincl;
58  delete ws;
59  delete mat;
60  delete incl;
61}
62
63G4HadFinalState* G4InclLightIonInterface::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
64{
65  G4int maxTries = 200;
66
67  G4int particleI, n = 0;
68
69  G4int bulletType = 0;
70
71  // Print diagnostic messages: 0 = silent, 1 and 2 = verbose
72  verboseLevel = 0;
73
74  // Increase the event number:
75  eventNumber++;
76
77  if (verboseLevel > 1) {
78    G4cout << " >>> G4InclLightIonInterface::ApplyYourself called" << G4endl;
79  }
80
81  if(verboseLevel > 1) {
82    G4cout <<"G4InclLightIonInterface: Now processing INCL4 event number:" << eventNumber << G4endl;
83  }
84
85  // INCL4 needs the energy in units MeV
86  G4double bulletE = aTrack.GetKineticEnergy() / MeV;
87
88  G4double targetA = theNucleus.GetN();
89  G4double targetZ = theNucleus.GetZ();
90
91  G4double eKin;
92  G4double momx = 0.0, momy = 0.0, momz = 0.0;
93  G4DynamicParticle *cascadeParticle = 0;
94  G4ParticleDefinition *aParticleDefinition = 0;
95
96  // INCL assumes the projectile particle is going in the direction of
97  // the Z-axis. Here we construct proper rotation to convert the
98  // momentum vectors of the outcoming particles to the original
99  // coordinate system.
100  G4LorentzVector projectileMomentum = aTrack.Get4Momentum();
101  G4LorentzRotation toZ;
102  toZ.rotateZ(-projectileMomentum.phi());
103  toZ.rotateY(-projectileMomentum.theta());
104  G4LorentzRotation toLabFrame = toZ.inverse();
105
106  theResult.Clear(); // Make sure the output data structure is clean.
107
108  // Map Geant4 particle types to corresponding INCL4 types.
109  enum bulletParticleType {nucleus = 0, proton = 1, neutron = 2, pionPlus = 3, pionZero = 4, 
110                           pionMinus = 5, deuteron = 6, triton = 7, he3 = 8, he4 = 9};
111
112  // Coding particles for use with INCL4 and ABLA
113  if (aTrack.GetDefinition() == G4Deuteron::Deuteron()   ) bulletType = deuteron;
114  if (aTrack.GetDefinition() == G4Triton::Triton()       ) bulletType = triton;
115  if (aTrack.GetDefinition() == G4He3::He3()             ) bulletType = he3;
116  if (aTrack.GetDefinition() == G4Alpha::Alpha()         ) bulletType = he4;
117
118  for(int i = 0; i < 15; i++) {
119    calincl->f[i] = 0.0; // Initialize INCL input data
120  }
121
122  // Check wheter the input is acceptable.
123  if((bulletType != 0) && ((targetA != 1) && (targetZ != 1))) {
124    calincl->f[0] = targetA;    // Target mass number
125    calincl->f[1] = targetZ;    // Charge number
126    calincl->f[6] = bulletType; // Type
127    calincl->f[2] = bulletE;    // Energy [MeV]
128    calincl->f[5] = 1.0;        // Time scaling
129    calincl->f[4] = 45.0;       // Nuclear potential
130
131    ws->nosurf = -2;  // Nucleus surface, -2 = Woods-Saxon
132    ws->xfoisa = 8;
133    ws->npaulstr = 0;
134
135    int nTries = 0;
136    varntp->ntrack = 0;
137
138    mat->nbmat = 1;
139    mat->amat[0] = int(calincl->f[0]);
140    mat->zmat[0] = int(calincl->f[1]);
141
142    incl->initIncl(true);
143
144    while((varntp->ntrack <= 0) && (nTries < maxTries)) { // Loop until we produce real cascade
145      nTries++;
146      if(verboseLevel > 1) {
147        G4cout <<"G4InclLightIonInterface: Try number = " << nTries << G4endl; 
148      }
149      incl->processEventIncl();
150
151      if(verboseLevel > 1) {
152        G4cout <<"G4InclLightIonInterface: number of tracks = " <<  varntp->ntrack <<G4endl;
153      }
154    }
155
156    if(verboseLevel > 1) {
157      /**
158       * Diagnostic output
159       */
160      G4cout <<"G4InclLightIonInterface: Bullet type: " << bulletType << G4endl;
161      G4cout <<"G4Incl4AblaCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
162
163      G4cout <<"G4InclLightIonInterface: Target A:  " << targetA << G4endl;
164      G4cout <<"G4InclLightIonInterface: Target Z:  " << targetZ << G4endl;
165
166      if(verboseLevel > 3) {
167        diagdata <<"G4InclLightIonInterface: Bullet type: " << bulletType << G4endl;
168        diagdata <<"G4InclLightIonInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
169       
170        diagdata <<"G4InclLightIonInterface: Target A:  " << targetA << G4endl;
171        diagdata <<"G4InclLightIonInterface: Target Z:  " << targetZ << G4endl;
172      }
173
174      for(particleI = 0; particleI < varntp->ntrack; particleI++) {
175        G4cout << n                       << " " << calincl->f[6]             << " " << calincl->f[2] << " ";
176        G4cout << varntp->massini         << " " << varntp->mzini             << " ";
177        G4cout << varntp->exini           << " " << varntp->mulncasc          << " " << varntp->mulnevap          << " " << varntp->mulntot << " ";
178        G4cout << varntp->bimpact         << " " << varntp->jremn             << " " << varntp->kfis              << " " << varntp->estfis << " ";
179        G4cout << varntp->izfis           << " " << varntp->iafis             << " " << varntp->ntrack            << " " << varntp->itypcasc[particleI] << " ";
180        G4cout << varntp->avv[particleI]  << " " << varntp->zvv[particleI]    << " " << varntp->enerj[particleI]  << " ";
181        G4cout << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
182        // For diagnostic output
183        if(verboseLevel > 3) {
184          diagdata << n                       << " " << calincl->f[6]             << " " << calincl->f[2] << " ";
185          diagdata << varntp->massini         << " " << varntp->mzini             << " ";
186          diagdata << varntp->exini           << " " << varntp->mulncasc          << " " << varntp->mulnevap          << " " << varntp->mulntot << " ";
187          diagdata << varntp->bimpact         << " " << varntp->jremn             << " " << varntp->kfis              << " " << varntp->estfis << " ";
188          diagdata << varntp->izfis           << " " << varntp->iafis             << " " << varntp->ntrack            << " ";
189          diagdata                                                                       << varntp->itypcasc[particleI] << " ";
190          diagdata << varntp->avv[particleI]  << " " << varntp->zvv[particleI]    << " " << varntp->enerj[particleI]  << " ";
191          diagdata << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
192        }
193      }
194    }
195
196    // Check whether a valid cascade was produced.
197    // If not return the original bullet particle with the same momentum.
198    if(varntp->ntrack <= 0) {
199      if(verboseLevel > 1) {
200        G4cout <<"WARNING G4InclLightIonInterface: No cascade. Returning original particle with original momentum." << G4endl;
201        G4cout <<"\t Reached maximum trials of 200 to produce inelastic scattering." << G4endl;
202      }
203
204      theResult.SetStatusChange(stopAndKill);
205     
206      if(bulletType == proton) {
207        aParticleDefinition = G4Proton::ProtonDefinition();
208      }
209      if(bulletType == neutron) {
210        aParticleDefinition = G4Neutron::NeutronDefinition();
211      }
212      if(bulletType == pionPlus) {
213        aParticleDefinition = G4PionPlus::PionPlusDefinition();
214      }
215      if(bulletType == pionZero) {
216        aParticleDefinition = G4PionZero::PionZeroDefinition();
217      }
218      if(bulletType == pionMinus) {
219        aParticleDefinition = G4PionMinus::PionMinusDefinition();
220      }
221
222      cascadeParticle = new G4DynamicParticle();
223      cascadeParticle->SetDefinition(aParticleDefinition);
224      cascadeParticle->Set4Momentum(aTrack.Get4Momentum());
225      theResult.AddSecondary(cascadeParticle);
226    }
227
228    // Convert INCL4 output to Geant4 compatible data structures.
229    // Elementary particles are converted to G4DynamicParticle.
230    theResult.SetStatusChange(stopAndKill);
231   
232    for(particleI = 0; particleI < varntp->ntrack; particleI++) { // Loop through the INCL4+ABLA output.
233      // Get energy/momentum and construct momentum vector in INCL4 coordinates.
234      momx = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::cos(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
235      momy = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::sin(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
236      momz = varntp->plab[particleI]*std::cos(varntp->tetlab[particleI]*CLHEP::pi/180.0)*MeV;
237
238      eKin = varntp->enerj[particleI] * MeV;
239
240      G4ThreeVector momDirection(momx, momy, momz); // Direction of the particle.
241      momDirection = momDirection.unit();
242      if(verboseLevel > 2) {
243        G4cout <<"G4InclLightIonInterface: " << G4endl;
244        G4cout <<"A    = " << varntp->avv[particleI] << " Z = "  << varntp->zvv[particleI] << G4endl;
245        G4cout <<"eKin = " << eKin                   << " MeV"   << G4endl;
246        G4cout <<"px   = " << momDirection.x()       << " py = " << momDirection.y()       <<" pz = " << momDirection.z() << G4endl;
247      }
248
249      G4int particleIdentified = 0; // Check particle ID.
250
251      if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 1)) { // Proton
252        cascadeParticle = 
253          new G4DynamicParticle(G4Proton::ProtonDefinition(), momDirection, eKin);
254        particleIdentified++;
255      }
256
257      if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 0)) { // Neutron
258        cascadeParticle = 
259          new G4DynamicParticle(G4Neutron::NeutronDefinition(), momDirection, eKin);
260        particleIdentified++;
261      }
262
263      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 1)) { // PionPlus
264        cascadeParticle = 
265          new G4DynamicParticle(G4PionPlus::PionPlusDefinition(), momDirection, eKin);
266        particleIdentified++;
267      }
268
269      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 0)) { // PionZero
270        cascadeParticle = 
271          new G4DynamicParticle(G4PionZero::PionZeroDefinition(), momDirection, eKin);
272        particleIdentified++;
273      }
274
275      if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == -1)) { // PionMinus
276        cascadeParticle = 
277          new G4DynamicParticle(G4PionMinus::PionMinusDefinition(), momDirection, eKin);
278        particleIdentified++;
279      }
280
281      if((varntp->avv[particleI] > 1) && (varntp->zvv[particleI] >= 1)) { // Nucleus fragment
282        G4ParticleDefinition * aIonDef = 0;
283        G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
284
285        G4int A = G4int(varntp->avv[particleI]);
286        G4int Z = G4int(varntp->zvv[particleI]);
287        G4double excitationE = G4double(varntp->exini) * MeV;
288
289        if(verboseLevel > 1) {
290          G4cout <<"Finding ion: A = " << A << " Z = " << Z << " E* = " << excitationE/MeV << G4endl;
291        }
292        aIonDef = theTableOfParticles->GetIon(Z, A, excitationE);
293       
294        if(aIonDef == 0) {
295          if(verboseLevel > 1) {
296            G4cout <<"G4InclLightIonInterface: " << G4endl;
297            G4cout <<"FATAL ERROR: aIonDef = 0" << G4endl;
298            G4cout <<"A = " << A << " Z = " << Z << " E* = " << excitationE << G4endl;
299          }
300        }
301
302        if(aIonDef != 0) { // If the ion was identified add it to output.
303          cascadeParticle =
304            new G4DynamicParticle(aIonDef, momDirection, eKin);
305          particleIdentified++;
306        }
307      }
308       
309      if(particleIdentified == 1) { // Particle identified properly.
310        cascadeParticle->Set4Momentum(cascadeParticle->Get4Momentum()*=toLabFrame);
311        theResult.AddSecondary(cascadeParticle); // Put data into G4HadFinalState.
312      }
313      else { // Particle identification failed.
314        if(particleIdentified > 1) { // Particle was identified as more than one particle type.
315          if(verboseLevel > 1) {
316            G4cout <<"G4InclLightIonInterface: One outcoming particle was identified as";
317            G4cout <<"more than one particle type. This is probably due to a bug in the interface." << G4endl;
318            G4cout <<"Particle A:" << varntp->avv[particleI] << "Z: " << varntp->zvv[particleI] << G4endl;
319            G4cout << "(particleIdentified =" << particleIdentified << ")"  << G4endl;
320          }
321        }
322      }
323    }
324
325    varntp->ntrack = 0; // Clean up the number of generated particles in the event.
326  }
327  /**
328   * Report unsupported features.
329   * (Check bullet, target, energy range)
330   */
331  else { // If the bullet type was not recognized by the interface, it will be returned back without any interaction.
332    theResult.SetStatusChange(stopAndKill);
333
334    G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
335    cascadeParticle = new G4DynamicParticle(theTableOfParticles->FindParticle(aTrack.GetDefinition()), aTrack.Get4Momentum());
336
337    theResult.AddSecondary(cascadeParticle);
338
339    if(verboseLevel > 1) {
340      G4cout <<"G4InclLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
341    }
342    if(verboseLevel > 3) {
343      diagdata <<"G4InclLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
344    }
345
346    if(bulletType == 0) {
347      if(verboseLevel > 1) {
348        G4cout <<"G4InclLightIonInterface: Unknown bullet type" << G4endl;     
349        G4cout <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
350      }     
351      if(verboseLevel > 3) {
352        diagdata <<"G4InclLightIonInterface: Unknown bullet type" << G4endl;     
353        diagdata <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
354      }
355    }
356
357    if((targetA == 1) && (targetZ == 1)) { // Unsupported target
358      if(verboseLevel > 1) {
359        G4cout <<"Unsupported target: " << G4endl;
360        G4cout <<"Target A: " << targetA << G4endl;
361        G4cout <<"TargetZ: " << targetZ << G4endl;
362      }
363      if(verboseLevel > 3) {
364        diagdata <<"Unsupported target: " << G4endl;
365        diagdata <<"Target A: " << targetA << G4endl;
366       diagdata <<"TargetZ: " << targetZ << G4endl;
367      }
368    }
369
370    if(bulletE < 100) { // INCL does not support E < 100 MeV.
371      if(verboseLevel > 1) {
372        G4cout <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
373        G4cout <<"WARNING: Returning the original bullet with original energy back to Geant4." << G4endl;
374      }
375      if(verboseLevel > 3) {
376        diagdata <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
377      }
378    }
379    if(verboseLevel > 3) {
380      diagdata <<"WARNING: returning the original bullet with original energy back to Geant4." << G4endl;
381    }
382  }
383
384  return &theResult;
385} 
386
387G4ReactionProductVector* G4InclLightIonInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) {
388  return 0;
389}
390
391
Note: See TracBrowser for help on using the repository browser.