source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4StringChipsInterface.cc @ 1315

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

update CVS release candidate geant4.9.3.01

File size: 17.3 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// !!! Was used in QBBC PL, NOW it is not. Must be absolete !!!
27// ============================================================
28
29//#define debug
30//#define pdebug
31
32#include "G4StringChipsInterface.hh"
33#include "globals.hh"
34#include <utility>
35#include <list>
36#include "G4KineticTrackVector.hh"
37#include "G4Nucleon.hh"
38#include "G4LorentzRotation.hh"
39
40G4StringChipsInterface::G4StringChipsInterface()
41{
42#ifdef CHIPSdebug
43  G4cout << "Please enter the energy loss per fermi in GeV"<<G4endl;
44  G4cin >> theEnergyLossPerFermi;
45#endif
46#ifdef debug
47  G4cout<<"G4StringChipsInterface::Constructor is called"<<G4endl;
48#endif
49  theEnergyLossPerFermi = 0.5*GeV;
50  // theEnergyLossPerFermi = 1.*GeV;
51}
52
53G4HadFinalState* G4StringChipsInterface::
54ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
55{
56#ifdef debug
57  G4cout<<"G4StringChipsInterface::ApplyYourself is called"<<G4endl;
58#endif
59  return theModel.ApplyYourself(aTrack, theNucleus);
60}
61
62G4ReactionProductVector* G4StringChipsInterface::
63Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus)
64{
65#ifdef debug
66  G4cout<<"G4StringChipsInterface::Propagate is called"<<G4endl;
67#endif
68  // Protection for non physical conditions
69 
70  if(theSecondaries->size() == 1) 
71    throw G4HadronicException(__FILE__, __LINE__, "G4StringChipsInterface: Only one particle from String models!");
72 
73  // Calculate the mean energy lost
74  std::pair<G4double, G4double> theImpact = theNucleus->RefetchImpactXandY();
75  G4double impactX = theImpact.first;
76  G4double impactY = theImpact.second;
77  G4double inpactPar2 = impactX*impactX + impactY*impactY;
78 
79  G4double radius2 = theNucleus->GetNuclearRadius(5*perCent);
80  radius2 *= radius2;
81  G4double pathlength = 0;
82  if(radius2 - inpactPar2>0) pathlength = 2.*std::sqrt(radius2 - inpactPar2);
83  G4double theEnergyLostInFragmentation = theEnergyLossPerFermi*pathlength/fermi;
84 
85  // now select all particles in range
86  std::list<std::pair<G4double, G4KineticTrack *> > theSorted;
87  std::list<std::pair<G4double, G4KineticTrack *> >::iterator current;
88  for(unsigned int secondary = 0; secondary<theSecondaries->size(); secondary++)
89  {
90    G4LorentzVector a4Mom = theSecondaries->operator[](secondary)->Get4Momentum();
91
92#ifdef CHIPSdebug
93    G4cout <<"ALL STRING particles "<<theSecondaries->operator[](secondary)->GetDefinition()->GetPDGCharge()<<" "
94           << theSecondaries->operator[](secondary)->GetDefinition()->GetPDGEncoding()<<" "
95           << a4Mom <<G4endl; 
96#endif
97
98    G4double toSort = a4Mom.rapidity();
99    std::pair<G4double, G4KineticTrack *> it;
100    it.first = toSort;
101    it.second = theSecondaries->operator[](secondary);
102    G4bool inserted = false;
103    for(current = theSorted.begin(); current!=theSorted.end(); current++)
104    {
105      if((*current).first > toSort)
106      {
107               theSorted.insert(current, it);
108               inserted = true;
109               break;
110      }
111    }
112    if(!inserted)
113    {
114      theSorted.push_front(it);
115    }
116  }
117 
118  G4LorentzVector proj4Mom(0.,0.,0.,0.);
119  // @@ Use the G4QContent class, which is exactly (nD,nU,nS,nAD,nAU,nAS) !
120  // The G4QContent class is a basic clas in CHIPS (not PDG Code as in GEANT4),
121  // so in CHIPS on can has a hadronic obgect (Quasmon) with any Quark Content.
122  // As a simple extantion for the hadron (which is a special case for Hadron)
123  // there is a clas G4QChipolino, which is a Quasmon, which can decay in two
124  // hadrons. In future the three-hadron class can be added...
125  G4int nD  = 0;
126  G4int nU  = 0;
127  G4int nS  = 0;
128  G4int nAD = 0;
129  G4int nAU = 0;
130  G4int nAS = 0;
131  std::list<std::pair<G4double, G4KineticTrack *> >::iterator firstEscaping = theSorted.begin();
132  G4double runningEnergy = 0;
133  G4int particleCount = 0;
134  G4LorentzVector theLow = (*(theSorted.begin())).second->Get4Momentum();
135  G4LorentzVector theHigh;
136
137#ifdef CHIPSdebug
138  G4cout << "CHIPS ENERGY LOST "<<theEnergyLostInFragmentation<<G4endl;
139  G4cout << "sorted rapidities event start"<<G4endl;
140#endif
141
142  G4ReactionProductVector * theResult = new G4ReactionProductVector;
143  G4ReactionProduct * theSec;
144  G4KineticTrackVector * secondaries;
145#ifdef pdebug
146  G4cout<<"G4StringChipsInterface::Propagate: Absorption"<<G4endl; 
147#endif
148 
149  // first decay and add all escaping particles.
150  for(current = theSorted.begin(); current!=theSorted.end(); current++)
151  {
152    firstEscaping = current;
153    if((*current).second->GetDefinition()->GetQuarkContent(3)!=0 ||
154       (*current).second->GetDefinition()->GetAntiQuarkContent(3) !=0)
155    {
156      G4KineticTrack * aResult = (*current).second;
157      G4ParticleDefinition * pdef=aResult->GetDefinition();
158      secondaries = NULL;
159      if ( pdef->GetPDGWidth() > 0 && pdef->GetPDGLifeTime() < 5E-17*s )
160      {
161        secondaries = aResult->Decay();
162      }
163      if ( secondaries == NULL )
164      {
165        theSec = new G4ReactionProduct(aResult->GetDefinition());
166        G4LorentzVector current4Mom = aResult->Get4Momentum();
167        theSec->SetTotalEnergy(current4Mom.t());
168        theSec->SetMomentum(current4Mom.vect());
169        theResult->push_back(theSec);
170      } 
171      else
172      {
173        for (unsigned int aSecondary=0; aSecondary<secondaries->size(); aSecondary++)
174        {
175          theSec = new G4ReactionProduct(secondaries->operator[](aSecondary)->GetDefinition());
176          G4LorentzVector current4Mom = secondaries->operator[](aSecondary)->Get4Momentum();
177          theSec->SetTotalEnergy(current4Mom.t());
178          theSec->SetMomentum(current4Mom.vect());
179          theResult->push_back(theSec);
180        }
181        std::for_each(secondaries->begin(), secondaries->end(), DeleteKineticTrack());
182        delete secondaries;
183      }
184      continue;
185    }
186    runningEnergy += (*current).second->Get4Momentum().t();
187    if(runningEnergy > theEnergyLostInFragmentation) break;
188   
189#ifdef CHIPSdebug
190    G4cout <<"ABSORBED STRING particles "<<current->second->GetDefinition()->GetPDGCharge()<<" "
191           << current->second->GetDefinition()->GetPDGEncoding()<<" "
192                  << current->second->Get4Momentum() <<G4endl; 
193#endif
194#ifdef pdebug
195    G4cout<<"G4StringChipsInterface::Propagate:C="
196          <<current->second->GetDefinition()->GetPDGCharge()<<", PDG="
197          << current->second->GetDefinition()->GetPDGEncoding()<<", 4M="
198                 << current->second->Get4Momentum() <<G4endl; 
199#endif
200
201    // projectile 4-momentum needed in constructor of quasmon
202    particleCount++;
203    theHigh = (*current).second->Get4Momentum(); 
204    proj4Mom += (*current).second->Get4Momentum(); 
205
206#ifdef CHIPSdebug
207    G4cout << "sorted rapidities "<<current->second->Get4Momentum().rapidity()<<G4endl;
208#endif
209   
210     // projectile quark contents needed for G4QContent construction (@@ ? -> G4QContent class)
211    nD  += (*current).second->GetDefinition()->GetQuarkContent(1);
212    nU  += (*current).second->GetDefinition()->GetQuarkContent(2);
213    nS  += (*current).second->GetDefinition()->GetQuarkContent(3);
214    nAD += (*current).second->GetDefinition()->GetAntiQuarkContent(1);
215    nAU += (*current).second->GetDefinition()->GetAntiQuarkContent(2);
216    nAS += (*current).second->GetDefinition()->GetAntiQuarkContent(3);
217  }
218  // construct G4QContent
219
220#ifdef CHIPSdebug
221  G4cout << "Quark content: d="<<nD<<", u="<<nU<< ", s="<< nS << G4endl;
222  G4cout << "Anti-quark content: anit-d="<<nAD<<", anti-u="<<nAU<< ", anti-s="<< nAS << G4endl;
223#endif
224
225  G4QContent theProjectiles(nD, nU, nS, nAD, nAU, nAS);
226
227#ifdef CHIPSdebug
228  G4cout << "G4QContent is constructed"<<endl;
229#endif
230 
231  // target properties needed in constructor of quasmon
232  // remove all hit nucleons to get Target code
233  theNucleus->StartLoop();
234  G4Nucleon * aNucleon;
235  G4int resA = 0;
236  G4int resZ = 0;
237  G4ThreeVector hitMomentum(0,0,0);
238  G4double hitMass = 0;
239  G4int hitCount = 0;
240  while((aNucleon = theNucleus->GetNextNucleon()))
241  {
242    if(!aNucleon->AreYouHit())
243    {
244      resA++;
245      resZ+=G4int (aNucleon->GetDefinition()->GetPDGCharge());
246    }
247    else
248    {
249      hitMomentum += aNucleon->GetMomentum().vect();
250      hitMass += aNucleon->GetMomentum().m();
251      hitCount ++;
252    }
253  }
254  G4int targetPDGCode = 90000000 + 1000*resZ + (resA-resZ);
255  G4double targetMass = theNucleus->GetMass();
256  targetMass -= hitMass;
257  G4double targetEnergy = std::sqrt(hitMomentum.mag2()+targetMass*targetMass);
258  // !! @@ Target should be at rest: hitMomentum=(0,0,0) @@ !! M.K.
259  G4LorentzVector targ4Mom(-1.*hitMomentum, targetEnergy);
260 
261  // construct the quasmon
262  //  G4int nop = 122; // clusters up to Alpha cluster
263  G4int nop = 152; // V.Ivanchenko set the same parameter as for all CHIPS models
264  G4double fractionOfSingleQuasiFreeNucleons = 0.5; // It is A-dependent (C=.85, U=.40)
265  G4double fractionOfPairedQuasiFreeNucleons = 0.05;
266  G4double clusteringCoefficient = 5.;
267  G4double temperature = 180.;
268  G4double halfTheStrangenessOfSee = 0.3; // = s/d = s/u
269  G4double etaToEtaPrime = 0.3;
270
271  G4QNucleus::SetParameters(fractionOfSingleQuasiFreeNucleons,
272                            fractionOfPairedQuasiFreeNucleons,
273                            clusteringCoefficient);
274  G4Quasmon::SetParameters(temperature,
275                           halfTheStrangenessOfSee,
276                           etaToEtaPrime);
277
278#ifdef CHIPSdebug
279  G4cout << "G4QNucleus parameters "<< fractionOfSingleQuasiFreeNucleons << " "
280         << fractionOfPairedQuasiFreeNucleons << " "<< clusteringCoefficient << G4endl;
281  G4cout << "G4Quasmon parameters "<< temperature << " "<< halfTheStrangenessOfSee << " "
282         <<etaToEtaPrime << G4endl;
283  G4cout << "The Target PDG code = "<<targetPDGCode<<G4endl;
284  G4cout << "The projectile momentum = "<<1./MeV*proj4Mom<<G4endl;
285  G4cout << "The target momentum = "<<1./MeV*targ4Mom<<G4endl;
286#endif
287
288 
289  // Chips expects all in target rest frame, along z.
290  // G4QCHIPSWorld aWorld(nop);              // Create CHIPS World of nop particles
291  G4QCHIPSWorld::Get()->GetParticles(nop);
292  G4QHadronVector projHV;
293  // target rest frame
294  proj4Mom.boost(-1.*targ4Mom.boostVector());
295  // now go along z
296  G4LorentzRotation toZ;
297  toZ.rotateZ(-1*proj4Mom.phi());
298  toZ.rotateY(-1*proj4Mom.theta());
299  proj4Mom = toZ*proj4Mom;
300  G4LorentzRotation toLab(toZ.inverse());
301
302#ifdef CHIPSdebug
303  G4cout << "a Boosted projectile vector along z"<<proj4Mom<<" "<<proj4Mom.mag()<<G4endl;
304#endif
305 
306  G4QHadron* iH = new G4QHadron(theProjectiles, 1./MeV*proj4Mom);
307  projHV.push_back(iH);
308
309  // now call chips with this info in place
310  G4QHadronVector * output = 0;
311  if (particleCount!=0)
312  {
313    G4QEnvironment* pan= new G4QEnvironment(projHV, targetPDGCode);
314    try
315    {
316      output = pan->Fragment();
317    }
318    catch(G4HadronicException & aR)
319    {
320      G4cerr << "Exception thrown passing through G4ChiralInvariantPhaseSpace "<<G4endl;
321      G4cerr << " targetPDGCode = "<< targetPDGCode <<G4endl;
322      G4cerr << " Dumping the information in the pojectile list"<<G4endl;
323      for(size_t i=0; i< projHV.size(); i++)
324      {
325               G4cerr <<"  Incoming 4-momentum and PDG code of "<<i<<"'th hadron: "
326               <<" "<< projHV[i]->Get4Momentum()<<" "<<projHV[i]->GetPDGCode()<<G4endl;
327      }
328      throw;
329    }
330    std::for_each(projHV.begin(), projHV.end(), DeleteQHadron());
331    projHV.clear();
332    delete pan;
333  }
334  else output = new G4QHadronVector;
335   
336  // Fill the result.
337#ifdef CHIPSdebug
338  G4cout << "NEXT EVENT"<<endl;
339#endif
340 
341  // first decay and add all escaping particles.
342  for(current = firstEscaping; current!=theSorted.end(); current++)
343  {
344    G4KineticTrack * aResult = (*current).second;
345    G4ParticleDefinition * pdef=aResult->GetDefinition();
346    secondaries = NULL;
347    if ( pdef->GetPDGWidth() > 0 && pdef->GetPDGLifeTime() < 5E-17*s )
348    {
349      secondaries = aResult->Decay();
350    }
351    if ( secondaries == NULL )
352    {
353      theSec = new G4ReactionProduct(aResult->GetDefinition());
354      G4LorentzVector current4Mom = aResult->Get4Momentum();
355      current4Mom = toLab*current4Mom;
356      current4Mom.boost(targ4Mom.boostVector());
357      theSec->SetTotalEnergy(current4Mom.t());
358      theSec->SetMomentum(current4Mom.vect());
359      theResult->push_back(theSec);
360    } 
361    else
362    {
363      for (unsigned int aSecondary=0; aSecondary<secondaries->size(); aSecondary++)
364      {
365        theSec = new G4ReactionProduct(secondaries->operator[](aSecondary)->GetDefinition());
366        G4LorentzVector current4Mom = secondaries->operator[](aSecondary)->Get4Momentum();
367        current4Mom = toLab*current4Mom;
368        current4Mom.boost(targ4Mom.boostVector());
369        theSec->SetTotalEnergy(current4Mom.t());
370        theSec->SetMomentum(current4Mom.vect());
371        theResult->push_back(theSec);
372      }
373      std::for_each(secondaries->begin(), secondaries->end(), DeleteKineticTrack());
374      delete secondaries;
375    }
376  }
377  std::for_each(theSecondaries->begin(), theSecondaries->end(), DeleteKineticTrack());
378  delete theSecondaries;
379   
380  // now add the quasmon output
381#ifdef CHIPSdebug
382  G4cout << "Number of particles from string"<<theResult->size()<<G4endl;
383  G4cout << "Number of particles from chips"<<output->size()<<G4endl;
384#endif
385
386  for(unsigned int particle = 0; particle < output->size(); particle++)
387  {
388    if(output->operator[](particle)->GetNFragments() != 0) 
389    {
390      delete output->operator[](particle);
391      continue;
392    }
393    theSec = new G4ReactionProduct; 
394    G4int pdgCode = output->operator[](particle)->GetPDGCode();
395    G4ParticleDefinition * theDefinition;
396    // Note that I still have to take care of strange nuclei
397    // For this I need the mass calculation, and a changed interface
398    // for ion-table ==> work for Hisaya @@@@@@@
399    // Then I can sort out the pdgCode. I also need a decau process
400    // for strange nuclei; may be another chips interface
401    if(pdgCode>90000000) 
402    {
403      G4int aZ = (pdgCode-90000000)/1000;
404      if (aZ>1000) aZ=aZ%1000;  // patch for strange nuclei, to be repaired @@@@
405      G4int anN = pdgCode-90000000-1000*aZ;
406      if(anN>1000) anN=anN%1000; // patch for strange nuclei, to be repaired @@@@
407      if(pdgCode==91000000) theDefinition = G4Lambda::LambdaDefinition();
408      else if(pdgCode==92000000) theDefinition = G4Lambda::LambdaDefinition();
409      else if(pdgCode==93000000) theDefinition = G4Lambda::LambdaDefinition();
410      else if(pdgCode==94000000) theDefinition = G4Lambda::LambdaDefinition();
411      else if(pdgCode==95000000) theDefinition = G4Lambda::LambdaDefinition();
412      else if(pdgCode==96000000) theDefinition = G4Lambda::LambdaDefinition();
413      else if(pdgCode==97000000) theDefinition = G4Lambda::LambdaDefinition();
414      else if(pdgCode==98000000) theDefinition = G4Lambda::LambdaDefinition();
415      else if(aZ == 0 && anN == 1) theDefinition = G4Neutron::Neutron();
416      else theDefinition = G4ParticleTable::GetParticleTable()->FindIon(aZ,anN+aZ,0,aZ);
417    }   
418    else
419    {
420      theDefinition = G4ParticleTable::GetParticleTable()->FindParticle(output->operator[](particle)->GetPDGCode());
421    }
422#ifdef CHIPSdebug
423    G4cout << "Particle code produced = "<< pdgCode <<G4endl;
424#endif
425
426    theSec = new G4ReactionProduct(theDefinition);
427    G4LorentzVector current4Mom = output->operator[](particle)->Get4Momentum();
428    current4Mom = toLab*current4Mom;
429    current4Mom.boost(targ4Mom.boostVector());
430    theSec->SetTotalEnergy(current4Mom.t());
431    theSec->SetMomentum(current4Mom.vect());
432    theResult->push_back(theSec);
433   
434#ifdef CHIPSdebug
435    G4cout <<"CHIPS particles "<<theDefinition->GetPDGCharge()<<" "
436           << theDefinition->GetPDGEncoding()<<" "
437           << current4Mom <<G4endl; 
438#endif
439
440    delete output->operator[](particle);
441  }
442  delete output;
443#ifdef CHIPSdebug
444  G4cout << "Number of particles"<<theResult->size()<<G4endl;
445  G4cout << G4endl;
446  // @@ G4QContent has even the out option!
447  G4cout << "QUASMON preparation info "
448         << 1./MeV*proj4Mom<<" "
449         << 1./MeV*targ4Mom<<" "
450         << nD<<" "<<nU<<" "<<nS<<" "<<nAD<<" "<<nAU<<" "<<nAS<<" "
451         << hitCount<<" "
452         << particleCount<<" "
453         << theLow<<" "
454         << theHigh<<" "
455         << G4endl;
456#endif
457
458  return theResult;
459} 
Note: See TracBrowser for help on using the repository browser.