source: trunk/source/processes/hadronic/models/cascade/cascade/src/G4PreCompoundCascadeInterface.cc @ 1315

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 10.9 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// $Id: G4PreCompoundCascadeInterface.cc,v 1.13 2010/05/21 18:07:30 mkelsey Exp $
26// Geant4 tag: $Name: geant4-09-04-beta-cand-01 $
27//
28// 20100114  M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29// 20100413  M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
30// 20100414  M. Kelsey -- Check for K0L/K0S before using G4InuclElemPart::type
31// 20100419  M. Kelsey -- Access G4CollisionOutput lists by const-ref, and
32//              const_iterator
33// 20100428  M. Kelsey -- Use G4InuclParticleNames enum
34// 20100429  M. Kelsey -- Change "case gamma:" to "case photon:"
35// 20100517  M. Kelsey -- Follow new ctors for G4*Collider family.
36// 20100520  M. Kelsey -- Add missing name string to ctor, follow code changes
37//              from G4CascadeInterface.
38
39#include "G4PreCompoundCascadeInterface.hh"
40#include "globals.hh"
41#include "G4CollisionOutput.hh"
42#include "G4DynamicParticle.hh"
43#include "G4InuclElementaryParticle.hh"
44#include "G4InuclNuclei.hh"
45#include "G4InuclParticle.hh"
46#include "G4InuclParticleNames.hh"
47#include "G4KaonZeroShort.hh"
48#include "G4KaonZeroLong.hh"
49#include "G4LorentzRotation.hh"
50#include "G4Nucleus.hh"
51#include "G4ParticleDefinition.hh"
52#include "G4PreCompoundInuclCollider.hh"
53#include "G4Track.hh"
54#include "G4V3DNucleus.hh"
55
56using namespace G4InuclParticleNames;
57
58
59typedef std::vector<G4InuclElementaryParticle>::const_iterator particleIterator;
60typedef std::vector<G4InuclNuclei>::const_iterator nucleiIterator;
61
62G4PreCompoundCascadeInterface::G4PreCompoundCascadeInterface(const G4String& nam)
63  :G4VIntraNuclearTransportModel(nam), verboseLevel(0)  {
64
65  if (verboseLevel > 3) {
66    G4cout << " >>> G4PreCompoundCascadeInterface::G4PreCompoundCascadeInterface" << G4endl;
67  }
68}
69   
70G4ReactionProductVector* G4PreCompoundCascadeInterface::Propagate(G4KineticTrackVector* , 
71                                                                  G4V3DNucleus* ) {
72  return 0;
73}
74
75// #define debug_G4PreCompoundCascadeInterface
76
77G4HadFinalState* G4PreCompoundCascadeInterface::ApplyYourself(const G4HadProjectile& aTrack, 
78                                                              G4Nucleus& theNucleus) {
79#ifdef debug_G4PreCompoundCascadeInterface
80  static G4int counter(0);
81  counter++;
82  G4cerr << "Reaction number "<< counter << " "<<aTrack.GetDynamicParticle()->GetDefinition()->GetParticleName()<<" "<< aTrack.GetDynamicParticle()->GetKineticEnergy()<<G4endl;
83#endif
84
85  theResult.Clear();
86
87  if (verboseLevel > 3) {
88    G4cout << " >>> G4PreCompoundCascadeInterface::ApplyYourself" << G4endl;
89  };
90
91  G4double eInit     = 0.0;
92  G4double eTot      = 0.0;
93  G4double sumBaryon = 0.0;
94  G4double sumEnergy = 0.0;
95
96  // Make conversion between native Geant4 and Bertini cascade classes.
97  // NOTE: Geant4 units are MeV = 1 and GeV = 1000. Cascade code by default use GeV = 1.
98
99  G4int bulletType;
100  if (aTrack.GetDefinition() == G4KaonZeroLong::KaonZeroLong() ||
101      aTrack.GetDefinition() == G4KaonZeroShort::KaonZeroShort() )
102    bulletType = (G4UniformRand() > 0.5) ? kaonZero : kaonZeroBar;
103  else 
104    bulletType = G4InuclElementaryParticle::type(aTrack.GetDefinition());
105
106  // Code momentum and energy.
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  G4LorentzVector momentumBullet(0., 0., aTrack.GetTotalMomentum()/GeV,
114                                 aTrack.GetTotalEnergy()/GeV);
115
116  G4InuclElementaryParticle *  bullet =
117    new G4InuclElementaryParticle(momentumBullet, bulletType); 
118
119  sumEnergy = bullet->getKineticEnergy(); // In GeV
120  sumBaryon += bullet->baryon();
121
122  // Set target
123  G4InuclNuclei*   target  = 0;
124  G4InuclParticle* targetH = 0;
125
126  G4double theNucleusA = theNucleus.GetN();
127
128  if ( !(G4int(theNucleusA) == 1) ) {
129    target  = new G4InuclNuclei(theNucleusA, theNucleus.GetZ());
130    eInit = bullet->getEnergy() + target->getEnergy();
131
132    sumBaryon += theNucleusA;
133
134    if (verboseLevel > 2) {
135      G4cout << "Bullet:  " << G4endl; 
136      bullet->printParticle();
137    }
138    if (verboseLevel > 2) {
139      G4cout << "Target:  " << G4endl; 
140      target->printParticle();
141    }
142  }
143
144  G4CollisionOutput output;
145
146  // Colliders initialisation
147  G4PreCompoundInuclCollider*  collider = new G4PreCompoundInuclCollider;
148
149  G4int  maxTries = 10; // maximum tries for inelastic collision to avoid infinite loop
150  G4int  nTries   = 0;  // try counter
151
152  if (G4int(theNucleusA) == 1) { // special treatment for target H(1,1) (proton)
153
154    targetH = new G4InuclElementaryParticle(1);
155
156    G4float cutElastic[32];
157
158    cutElastic[proton   ]   = 1.0; // 1 GeV
159    cutElastic[neutron  ]   = 1.0;
160    cutElastic[lambda]      = 1.0;
161    cutElastic[sigmaPlus]   = 1.0;
162    cutElastic[sigmaZero]   = 1.0;
163    cutElastic[sigmaMinus]  = 1.0;
164    cutElastic[xiZero]      = 1.0;
165    cutElastic[xiMinus]     = 1.0;
166
167    cutElastic[pionPlus ]   = 0.6; // 0.6 GeV
168
169    cutElastic[kaonPlus ]   = 0.5; // 0.5 GeV
170    cutElastic[kaonMinus]   = 0.5;
171    cutElastic[kaonZero]    = 0.5;
172    cutElastic[kaonZeroBar] = 0.5;
173
174    cutElastic[pionMinus]   = 0.2; // 0.2 GeV
175    cutElastic[pionZero ]   = 0.2;
176
177
178    if (momentumBullet.z() > cutElastic[bulletType]) { // inelastic collision possible
179
180      do {   // we try to create inelastic interaction
181        output.reset();
182        collider->collide(bullet, targetH, output);
183        nTries++;
184      } while(
185              (nTries < maxTries) &&
186              (output.getOutgoingParticles().size() == 2 && // elastic: bullet + p = H(1,1) coming out
187               (output.getOutgoingParticles().begin()->type() == bulletType || output.getOutgoingParticles().begin()->type() == proton)
188               )
189              );
190
191    } else { // only elastic collision is energetically possible
192      collider->collide(bullet, targetH, output);
193    }
194
195    sumBaryon += 1;
196
197    eInit = bullet->getEnergy() + target->getEnergy();
198
199    if (verboseLevel > 2) {
200      G4cout << "Target:  " << G4endl;
201      targetH->printParticle();
202    }
203
204  } else {  // treat all other targets excepet H(1,1)
205
206    do  // we try to create inelastic interaction
207      {
208        output.reset();
209        collider->collide(bullet, target, output);
210        nTries++;
211      } while (
212               (nTries < maxTries) &&
213               output.getOutgoingParticles().size() == 1 &&     // we retry when elastic collision happened
214               output.getNucleiFragments().size() == 1 &&           
215               output.getOutgoingParticles().begin()->type() == bullet->type() &&
216               output.getNucleiFragments().begin()->getA() == target->getA() && 
217               output.getNucleiFragments().begin()->getZ() == target->getZ() 
218               );
219  }
220 
221  if (verboseLevel > 1) {
222      G4cout << " Cascade output: " << G4endl;
223      output.printCollisionOutput();
224  }
225 
226  // Rotate event to put Z axis along original projectile direction
227  output.rotateEvent(toLabFrame);
228
229  // Convert cascade data to use hadronics interface
230  const std::vector<G4InuclNuclei>& nucleiFragments = output.getNucleiFragments();
231  const std::vector<G4InuclElementaryParticle>& particles = output.getOutgoingParticles();
232
233  theResult.SetStatusChange(stopAndKill);
234
235  // Get outcoming particles
236  G4DynamicParticle* cascadeParticle = 0;
237  if (!particles.empty()) { 
238    particleIterator ipart = particles.begin();
239    for (; ipart != particles.end(); ipart++) {
240      G4int outgoingType = ipart->type();
241
242      eTot += ipart->getEnergy();
243      sumBaryon -= ipart->baryon();
244      sumEnergy -= ipart->getKineticEnergy();
245
246      if (!ipart->valid() || ipart->quasi_deutron()) {
247        G4cerr << " ERROR: G4PreCompoundCascadeInterface::Propagate incompatible"
248               << " particle type " << ipart->type() << G4endl;
249        continue;
250      }
251
252      // Copy local G4DynPart to public output (handle kaon mixing specially)
253      if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
254        G4ThreeVector momDir = ipart->getMomentum().vect().unit();
255        G4double ekin = ipart->getKineticEnergy();
256
257        G4ParticleDefinition* pd = G4KaonZeroShort::Definition();
258        if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
259
260        cascadeParticle = new G4DynamicParticle(pd, momDir, ekin);
261      } else {
262        cascadeParticle = new G4DynamicParticle(ipart->getDynamicParticle());
263      }
264 
265      theResult.AddSecondary(cascadeParticle); 
266    }
267  }
268
269  // get nuclei fragments
270  G4DynamicParticle * aFragment = 0;
271  if (!nucleiFragments.empty()) { 
272    nucleiIterator ifrag = nucleiFragments.begin();
273    for (; ifrag != nucleiFragments.end(); ifrag++) {
274      eTot += ifrag->getEnergy();
275      sumBaryon -= ifrag->getA();
276      sumEnergy -= ifrag->getKineticEnergy();
277
278      if (verboseLevel > 2) {
279        G4cout << " Nuclei fragment: " << G4endl;
280        ifrag->printParticle();
281      }
282
283      // Copy local G4DynPart to public output
284      aFragment =  new G4DynamicParticle(ifrag->getDynamicParticle());
285      theResult.AddSecondary(aFragment); 
286    }
287  }
288
289  // Report violations of energy, baryon conservation
290  if (verboseLevel > 2) {
291    if (sumBaryon != 0) {
292      G4cout << "ERROR: no baryon number conservation, sum of baryons = "
293             << sumBaryon << G4endl;
294    }
295
296    if (sumEnergy > 0.01 ) {
297      G4cout << "Kinetic energy conservation violated by "
298             << sumEnergy << " GeV" << G4endl;
299    }
300     
301    G4cout << "Total energy conservation at level ~"
302           << (eInit - eTot) * GeV << " MeV" << G4endl;
303   
304    if (sumEnergy < -5.0e-5 ) { // 0.05 MeV
305      G4cout << "FATAL ERROR: energy created  "
306             << sumEnergy * GeV << " MeV" << G4endl;
307    }
308  }
309
310  delete bullet;
311  delete collider;
312
313  if(target != 0) delete target;
314  if(targetH != 0) delete targetH;
315
316  return &theResult;
317  }
Note: See TracBrowser for help on using the repository browser.