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

Last change on this file since 1245 was 1228, checked in by garnier, 15 years ago

update geant4.9.3 tag

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