source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QStringChipsParticleLevelInterface.cc @ 1228

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

update geant4.9.3 tag

File size: 26.6 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// Short description: Interface of QGSC to CHIPS (EnergyFlow)
27// ----------------------------------------------------------------------------
28
29//#define debug
30//#define pdebug
31
32#include "G4QStringChipsParticleLevelInterface.hh"
33#include "globals.hh"
34#include <utility>
35#include <list>
36#include <vector>
37#include "G4KineticTrackVector.hh"
38#include "G4Nucleon.hh"
39#include "G4Proton.hh"
40#include "G4Neutron.hh"
41#include "G4LorentzRotation.hh"
42#include "G4HadronicException.hh"
43// #define CHIPSdebug
44// #define CHIPSdebug_1
45
46G4QStringChipsParticleLevelInterface::G4QStringChipsParticleLevelInterface()
47{
48#ifdef debug
49  G4cout<<"G4QStringChipsParticleLevelInterface::Constructor is called"<<G4endl;
50#endif
51  theEnergyLossPerFermi = 1.*GeV;
52  nop = 152;                               // clusters (A<6)
53  fractionOfSingleQuasiFreeNucleons = 0.5; // It is A-dependent (C=.85, U=.40)
54  fractionOfPairedQuasiFreeNucleons = 0.05;
55  clusteringCoefficient = 5.;
56  temperature = 180.;
57  halfTheStrangenessOfSee = 0.3; // = s/d = s/u
58  etaToEtaPrime = 0.3;
59  fusionToExchange = 1.;
60  theInnerCoreDensityCut = 50.;
61 
62  if(getenv("ChipsParameterTuning"))
63  {
64    G4cout << "Please enter the energy loss per fermi in GeV"<<G4endl;
65    G4cin >> theEnergyLossPerFermi;
66    theEnergyLossPerFermi *= GeV;
67    G4cout << "Please enter nop"<<G4endl;
68    G4cin >> nop;
69    G4cout << "Please enter the fractionOfSingleQuasiFreeNucleons"<<G4endl;
70    G4cin >> fractionOfSingleQuasiFreeNucleons;
71    G4cout << "Please enter the fractionOfPairedQuasiFreeNucleons"<<G4endl;
72    G4cin >> fractionOfPairedQuasiFreeNucleons;
73    G4cout << "Please enter the clusteringCoefficient"<<G4endl;
74    G4cin >> clusteringCoefficient;
75    G4cout << "Please enter the temperature"<<G4endl;
76    G4cin >> temperature;
77    G4cout << "Please enter halfTheStrangenessOfSee"<<G4endl;
78    G4cin >> halfTheStrangenessOfSee;
79    G4cout << "Please enter the etaToEtaPrime"<<G4endl;
80    G4cin >> etaToEtaPrime;
81    G4cout << "Please enter the fusionToExchange"<<G4endl;
82    G4cin >> fusionToExchange;
83    G4cout << "Please enter cut-off for calculating the nuclear radius in percent"<<G4endl;
84    G4cin >> theInnerCoreDensityCut;
85  }
86}
87
88G4HadFinalState* G4QStringChipsParticleLevelInterface::
89ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
90{
91#ifdef debug
92  G4cout<<"G4QStringChipsParticleLevelInterface::ApplyYourself is called"<<G4endl;
93#endif
94  return theModel.ApplyYourself(aTrack, theNucleus);
95}
96
97G4ReactionProductVector* G4QStringChipsParticleLevelInterface::
98Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus)
99{
100  static const G4double mProt=G4Proton::Proton()->GetPDGMass();
101  static const G4double mNeut=G4Neutron::Neutron()->GetPDGMass();
102  static const G4double mLamb=G4Lambda::Lambda()->GetPDGMass();
103#ifdef debug
104  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate is called"<<G4endl;
105#endif
106  // Protection for non physical conditions
107 
108  if(theSecondaries->size() == 1) 
109  {
110    G4ReactionProductVector* theFastResult = new G4ReactionProductVector;
111    G4ReactionProduct* theFastSec;
112    theFastSec = new G4ReactionProduct((*theSecondaries)[0]->GetDefinition());
113    G4LorentzVector current4Mom = (*theSecondaries)[0]->Get4Momentum();
114    theFastSec->SetTotalEnergy(current4Mom.t());
115    theFastSec->SetMomentum(current4Mom.vect());
116    theFastResult->push_back(theFastSec);
117    return theFastResult;
118    //throw G4HadronicException(__FILE__,__LINE__,
119    //      "G4QStringChipsParticleLevelInterface: Only one particle from String models!");
120  }
121 
122  // target properties needed in constructor of quasmon, and for boosting to
123  // target rest frame
124  // remove all nucleons already involved in STRING interaction, to make the ResidualTarget
125  theNucleus->StartLoop();
126  G4Nucleon * aNucleon;
127  G4int resA = 0;
128  G4int resZ = 0;
129  G4ThreeVector hitMomentum(0,0,0);
130  G4double hitMass = 0;
131  unsigned int hitCount = 0;
132  while((aNucleon = theNucleus->GetNextNucleon()))
133  {
134    if(!aNucleon->AreYouHit())
135    {
136      resA++;                                                  // Collect A of the ResidNuc
137      resZ+=G4int (aNucleon->GetDefinition()->GetPDGCharge()); // Collect Z of the ResidNuc
138    }
139    else
140    {
141      hitMomentum += aNucleon->GetMomentum().vect();           // Sum 3-mom of StringHadr's
142      hitMass += aNucleon->GetMomentum().m();                  // Sum masses of StringHadrs
143      hitCount ++;                                             // Calculate STRING hadrons
144    }
145  }
146  G4int targetPDGCode = 90000000 + 1000*resZ + (resA-resZ);    // PDG of theResidualNucleus
147  G4double targetMass=mNeut;
148  if (!resZ)                                                   // Nucleus of only neutrons
149  {
150    if (resA>1) targetMass*=resA;
151  }
152  else targetMass=G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ)
153                                                                            ->GetPDGMass();
154  G4double targetEnergy = std::sqrt(hitMomentum.mag2()+targetMass*targetMass);
155  // !! @@ Target should be at rest: hitMomentum=(0,0,0) @@ !! M.K. (go to this system)
156  G4LorentzVector targ4Mom(-1.*hitMomentum, targetEnergy);
157 
158  // Calculate the mean energy lost
159  std::pair<G4double, G4double> theImpact = theNucleus->RefetchImpactXandY();
160  G4double impactX = theImpact.first;
161  G4double impactY = theImpact.second;
162  G4double impactPar2 = impactX*impactX + impactY*impactY;
163 
164  G4double radius2 = theNucleus->GetNuclearRadius(theInnerCoreDensityCut*perCent);
165  radius2 *= radius2;
166  G4double pathlength = 0.;
167#ifdef pdebug
168  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: r="<<std::sqrt(radius2)/fermi
169        <<", b="<<std::sqrt(impactPar2)/fermi<<", R="<<theNucleus->GetOuterRadius()/fermi
170        <<", b/r="<<std::sqrt(impactPar2/radius2)<<G4endl; 
171#endif
172  if(radius2 - impactPar2>0) pathlength = 2.*std::sqrt(radius2 - impactPar2);
173  G4double theEnergyLostInFragmentation = theEnergyLossPerFermi*pathlength/fermi;
174 
175  // now select all particles in range
176  std::list<std::pair<G4double, G4KineticTrack *> > theSorted;         // Output
177  std::list<std::pair<G4double, G4KineticTrack *> >::iterator current; // Input
178  for(unsigned int secondary = 0; secondary<theSecondaries->size(); secondary++)
179  {
180    G4LorentzVector a4Mom = theSecondaries->operator[](secondary)->Get4Momentum();
181#ifdef CHIPSdebug
182    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: ALL STRING particles "
183          << theSecondaries->operator[](secondary)->GetDefinition()->GetPDGCharge()<<" "
184          << theSecondaries->operator[](secondary)->GetDefinition()->GetPDGEncoding()<<" "
185          << a4Mom <<G4endl; 
186#endif
187#ifdef pdebug
188    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: in C="
189          <<theSecondaries->operator[](secondary)->GetDefinition()->GetPDGCharge()<<",PDG="
190          <<theSecondaries->operator[](secondary)->GetDefinition()->GetPDGEncoding()
191          <<",4M="<<a4Mom<<", current nS="<<theSorted.size()<<G4endl; 
192#endif
193    G4double toSort = a4Mom.rapidity();          // Rapidity is used for the ordering (?!)
194    std::pair<G4double, G4KineticTrack *> it;
195    it.first = toSort;
196    it.second = theSecondaries->operator[](secondary);
197    G4bool inserted = false;
198    for(current = theSorted.begin(); current!=theSorted.end(); current++)
199    {
200      if((*current).first > toSort)        // The current is smaller then existing
201      {
202        theSorted.insert(current, it);     // It shifts the others up
203        inserted = true;
204        break;
205      }
206    }
207    if(!inserted) theSorted.push_back(it); // It is bigger than any previous
208  }
209 
210  G4LorentzVector proj4Mom(0.,0.,0.,0.);
211  G4int nD  = 0;
212  G4int nU  = 0;
213  G4int nS  = 0;
214  G4int nAD = 0;
215  G4int nAU = 0;
216  G4int nAS = 0;
217  std::list<std::pair<G4double,G4KineticTrack*> >::iterator firstEscape=theSorted.begin();
218  G4double runningEnergy = 0;
219  G4int particleCount = 0;
220  G4LorentzVector theLow = (*(theSorted.begin())).second->Get4Momentum();
221  G4LorentzVector theHigh;
222
223#ifdef CHIPSdebug
224  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: CHIPS ENERGY LOST "
225        <<theEnergyLostInFragmentation<<". Sorted rapidities event start"<<G4endl;
226#endif
227#ifdef pdebug
228  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: total CHIPS energy = "
229        <<theEnergyLostInFragmentation<<". Start rapidity sorting nS="<<theSorted.size()
230        <<G4endl;
231#endif
232
233  G4QHadronVector projHV;
234  std::vector<G4QContent> theContents;
235  std::vector<G4LorentzVector*> theMomenta;
236  G4ReactionProductVector* theResult = new G4ReactionProductVector;
237  G4ReactionProduct* theSec;
238  G4KineticTrackVector* secondaries;
239#ifdef pdebug
240  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: Absorption nS="
241        <<theSorted.size()<<G4endl; 
242#endif
243  G4bool EscapeExists = false;
244  for(current = theSorted.begin(); current!=theSorted.end(); current++)
245  {
246#ifdef pdebug
247    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: nq="
248          <<(*current).second->GetDefinition()->GetQuarkContent(3)<<", naq="
249          <<(*current).second->GetDefinition()->GetAntiQuarkContent(3)<<", PDG="
250          <<(*current).second->GetDefinition()->GetPDGEncoding()<<",4M="
251          <<(*current).second->Get4Momentum()<<G4endl; 
252#endif
253    firstEscape = current;              // Remember to make decays for the rest
254    G4KineticTrack* aResult = (*current).second;
255    // This is an old (H.P.) solution, which makes an error in En/Mom conservation
256    //
257    // @@ Now it does not include strange particle for the absorption in nuclei (?!) M.K.
258    //if((*current).second->GetDefinition()->GetQuarkContent(3)!=0 ||
259    //   (*current).second->GetDefinition()->GetAntiQuarkContent(3) !=0) // Strange quarks
260    //{
261    //  G4ParticleDefinition* pdef = aResult->GetDefinition();
262    //  secondaries = NULL;
263    //  if ( pdef->GetPDGWidth() > 0 && pdef->GetPDGLifeTime() < 5E-17*s )
264    //     secondaries = aResult->Decay(); // @@ Decay of only strange resonances (?!) M.K.
265    //  if ( secondaries == NULL )           // No decay
266    //  {
267    //    theSec = new G4ReactionProduct(aResult->GetDefinition());
268    //    G4LorentzVector current4Mom = aResult->Get4Momentum();
269    //    current4Mom.boost(targ4Mom.boostVector()); // boost from the targetAtRes system
270    //    theSec->SetTotalEnergy(current4Mom.t());
271    //    theSec->SetMomentum(current4Mom.vect());
272    //    theResult->push_back(theSec);
273    //  }
274    //  else                                 // The decay happened
275    //  {
276    //   for (unsigned int aSecondary=0; aSecondary<secondaries->size(); aSecondary++)
277    //   {
278    //     theSec =
279    //        new G4ReactionProduct(secondaries->operator[](aSecondary)->GetDefinition());
280    //     G4LorentzVector current4Mom=secondaries->operator[](aSecondary)->Get4Momentum();
281    //     current4Mom.boost(targ4Mom.boostVector());
282    //     theSec->SetTotalEnergy(current4Mom.t());
283    //     theSec->SetMomentum(current4Mom.vect());
284    //     theResult->push_back(theSec);
285    //   }
286    //   std::for_each(secondaries->begin(), secondaries->end(), DeleteKineticTrack());
287    //   delete secondaries;
288    //  }
289    //}
290    //
291    //runningEnergy += (*current).second->Get4Momentum().t();
292    //if((*current).second->GetDefinition() == G4Proton::Proton())
293    //                                   runningEnergy-=G4Proton::Proton()->GetPDGMass();
294    //if((*current).second->GetDefinition() == G4Neutron::Neutron())
295    //                                   runningEnergy-=G4Neutron::Neutron()->GetPDGMass();
296    //if((*current).second->GetDefinition() == G4Lambda::Lambda())
297    //                                   runningEnergy-=G4Lambda::Lambda()->GetPDGMass();
298    //
299    // New solution starts from here (M.Kossov March 2006) [Strange particles included]
300    runningEnergy += aResult->Get4Momentum().t();
301    G4double charge=aResult->GetDefinition()->GetPDGCharge(); // Charge of the particle
302    G4int strang=aResult->GetDefinition()->GetQuarkContent(3);// Its strangeness
303    G4int baryn=aResult->GetDefinition()->GetBaryonNumber();  // Its baryon number
304    if     (baryn>0 && charge>0 && strang<1) runningEnergy-=mProt;  // For positive baryons
305    else if(baryn>0 && strang<1) runningEnergy-=mNeut;              // For neut/neg baryons
306    else if(baryn>0) runningEnergy-=mLamb;                          // For strange baryons
307    else if(baryn<0) runningEnergy+=mProt;                          // For anti-particles
308    // ------------ End of the new solution
309#ifdef CHIPSdebug
310    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: sorted rapidities "
311                                    <<(*current).second->Get4Momentum().rapidity()<<G4endl;
312#endif
313
314#ifdef pdebug
315    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: E="<<runningEnergy<<", EL="
316                                                    <<theEnergyLostInFragmentation<<G4endl;
317#endif
318    if(runningEnergy > theEnergyLostInFragmentation)
319    {
320      EscapeExists = true;
321      break;
322    }
323#ifdef CHIPSdebug
324    G4cout <<"G4QStringChipsParticleLevelInterface::Propagate: ABSORBED STRING particles "
325           <<(*current).second->GetDefinition()->GetPDGCharge()<<" "
326           << (*current).second->GetDefinition()->GetPDGEncoding()<<" "
327           << (*current).second->Get4Momentum() <<G4endl; 
328#endif
329#ifdef pdebug
330    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate:C="
331          <<current->second->GetDefinition()->GetPDGCharge()<<", PDG="
332          <<current->second->GetDefinition()->GetPDGEncoding()<<", 4M="
333          <<current->second->Get4Momentum()<<G4endl; 
334#endif
335
336    // projectile 4-momentum in target rest frame needed in constructor of QHadron
337    particleCount++;
338    theHigh = (*current).second->Get4Momentum(); 
339    proj4Mom = (*current).second->Get4Momentum(); 
340    proj4Mom.boost(-1.*targ4Mom.boostVector());   // Back to the system of nucleusAtRest
341    nD = (*current).second->GetDefinition()->GetQuarkContent(1);
342    nU = (*current).second->GetDefinition()->GetQuarkContent(2);
343    nS = (*current).second->GetDefinition()->GetQuarkContent(3);
344    nAD = (*current).second->GetDefinition()->GetAntiQuarkContent(1);
345    nAU = (*current).second->GetDefinition()->GetAntiQuarkContent(2);
346    nAS = (*current).second->GetDefinition()->GetAntiQuarkContent(3);
347    G4QContent aProjectile(nD, nU, nS, nAD, nAU, nAS);
348
349#ifdef CHIPSdebug_1
350    G4cout <<G4endl;
351    G4cout <<"G4QStringChipsParticleLevelInterface::Propagate: Quark content: d="<<nD
352           <<", u="<<nU<<", s="<<nS<< "Anti-quark content: anit-d="<<nAD<<", anti-u="<<nAU
353           <<", anti-s="<<nAS<<". G4QContent is constructed"<<endl;
354#endif
355
356    theContents.push_back(aProjectile);
357    G4LorentzVector* aVec = new G4LorentzVector((1./MeV)*proj4Mom); // @@ MeV is basic
358
359#ifdef CHIPSdebug_1
360    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: projectile momentum = "
361          <<*aVec<<G4endl;
362    G4cout << G4endl;
363#endif
364   
365    theMomenta.push_back(aVec);
366  }
367  std::vector<G4QContent> theFinalContents;
368  std::vector<G4LorentzVector*> theFinalMomenta;
369  if(theContents.size()<hitCount || 1) // Looks like the "else" is closed by "|| 1"
370  {
371    for(unsigned int hp = 0; hp<theContents.size(); hp++)
372    {
373      G4QHadron* aHadron = new G4QHadron(theContents[hp], *(theMomenta[hp]) );
374      projHV.push_back(aHadron);
375    }
376  }
377  else                                 // Never come here (!?)
378  {
379    unsigned int hp;
380    for(hp=0; hp<hitCount; hp++)       // Initialize the arrays
381    {
382      G4QContent co(0, 0, 0, 0, 0, 0);
383      theFinalContents.push_back(co);
384      G4LorentzVector* mo = new G4LorentzVector(0,0,0,0);
385      theFinalMomenta.push_back(mo);
386    }
387    unsigned int running = 0;
388    while (running<theContents.size())
389    {
390      for(hp = 0; hp<hitCount; hp++)
391      {
392        theFinalContents[hp] +=theContents[running];
393        *(theFinalMomenta[hp])+=*(theMomenta[running]);
394        running++;
395        if(running == theContents.size()) break;
396      }
397    }
398    for(hp = 0; hp<hitCount; hp++)
399    {
400      G4QHadron* aHadron = new G4QHadron(theFinalContents[hp], *theFinalMomenta[hp]);
401      projHV.push_back(aHadron);
402    }
403  }
404  // construct the quasmon
405  size_t i;
406  for (i=0; i<theFinalMomenta.size(); i++) delete theFinalMomenta[i];
407  for (i=0; i<theMomenta.size();      i++) delete theMomenta[i];
408  theFinalMomenta.clear();
409  theMomenta.clear();
410
411  G4QNucleus::SetParameters(fractionOfSingleQuasiFreeNucleons,
412                            fractionOfPairedQuasiFreeNucleons,
413                            clusteringCoefficient, 
414                            fusionToExchange);
415  G4Quasmon::SetParameters(temperature, halfTheStrangenessOfSee, etaToEtaPrime);
416
417#ifdef CHIPSdebug
418  G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: G4QNucleus parameters "
419        <<fractionOfSingleQuasiFreeNucleons<<" "<<fractionOfPairedQuasiFreeNucleons
420        <<" "<<clusteringCoefficient<<G4endl;
421  G4cout<<"G4Quasmon parameters "<<temperature<<" "<<halfTheStrangenessOfSee<<" "
422        <<etaToEtaPrime << G4endl;
423  G4cout<<"The Target PDG code = "<<targetPDGCode<<G4endl;
424  G4cout<<"The projectile momentum = "<<1./MeV*proj4Mom<<G4endl;
425  G4cout<<"The target momentum = "<<1./MeV*targ4Mom<<G4endl;
426#endif
427
428  // now call chips with this info in place
429  G4QHadronVector* output = 0;
430  if (particleCount!=0 && resA!=0)
431  {
432    //  G4QCHIPSWorld aWorld(nop);              // Create CHIPS World of nop particles
433    G4QCHIPSWorld::Get()->GetParticles(nop);
434    G4QEnvironment* pan= new G4QEnvironment(projHV, targetPDGCode);
435#ifdef pdebug
436      G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: CHIPS fragmentation, rA="
437            <<resA<<", #AbsPt="<<particleCount<<G4endl; 
438#endif
439    try
440    {
441      output = pan->Fragment();                 // The main fragmentation member function
442    }
443    catch(G4HadronicException& aR)
444    {
445      G4cerr << "Exception thrown of G4QStringChipsParticleLevelInterface "<<G4endl;
446      G4cerr << " targetPDGCode = "<< targetPDGCode <<G4endl;
447      G4cerr << " The projectile momentum = "<<1./MeV*proj4Mom<<G4endl;
448      G4cerr << " The target momentum = "<<1./MeV*targ4Mom<<G4endl<<G4endl;
449      G4cerr << " Dumping the information in the pojectile list"<<G4endl;
450      for(size_t i=0; i< projHV.size(); i++)
451      {
452        G4cerr <<"  Incoming 4-momentum and PDG code of "<<i<<"'th hadron: "
453        <<" "<< projHV[i]->Get4Momentum()<<" "<<projHV[i]->GetPDGCode()<<G4endl;
454      }
455      throw;
456    }
457    // clean up particles
458    std::for_each(projHV.begin(), projHV.end(), DeleteQHadron());
459    projHV.clear();
460    delete pan;
461  }
462  else
463  {
464#ifdef pdebug
465    G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: NO CHIPS fragmentation, rA="
466          <<resA<<", #AbsPt="<<particleCount<<G4endl; 
467#endif
468    output = new G4QHadronVector;
469  }   
470  // Fill the result.
471#ifdef CHIPSdebug
472  G4cout << "NEXT EVENT, EscapeExists="<<EscapeExists<<G4endl;
473#endif
474
475  // first decay and add all escaping particles.
476  if (EscapeExists) for (current = firstEscape; current!=theSorted.end(); current++)
477  {
478    G4KineticTrack* aResult = (*current).second;
479    G4ParticleDefinition* pdef=aResult->GetDefinition();
480    secondaries = NULL;
481    if(pdef->GetPDGWidth() > 0 && pdef->GetPDGLifeTime() < 5E-17*s )
482    {
483      secondaries = aResult->Decay();  // @@ Uses standard Decay, which is now wrong!
484    }
485    if ( secondaries == NULL )
486    {
487      theSec = new G4ReactionProduct(aResult->GetDefinition());
488      G4LorentzVector current4Mom = aResult->Get4Momentum();
489      current4Mom.boost(targ4Mom.boostVector());
490      theSec->SetTotalEnergy(current4Mom.t());
491      theSec->SetMomentum(current4Mom.vect());
492#ifdef pdebug
493      G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: *OUT* QGS stable PDG="
494            <<aResult->GetDefinition()->GetPDGEncoding()<<",4M="<<current4Mom<<G4endl; 
495#endif
496      theResult->push_back(theSec);
497    } 
498    else
499    {
500      for (unsigned int aSecondary=0; aSecondary<secondaries->size(); aSecondary++)
501      {
502        theSec=new G4ReactionProduct(secondaries->operator[](aSecondary)->GetDefinition());
503        G4LorentzVector current4Mom = secondaries->operator[](aSecondary)->Get4Momentum();
504        current4Mom.boost(targ4Mom.boostVector());
505        theSec->SetTotalEnergy(current4Mom.t());
506        theSec->SetMomentum(current4Mom.vect());
507#ifdef pdebug
508        G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: *OUT* QGS decay PDG="
509              <<secondaries->operator[](aSecondary)->GetDefinition()->GetPDGEncoding()
510              <<",4M="<<current4Mom<<G4endl; 
511#endif
512        theResult->push_back(theSec);
513      }
514      std::for_each(secondaries->begin(), secondaries->end(), DeleteKineticTrack());
515      delete secondaries;
516    }
517  }
518  std::for_each(theSecondaries->begin(), theSecondaries->end(), DeleteKineticTrack());
519  delete theSecondaries;
520   
521  // now add the quasmon output
522  G4int maxParticle=output->size();
523#ifdef CHIPSdebug
524  G4cout << "Number of particles from string"<<theResult->size()<<G4endl;
525  G4cout << "Number of particles from chips"<<maxParticle<<G4endl;
526#endif
527#ifdef pdebug
528  G4cout << "Number of particles from QGS="<<theResult->size()<<G4endl;
529  G4cout << "Number of particles from CHIPS="<<maxParticle<<G4endl;
530#endif
531  if(maxParticle) for(G4int particle = 0; particle < maxParticle; particle++)
532  {
533    if(output->operator[](particle)->GetNFragments() != 0) 
534    {
535      delete output->operator[](particle);
536      continue;
537    }
538    G4int pdgCode = output->operator[](particle)->GetPDGCode();
539
540
541#ifdef CHIPSdebug
542    G4cerr << "PDG code of chips particle = "<<pdgCode<<G4endl;
543#endif
544
545    G4ParticleDefinition * theDefinition;
546    // Note that I still have to take care of strange nuclei
547    // For this I need the mass calculation, and a changed interface
548    // for ion-table ==> work for Hisaya @@@@@@@
549    // Then I can sort out the pdgCode. I also need a decau process
550    // for strange nuclei; may be another chips interface
551    if(pdgCode>90000000) 
552    {
553      G4int aZ = (pdgCode-90000000)/1000;
554      if (aZ>1000) aZ=aZ%1000;  // patch for strange nuclei, to be repaired @@@@
555      G4int anN = pdgCode-90000000-1000*aZ;
556      if(anN>1000) anN=anN%1000; // patch for strange nuclei, to be repaired @@@@
557      if(pdgCode==91000000) theDefinition = G4Lambda::LambdaDefinition();
558      else if(pdgCode==92000000) theDefinition = G4Lambda::LambdaDefinition();
559      else if(pdgCode==93000000) theDefinition = G4Lambda::LambdaDefinition();
560      else if(pdgCode==94000000) theDefinition = G4Lambda::LambdaDefinition();
561      else if(pdgCode==95000000) theDefinition = G4Lambda::LambdaDefinition();
562      else if(pdgCode==96000000) theDefinition = G4Lambda::LambdaDefinition();
563      else if(pdgCode==97000000) theDefinition = G4Lambda::LambdaDefinition();
564      else if(pdgCode==98000000) theDefinition = G4Lambda::LambdaDefinition();
565      else if(aZ == 0 && anN == 1) theDefinition = G4Neutron::Neutron();
566      else theDefinition = G4ParticleTable::GetParticleTable()->FindIon(aZ,anN+aZ,0,aZ);
567    }   
568    else theDefinition = G4ParticleTable::GetParticleTable()->FindParticle(pdgCode);
569
570#ifdef CHIPSdebug
571    G4cout << "Particle code produced = "<< pdgCode <<G4endl;
572#endif
573
574    if(theDefinition)
575    {
576      theSec = new G4ReactionProduct(theDefinition);
577      G4LorentzVector current4Mom = output->operator[](particle)->Get4Momentum();
578      current4Mom.boost(targ4Mom.boostVector());
579      theSec->SetTotalEnergy(current4Mom.t());
580      theSec->SetMomentum(current4Mom.vect());
581#ifdef pdebug
582      G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: *OUT* CHIPS PDG="
583              <<theDefinition->GetPDGEncoding()<<",4M="<<current4Mom<<G4endl; 
584#endif
585      theResult->push_back(theSec);
586    }
587    else
588    {
589      G4cerr << G4endl<<"WARNING: "<<G4endl;
590      G4cerr << "Getting unknown pdgCode from chips in ParticleLevelInterface"<<G4endl;
591      G4cerr << "skipping particle with pdgCode = "<<pdgCode<<G4endl<<G4endl;
592    }
593   
594#ifdef CHIPSdebug
595    G4cout <<"CHIPS particles "<<theDefinition->GetPDGCharge()<<" "
596           << theDefinition->GetPDGEncoding()<<" "
597    << output->operator[](particle)->Get4Momentum() <<G4endl; 
598#endif
599
600    delete output->operator[](particle);
601  }
602  else
603  {
604    if(resA>0)
605    {
606      G4ParticleDefinition* theDefinition = G4Neutron::Neutron();
607      if(resA==1) // The residual nucleus at rest must be added to conserve BaryN & Charge
608      {
609        if(resZ == 1) theDefinition = G4Proton::Proton();
610      }
611      else theDefinition = G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ);
612      theSec = new G4ReactionProduct(theDefinition);
613      theSec->SetTotalEnergy(theDefinition->GetPDGMass());
614      theSec->SetMomentum(G4ThreeVector(0.,0.,0.));
615      theResult->push_back(theSec);
616      if(!resZ && resA>0) for(G4int ni=1; ni<resA; ni++) 
617      {
618        theSec = new G4ReactionProduct(theDefinition);
619        theSec->SetTotalEnergy(theDefinition->GetPDGMass());
620        theSec->SetMomentum(G4ThreeVector(0.,0.,0.));
621        theResult->push_back(theSec);
622      }
623    }
624  }
625  delete output;
626
627#ifdef CHIPSdebug
628  G4cout << "Number of particles"<<theResult->size()<<G4endl;
629  G4cout << G4endl;
630  G4cout << "QUASMON preparation info "
631         << 1./MeV*proj4Mom<<" "
632  << 1./MeV*targ4Mom<<" "
633  << nD<<" "<<nU<<" "<<nS<<" "<<nAD<<" "<<nAU<<" "<<nAS<<" "
634  << hitCount<<" "
635  << particleCount<<" "
636  << theLow<<" "
637  << theHigh<<" "
638  << G4endl;
639#endif
640
641  return theResult;
642} 
Note: See TracBrowser for help on using the repository browser.