source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPEnAngCorrelation.cc @ 1196

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

update processes

File size: 4.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//
26// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30#include "G4NeutronHPEnAngCorrelation.hh"
31#include "G4LorentzRotation.hh"
32#include "G4LorentzVector.hh"
33
34G4ReactionProduct * G4NeutronHPEnAngCorrelation::SampleOne(G4double anEnergy)
35{ 
36  G4ReactionProduct * result = new G4ReactionProduct;
37 
38  // do we have an appropriate distribution
39  if(nProducts!=1) throw G4HadronicException(__FILE__, __LINE__, "More than one product in SampleOne");
40 
41  // get the result
42  G4ReactionProductVector * temp=0;
43  G4int i=0;
44  while(temp == 0) temp = theProducts[i++].Sample(anEnergy);
45 
46  // is the multiplicity correct
47  if(temp->size()!=1) throw G4HadronicException(__FILE__, __LINE__, "SampleOne: Yield not correct");
48 
49  // fill result
50  result = temp->operator[](0);
51 
52  // some garbage collection
53  delete temp;
54 
55  // return result
56  return result;
57}
58
59G4ReactionProductVector * G4NeutronHPEnAngCorrelation::Sample(G4double anEnergy)
60{
61  G4ReactionProductVector * result = new G4ReactionProductVector;
62  G4int i;
63  G4ReactionProductVector * it;
64  G4ReactionProduct theCMS;
65  G4LorentzRotation toZ;
66  if(frameFlag==2)
67  {
68    // simplify and double check @
69    G4ThreeVector the3Neutron = theNeutron.GetMomentum();
70    G4double nEnergy = theNeutron.GetTotalEnergy();
71    G4ThreeVector the3Target = theTarget.GetMomentum();
72    G4double tEnergy = theTarget.GetTotalEnergy();
73    G4double totE = nEnergy+tEnergy;
74    G4ThreeVector the3CMS = the3Target+the3Neutron;
75    theCMS.SetMomentum(the3CMS);
76    G4double cmsMom = std::sqrt(the3CMS*the3CMS);
77    G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
78    theCMS.SetMass(sqrts);
79    theCMS.SetTotalEnergy(totE);
80    G4ReactionProduct aNeutron;
81    aNeutron.Lorentz(theNeutron, theCMS);
82    anEnergy = aNeutron.GetKineticEnergy();
83
84    G4LorentzVector Ptmp (aNeutron.GetMomentum(), aNeutron.GetTotalEnergy());
85    toZ.rotateZ(-1*Ptmp.phi());
86    toZ.rotateY(-1*Ptmp.theta());
87
88  }
89  theTotalMeanEnergy=0;
90  G4LorentzRotation toLab(toZ.inverse());
91  for(i=0; i<nProducts; i++)
92  {
93    it = theProducts[i].Sample(anEnergy);
94    G4double aMeanEnergy = theProducts[i].MeanEnergyOfThisInteraction();
95    if(aMeanEnergy>0)
96    {
97      theTotalMeanEnergy += aMeanEnergy;
98    }
99    else
100    {
101      theTotalMeanEnergy = anEnergy/nProducts+theProducts[i].GetQValue();
102    }
103    if(it!=0)
104    {
105      for(unsigned int ii=0; ii<it->size(); ii++)
106      {
107        G4LorentzVector pTmp1 (it->operator[](ii)->GetMomentum(),
108                               it->operator[](ii)->GetTotalEnergy());
109        pTmp1 = toLab*pTmp1;
110        it->operator[](ii)->SetMomentum(pTmp1.vect());
111        it->operator[](ii)->SetTotalEnergy(pTmp1.e());
112        if(frameFlag==1) // target rest
113        {
114          it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theTarget);
115        }
116        else if(frameFlag==2) // CMS
117        {
118#ifdef G4_NHP_DEBUG
119          cout <<"G4NeutronHPEnAngCorrelation: "<<
120                 it->at(ii)->GetTotalEnergy()<<" "<<
121                 it->at(ii)->GetMomentum()<<G4endl;
122#endif
123          it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
124        }
125        else
126        {
127          throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPEnAngCorrelation::Sample: The frame of the finalstate is not specified");
128        }
129        result->push_back(it->operator[](ii));
130      }
131    delete it;
132    }
133  }   
134  return result;
135}
136
Note: See TracBrowser for help on using the repository browser.