source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPInelasticCompFS.cc @ 1340

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

update ti head

File size: 27.0 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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
31// 070606 bug fix and migrate to enable to Partial cases by T. Koi
32// 080603 bug fix for Hadron Hyper News #932 by T. Koi
33// 080612 bug fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4,6
34// 080717 bug fix of calculation of residual momentum by T. Koi
35// 080801 protect negative avalable energy by T. Koi
36//        introduce theNDLDataA,Z which has A and Z of NDL data by T. Koi
37// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
38// 090514 Fix bug in IC electron emission case
39//        Contribution from Chao Zhang (Chao.Zhang@usd.edu) and Dongming Mei(Dongming.Mei@usd.edu)
40// 100406 "nothingWasKnownOnHadron=1" then sample mu isotropic in CM
41//        add two_body_reaction
42// 100909 add safty
43//
44#include "G4NeutronHPInelasticCompFS.hh"
45#include "G4Nucleus.hh"
46#include "G4NucleiProperties.hh"
47#include "G4He3.hh"
48#include "G4Alpha.hh"
49#include "G4Electron.hh"
50#include "G4NeutronHPDataUsed.hh"
51#include "G4ParticleTable.hh"
52
53void G4NeutronHPInelasticCompFS::InitGammas(G4double AR, G4double ZR)
54{
55  //   char the[100] = {""};
56  //   std::ostrstream ost(the, 100, std::ios::out);
57  //   ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
58  //   G4String * aName = new G4String(the);
59  //   std::ifstream from(*aName, std::ios::in);
60
61   std::ostringstream ost;
62   ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
63   G4String aName = ost.str();
64   std::ifstream from(aName, std::ios::in);
65
66   if(!from) return; // no data found for this isotope
67   //   std::ifstream theGammaData(*aName, std::ios::in);
68   std::ifstream theGammaData(aName, std::ios::in);
69   
70   theGammas.Init(theGammaData);
71   //   delete aName;
72}
73
74void G4NeutronHPInelasticCompFS::Init (G4double A, G4double Z, G4String & dirName, G4String & aFSType)
75{
76
77  gammaPath = "/Inelastic/Gammas/";
78    if(!getenv("G4NEUTRONHPDATA")) 
79       throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
80  G4String tBase = getenv("G4NEUTRONHPDATA");
81  gammaPath = tBase+gammaPath;
82  G4String tString = dirName;
83  G4bool dbool;
84  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), tString, aFSType, dbool);
85  G4String filename = aFile.GetName();
86  theBaseA = aFile.GetA();
87  theBaseZ = aFile.GetZ();
88   theNDLDataA = (int)aFile.GetA();
89   theNDLDataZ = aFile.GetZ();
90  if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
91  {
92    if(getenv("NeutronHPNamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
93    hasAnyData = false;
94    hasFSData = false; 
95    hasXsec = false;
96    return;
97  }
98  theBaseA = A;
99  theBaseZ = G4int(Z+.5);
100  std::ifstream theData(filename, std::ios::in);
101  if(!theData)
102  {
103    hasAnyData = false;
104    hasFSData = false; 
105    hasXsec = false;
106    theData.close();
107    return;
108  }
109  // here we go
110  G4int infoType, dataType, dummy;
111  G4int sfType, it;
112  hasFSData = false; 
113  while (theData >> infoType)
114  {
115    hasFSData = true; 
116    theData >> dataType;
117    theData >> sfType >> dummy;
118    it = 50;
119    if(sfType>=600||(sfType<100&&sfType>=50)) it = sfType%50;
120    if(dataType==3) 
121    {
122      theData >> dummy >> dummy;
123      theXsection[it] = new G4NeutronHPVector;
124      G4int total;
125      theData >> total;
126      theXsection[it]->Init(theData, total, eV);
127      //std::cout << theXsection[it]->GetXsec(1*MeV) << std::endl;
128    }
129    else if(dataType==4)
130    {
131      theAngularDistribution[it] = new G4NeutronHPAngular;
132      theAngularDistribution[it]->Init(theData);
133    }
134    else if(dataType==5)
135    {
136      theEnergyDistribution[it] = new G4NeutronHPEnergyDistribution;
137      theEnergyDistribution[it]->Init(theData); 
138    }
139    else if(dataType==6)
140    {
141      theEnergyAngData[it] = new G4NeutronHPEnAngCorrelation;
142      theEnergyAngData[it]->Init(theData);
143    }
144    else if(dataType==12)
145    {
146      theFinalStatePhotons[it] = new G4NeutronHPPhotonDist;
147      theFinalStatePhotons[it]->InitMean(theData);
148    }
149    else if(dataType==13)
150    {
151      theFinalStatePhotons[it] = new G4NeutronHPPhotonDist;
152      theFinalStatePhotons[it]->InitPartials(theData);
153    }
154    else if(dataType==14)
155    {
156      theFinalStatePhotons[it]->InitAngular(theData);
157    }
158    else if(dataType==15)
159    {
160      theFinalStatePhotons[it]->InitEnergies(theData);
161    }
162    else
163    {
164      throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4NeutronHPInelasticCompFS");
165    }
166  }
167  theData.close();
168}
169
170G4int G4NeutronHPInelasticCompFS::SelectExitChannel(G4double eKinetic)
171{
172
173// it = 0 has without Photon
174  G4double running[50];
175  running[0] = 0;
176  unsigned int i;
177  for(i=0; i<50; i++)
178  {
179    if(i!=0) running[i]=running[i-1];
180    if(theXsection[i] != 0) 
181    {
182      running[i] += std::max(0., theXsection[i]->GetXsec(eKinetic));
183    }
184  }
185  G4double random = G4UniformRand();
186  G4double sum = running[49];
187  G4int it = 50;
188  if(0!=sum)
189  {
190    G4int i0;
191    for(i0=0; i0<50; i0++)
192    {
193      it = i0;
194      if(random < running[i0]/sum) break;
195    }
196  }
197//debug:  it = 1;
198  return it;
199}
200
201
202                                                                                                       //n,p,d,t,he3,a
203void G4NeutronHPInelasticCompFS::CompositeApply(const G4HadProjectile & theTrack, G4ParticleDefinition * aDefinition)
204{
205
206// prepare neutron
207    theResult.Clear();
208    G4double eKinetic = theTrack.GetKineticEnergy();
209    const G4HadProjectile *incidentParticle = &theTrack;
210    G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
211    theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
212    theNeutron.SetKineticEnergy( eKinetic );
213
214// prepare target
215    G4int i;
216    for(i=0; i<50; i++)
217    { if(theXsection[i] != 0) { break; } } 
218
219    G4double targetMass=0;
220    G4double eps = 0.0001;
221    targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
222                   G4Neutron::Neutron()->GetPDGMass();
223//    if(theEnergyAngData[i]!=0)
224//        targetMass = theEnergyAngData[i]->GetTargetMass();
225//    else if(theAngularDistribution[i]!=0)
226//        targetMass = theAngularDistribution[i]->GetTargetMass();
227//    else if(theFinalStatePhotons[50]!=0)
228//        targetMass = theFinalStatePhotons[50]->GetTargetMass();
229    G4Nucleus aNucleus;
230    G4ReactionProduct theTarget; 
231    G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
232    theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
233
234// prepare the residual mass
235    G4double residualMass=0;
236    G4double residualZ = theBaseZ - aDefinition->GetPDGCharge();
237    G4double residualA = theBaseA - aDefinition->GetBaryonNumber()+1;
238    residualMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(residualA+eps), static_cast<G4int>(residualZ+eps)) ) /
239                     G4Neutron::Neutron()->GetPDGMass();
240
241// prepare energy in target rest frame
242    G4ReactionProduct boosted;
243    boosted.Lorentz(theNeutron, theTarget);
244    eKinetic = boosted.GetKineticEnergy();
245//    G4double momentumInCMS = boosted.GetTotalMomentum();
246 
247// select exit channel for composite FS class.
248    G4int it = SelectExitChannel( eKinetic );
249   
250// set target and neutron in the relevant exit channel
251    InitDistributionInitialState(theNeutron, theTarget, it);   
252
253    G4ReactionProductVector * thePhotons = 0;
254    G4ReactionProductVector * theParticles = 0;
255    G4ReactionProduct aHadron;
256    aHadron.SetDefinition(aDefinition); // what if only cross-sections exist ==> Na 23 11 @@@@   
257    G4double availableEnergy = theNeutron.GetKineticEnergy() + theNeutron.GetMass() - aHadron.GetMass() +
258                             (targetMass - residualMass)*G4Neutron::Neutron()->GetPDGMass();
259//080730c
260    if ( availableEnergy < 0 )
261    {
262       //G4cout << "080730c Adjust availavleEnergy " << G4endl;
263       availableEnergy = 0; 
264    }
265    G4int nothingWasKnownOnHadron = 0;
266    G4int dummy;
267    G4double eGamm = 0;
268    G4int iLevel=it-1;
269
270//  TK without photon has it = 0
271    if( 50 == it ) 
272    {
273
274//    TK Excitation level is not determined
275      iLevel=-1;
276      aHadron.SetKineticEnergy(availableEnergy*residualMass*G4Neutron::Neutron()->GetPDGMass()/
277                               (aHadron.GetMass()+residualMass*G4Neutron::Neutron()->GetPDGMass()));
278
279      //aHadron.SetMomentum(theNeutron.GetMomentum()*(1./theNeutron.GetTotalMomentum())*
280      //                  std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
281      //                            aHadron.GetMass()*aHadron.GetMass()));
282
283      //TK add safty 100909
284      G4double p2 = ( aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy() - aHadron.GetMass()*aHadron.GetMass() );
285      G4double p = 0.0;
286      if ( p2 > 0.0 ) p = std::sqrt( p ); 
287
288      aHadron.SetMomentum(theNeutron.GetMomentum()*(1./theNeutron.GetTotalMomentum())*p );
289
290    }
291    else
292    {
293      while( iLevel!=-1 && theGammas.GetLevel(iLevel) == 0 ) { iLevel--; }
294    }
295
296
297    if ( theAngularDistribution[it] != 0 ) // MF4
298    {
299      if(theEnergyDistribution[it]!=0) // MF5
300      {
301        aHadron.SetKineticEnergy(theEnergyDistribution[it]->Sample(eKinetic, dummy));
302        G4double eSecN = aHadron.GetKineticEnergy();
303        eGamm = eKinetic-eSecN;
304        for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
305        {
306          if(theGammas.GetLevelEnergy(iLevel)<eGamm) break;
307        }
308        G4double random = 2*G4UniformRand();
309        iLevel+=G4int(random);
310        if(iLevel>theGammas.GetNumberOfLevels()-1)iLevel = theGammas.GetNumberOfLevels()-1;
311      }
312      else
313      {
314        G4double eExcitation = 0;
315        if(iLevel>=0) eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();   
316        while (eKinetic-eExcitation < 0 && iLevel>0)
317        {
318          iLevel--;
319          eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();   
320        }
321       
322        if(getenv("InelasticCompFSLogging") && eKinetic-eExcitation < 0) 
323        {
324          throw G4HadronicException(__FILE__, __LINE__, "SEVERE: InelasticCompFS: Consistency of data not good enough, please file report");
325        }
326        if(eKinetic-eExcitation < 0) eExcitation = 0;
327        if(iLevel!= -1) aHadron.SetKineticEnergy(eKinetic - eExcitation);
328       
329      }
330      theAngularDistribution[it]->SampleAndUpdate(aHadron);
331
332      if( theFinalStatePhotons[it] == 0 )
333      {
334// TK comment Most n,n* eneter to this 
335        thePhotons = theGammas.GetDecayGammas(iLevel);
336        eGamm -= theGammas.GetLevelEnergy(iLevel);
337        if(eGamm>0) // @ ok for now, but really needs an efficient way of correllated sampling @
338        {
339          G4ReactionProduct * theRestEnergy = new G4ReactionProduct;
340          theRestEnergy->SetDefinition(G4Gamma::Gamma());
341          theRestEnergy->SetKineticEnergy(eGamm);
342          G4double costh = 2.*G4UniformRand()-1.;
343          G4double phi = twopi*G4UniformRand();
344          theRestEnergy->SetMomentum(eGamm*std::sin(std::acos(costh))*std::cos(phi), 
345                                     eGamm*std::sin(std::acos(costh))*std::sin(phi),
346                                     eGamm*costh);
347          if(thePhotons == 0) { thePhotons = new G4ReactionProductVector; }
348          thePhotons->push_back(theRestEnergy);
349        }
350      }
351    }
352    else if(theEnergyAngData[it] != 0) // MF6 
353    {
354      theParticles = theEnergyAngData[it]->Sample(eKinetic);
355    }
356    else
357    {
358      // @@@ what to do, if we have photon data, but no info on the hadron itself
359      nothingWasKnownOnHadron = 1;
360    }
361
362    //G4cout << "theFinalStatePhotons it " << it << G4endl;
363    //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
364    //G4cout << "theFinalStatePhotons it " << it << G4endl;
365    //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
366    //G4cout << "thePhotons " << thePhotons << G4endl;
367
368    if ( theFinalStatePhotons[it] != 0 ) 
369    {
370       // the photon distributions are in the Nucleus rest frame.
371       // TK residual rest frame
372      G4ReactionProduct boosted;
373      boosted.Lorentz(theNeutron, theTarget);
374      G4double anEnergy = boosted.GetKineticEnergy();
375      thePhotons = theFinalStatePhotons[it]->GetPhotons(anEnergy);
376      G4double aBaseEnergy = theFinalStatePhotons[it]->GetLevelEnergy();
377      G4double testEnergy = 0;
378      if(thePhotons!=0 && thePhotons->size()!=0)
379      { aBaseEnergy-=thePhotons->operator[](0)->GetTotalEnergy(); }
380      if(theFinalStatePhotons[it]->NeedsCascade())
381      {
382        while(aBaseEnergy>0.01*keV)
383        {
384          // cascade down the levels
385          G4bool foundMatchingLevel = false;
386          G4int closest = 2;
387          G4double deltaEold = -1;
388          for(G4int i=1; i<it; i++)
389          {
390            if(theFinalStatePhotons[i]!=0) 
391            {
392              testEnergy = theFinalStatePhotons[i]->GetLevelEnergy();
393            }
394            else
395            {
396              testEnergy = 0;
397            }
398            G4double deltaE = std::abs(testEnergy-aBaseEnergy);
399            if(deltaE<0.1*keV)
400            {
401              G4ReactionProductVector * theNext = 
402                theFinalStatePhotons[i]->GetPhotons(anEnergy);
403              thePhotons->push_back(theNext->operator[](0));
404              aBaseEnergy = testEnergy-theNext->operator[](0)->GetTotalEnergy();
405              delete theNext;
406              foundMatchingLevel = true;
407              break; // ===>
408            }
409            if(theFinalStatePhotons[i]!=0 && ( deltaE<deltaEold||deltaEold<0.) )
410            {
411              closest = i;
412              deltaEold = deltaE;     
413            }
414          } // <=== the break goes here.
415          if(!foundMatchingLevel)
416          {
417            G4ReactionProductVector * theNext = 
418               theFinalStatePhotons[closest]->GetPhotons(anEnergy);
419            thePhotons->push_back(theNext->operator[](0));
420            aBaseEnergy = aBaseEnergy-theNext->operator[](0)->GetTotalEnergy();
421            delete theNext;
422          }
423        } 
424      }
425    }
426    unsigned int i0;
427    if(thePhotons!=0)
428    {
429      for(i0=0; i0<thePhotons->size(); i0++)
430      {
431        // back to lab
432        thePhotons->operator[](i0)->Lorentz(*(thePhotons->operator[](i0)), -1.*theTarget);
433      }
434    }
435    //G4cout << "nothingWasKnownOnHadron " << nothingWasKnownOnHadron << G4endl;
436    if(nothingWasKnownOnHadron)
437    {
438//    TKDB 100405
439//    In this case, hadron should be isotropic in CM
440//    mu and p should be correlated
441//
442      G4double totalPhotonEnergy = 0.0;
443      if ( thePhotons != 0 )
444      {
445         unsigned int nPhotons = thePhotons->size();
446         unsigned int i0;
447         for ( i0=0; i0<nPhotons; i0++)
448         {
449            //thePhotons has energies at LAB system
450            totalPhotonEnergy += thePhotons->operator[](i0)->GetTotalEnergy();
451         }
452      }
453
454      //isotropic distribution in CM
455      G4double mu = 1.0 - 2 * G4UniformRand();
456
457      // need momentums in target rest frame;
458      G4LorentzVector target_in_LAB ( theTarget.GetMomentum() , theTarget.GetTotalEnergy() );
459      G4ThreeVector boostToTargetRest = -target_in_LAB.boostVector();
460      G4LorentzVector proj_in_LAB = incidentParticle->Get4Momentum();
461
462      G4DynamicParticle* proj = new G4DynamicParticle( G4Neutron::Neutron() , proj_in_LAB.boost( boostToTargetRest ) ); 
463      G4DynamicParticle* targ = new G4DynamicParticle( G4ParticleTable::GetParticleTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , totalPhotonEnergy )  , G4ThreeVector(0) );
464      G4DynamicParticle* hadron = new G4DynamicParticle( aHadron.GetDefinition() , G4ThreeVector(0) );  // will be fill momentum
465
466      two_body_reaction ( proj , targ , hadron , mu );
467
468      G4LorentzVector hadron_in_trag_rest = hadron->Get4Momentum();
469      G4LorentzVector hadron_in_LAB = hadron_in_trag_rest.boost ( -boostToTargetRest );
470      aHadron.SetMomentum( hadron_in_LAB.v() );
471      aHadron.SetKineticEnergy ( hadron_in_LAB.e() - hadron_in_LAB.m() );
472
473      delete proj;
474      delete targ; 
475      delete hadron;
476
477//TKDB 100405
478/*
479      G4double totalPhotonEnergy = 0;
480      if(thePhotons!=0)
481      {
482        unsigned int nPhotons = thePhotons->size();
483        unsigned int i0;
484        for(i0=0; i0<nPhotons; i0++)
485        {
486          totalPhotonEnergy += thePhotons->operator[](i0)->GetTotalEnergy();
487        }
488      }
489      availableEnergy -= totalPhotonEnergy;
490      residualMass += totalPhotonEnergy/G4Neutron::Neutron()->GetPDGMass();
491      aHadron.SetKineticEnergy(availableEnergy*residualMass*G4Neutron::Neutron()->GetPDGMass()/
492                               (aHadron.GetMass()+residualMass*G4Neutron::Neutron()->GetPDGMass()));
493      G4double CosTheta = 1.0 - 2.0*G4UniformRand();
494      G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
495      G4double Phi = twopi*G4UniformRand();
496      G4ThreeVector Vector(std::cos(Phi)*SinTheta, std::sin(Phi)*SinTheta, CosTheta);
497      //aHadron.SetMomentum(Vector* std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
498      //                                 aHadron.GetMass()*aHadron.GetMass()));
499      G4double p2 = aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()- aHadron.GetMass()*aHadron.GetMass();
500
501      G4double p = 0.0;
502      if ( p2 > 0.0 )
503         p = std::sqrt ( p2 );
504
505      aHadron.SetMomentum( Vector*p );
506*/
507
508    }
509
510// fill the result
511// Beware - the recoil is not necessarily in the particles...
512// Can be calculated from momentum conservation?
513// The idea is that the particles ar emitted forst, and the gammas only once the
514// recoil is on the residual; assumption is that gammas do not contribute to
515// the recoil.
516// This needs more design @@@
517
518    G4int nSecondaries = 2; // the hadron and the recoil
519    G4bool needsSeparateRecoil = false;
520    G4int totalBaryonNumber = 0;
521    G4int totalCharge = 0;
522    G4ThreeVector totalMomentum(0);
523    if(theParticles != 0) 
524    {
525      nSecondaries = theParticles->size();
526      G4ParticleDefinition * aDef;
527      unsigned int i0;
528      for(i0=0; i0<theParticles->size(); i0++)
529      {
530        aDef = theParticles->operator[](i0)->GetDefinition();
531        totalBaryonNumber+=aDef->GetBaryonNumber();
532        totalCharge+=G4int(aDef->GetPDGCharge()+eps);
533        totalMomentum += theParticles->operator[](i0)->GetMomentum();
534      } 
535      if(totalBaryonNumber!=G4int(theBaseA+eps+incidentParticle->GetDefinition()->GetBaryonNumber())) 
536      {
537        needsSeparateRecoil = true;
538        nSecondaries++;
539        residualA = G4int(theBaseA+eps+incidentParticle->GetDefinition()->GetBaryonNumber()
540                          -totalBaryonNumber);
541        residualZ = G4int(theBaseZ+eps+incidentParticle->GetDefinition()->GetPDGCharge()
542                          -totalCharge);
543      }
544    }
545   
546    G4int nPhotons = 0;
547    if(thePhotons!=0) { nPhotons = thePhotons->size(); }
548    nSecondaries += nPhotons;
549       
550    G4DynamicParticle * theSec;
551   
552    if( theParticles==0 )
553    {
554      theSec = new G4DynamicParticle;   
555      theSec->SetDefinition(aHadron.GetDefinition());
556      theSec->SetMomentum(aHadron.GetMomentum());
557      theResult.AddSecondary(theSec);   
558 
559        aHadron.Lorentz(aHadron, theTarget);
560        G4ReactionProduct theResidual;   
561        theResidual.SetDefinition(G4ParticleTable::GetParticleTable()
562                                  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0)); 
563        theResidual.SetKineticEnergy(aHadron.GetKineticEnergy()*aHadron.GetMass()/theResidual.GetMass());
564
565        //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #6
566        //theResidual.SetMomentum(-1.*aHadron.GetMomentum());
567        G4ThreeVector incidentNeutronMomentum = theNeutron.GetMomentum();
568        theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
569
570        theResidual.Lorentz(theResidual, -1.*theTarget);
571        G4ThreeVector totalPhotonMomentum(0,0,0);
572        if(thePhotons!=0)
573        {
574          for(i=0; i<nPhotons; i++)
575          {
576            totalPhotonMomentum += thePhotons->operator[](i)->GetMomentum();
577          }
578        }
579        theSec = new G4DynamicParticle;   
580        theSec->SetDefinition(theResidual.GetDefinition());
581        theSec->SetMomentum(theResidual.GetMomentum()-totalPhotonMomentum);
582        theResult.AddSecondary(theSec);   
583    }
584    else
585    {
586      for(i0=0; i0<theParticles->size(); i0++)
587      {
588        theSec = new G4DynamicParticle; 
589        theSec->SetDefinition(theParticles->operator[](i0)->GetDefinition());
590        theSec->SetMomentum(theParticles->operator[](i0)->GetMomentum());
591        theResult.AddSecondary(theSec); 
592        delete theParticles->operator[](i0); 
593      } 
594      delete theParticles;
595      if(needsSeparateRecoil && residualZ!=0)
596      {
597        G4ReactionProduct theResidual;   
598        theResidual.SetDefinition(G4ParticleTable::GetParticleTable()
599                                  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0)); 
600        G4double resiualKineticEnergy  = theResidual.GetMass()*theResidual.GetMass();
601                 resiualKineticEnergy += totalMomentum*totalMomentum;
602                 resiualKineticEnergy  = std::sqrt(resiualKineticEnergy) - theResidual.GetMass();
603//        cout << "Kinetic energy of the residual = "<<resiualKineticEnergy<<endl;
604        theResidual.SetKineticEnergy(resiualKineticEnergy);
605
606        //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4
607        //theResidual.SetMomentum(-1.*totalMomentum);
608        //G4ThreeVector incidentNeutronMomentum = theNeutron.GetMomentum();
609        //theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
610//080717 TK Comment still do NOT include photon's mometum which produce by thePhotons
611        theResidual.SetMomentum( theNeutron.GetMomentum() + theTarget.GetMomentum() - totalMomentum );
612
613        theSec = new G4DynamicParticle;   
614        theSec->SetDefinition(theResidual.GetDefinition());
615        theSec->SetMomentum(theResidual.GetMomentum());
616        theResult.AddSecondary(theSec); 
617      } 
618    }
619    if(thePhotons!=0)
620    {
621      for(i=0; i<nPhotons; i++)
622      {
623        theSec = new G4DynamicParticle;   
624        //Bug reported Chao Zhang (Chao.Zhang@usd.edu), Dongming Mei(Dongming.Mei@usd.edu) Feb. 25, 2009
625        //theSec->SetDefinition(G4Gamma::Gamma());
626        theSec->SetDefinition( thePhotons->operator[](i)->GetDefinition() );
627        //But never cause real effect at least with G4NDL3.13 TK
628        theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
629        theResult.AddSecondary(theSec); 
630        delete thePhotons->operator[](i);
631      }
632// some garbage collection
633      delete thePhotons;
634    }
635
636//080721
637   G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , 0.0 );
638   G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
639   G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
640   G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
641   adjust_final_state ( init_4p_lab ); 
642
643// clean up the primary neutron
644    theResult.SetStatusChange(stopAndKill);
645}
646
647
648
649#include "G4RotationMatrix.hh"
650void G4NeutronHPInelasticCompFS::two_body_reaction ( G4DynamicParticle* proj, G4DynamicParticle* targ, G4DynamicParticle* hadron, G4double mu ) 
651{
652
653// Target rest flame
654// 4vector in targ rest frame;
655// targ could have excitation energy (photon energy will be emiited) tricky but,,,
656
657   G4LorentzVector before = proj->Get4Momentum() + targ->Get4Momentum();
658
659   G4ThreeVector p3_proj = proj->GetMomentum();
660   G4ThreeVector d = p3_proj.unit();
661   G4RotationMatrix rot; 
662   G4RotationMatrix rot1; 
663   rot1.setPhi( pi/2 + d.phi() );
664   G4RotationMatrix rot2; 
665   rot2.setTheta( d.theta() );
666   rot=rot2*rot1;
667   proj->SetMomentum( rot*p3_proj );
668
669// Now proj only has pz component;
670
671// mu in CM system
672
673   //Valid only for neutron incidence
674   G4DynamicParticle* residual = new G4DynamicParticle ( G4ParticleTable::GetParticleTable()->GetIon ( (G4int)( targ->GetDefinition()->GetPDGCharge() - hadron->GetDefinition()->GetPDGCharge() ) , (G4int)(targ->GetDefinition()->GetBaryonNumber() - hadron->GetDefinition()->GetBaryonNumber()+1) , 0 ) , G4ThreeVector(0) ); 
675
676   G4double Q = proj->GetDefinition()->GetPDGMass() + targ->GetDefinition()->GetPDGMass() 
677              - ( hadron->GetDefinition()->GetPDGMass() + residual->GetDefinition()->GetPDGMass() );
678
679   // Non Relativistic Case
680   G4double A = targ->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass();
681   G4double AA = hadron->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass(); 
682   G4double E1 = proj->GetKineticEnergy();
683   G4double beta = std::sqrt ( A*(A+1-AA)/AA*(1+(1+A)/A*Q/E1) );
684   G4double gamma = AA/(A+1-AA)*beta;
685   G4double E3 = AA/std::pow((1+A),2)*(beta*beta+1+2*beta*mu)*E1;
686   G4double omega3 = (1+beta*mu)/std::sqrt(beta*beta+1+2*beta*mu);
687
688   G4double E4 = (A+1-AA)/std::pow((1+A),2)*(gamma*gamma+1-2*gamma*mu)*E1;
689   G4double omega4 = (1-gamma*mu)/std::sqrt(gamma*gamma+1-2*gamma*mu);
690
691   hadron->SetKineticEnergy ( E3 );
692   
693   G4double M = hadron->GetDefinition()->GetPDGMass();
694   G4double pmag = std::sqrt ((E3+M)*(E3+M)-M*M) ;
695   G4ThreeVector p ( 0 , pmag*std::sqrt(1-omega3*omega3), pmag*omega3 );
696
697   G4double M4 = residual->GetDefinition()->GetPDGMass();
698   G4double pmag4 = std::sqrt ((E4+M4)*(E4+M4)-M4*M4) ;
699   G4ThreeVector p4 ( 0 , -pmag4*std::sqrt(1-omega4*omega4), pmag4*omega4 );
700
701// Rotate to orginal target rest flame.
702   p *= rot.inverse();
703   hadron->SetMomentum( p );
704// Now hadron had 4 momentum in target rest flame
705
706// TypeA
707   p4 *= rot.inverse();
708   residual->SetMomentum ( p4 );
709
710//TypeB1
711   //residual->Set4Momentum ( p4_residual );
712//TypeB2
713   //residual->SetMomentum ( p4_residual.v() );
714
715// Type A make difference in Momenutum
716// Type B1 make difference in Mass of residual
717// Type B2 make difference in total energy.
718
719   delete residual;
720
721}
Note: See TracBrowser for help on using the repository browser.