source: trunk/source/processes/hadronic/models/de_excitation/evaporation/src/G4Evaporation.cc @ 843

Last change on this file since 843 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 10.5 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//
27// $Id: G4Evaporation.cc,v 1.7 2007/02/14 13:37:49 ahoward Exp $
28// GEANT4 tag $Name:  $
29//
30// Hadronic Process: Nuclear De-excitations
31// by V. Lara (Oct 1998)
32//
33// Alex Howard - added protection for negative probabilities in the sum, 14/2/07
34//
35#include "G4Evaporation.hh"
36#include "G4EvaporationFactory.hh"
37#include "G4EvaporationGEMFactory.hh"
38#include "G4HadronicException.hh"
39#include <numeric>
40
41G4Evaporation::G4Evaporation() 
42{
43  theChannelFactory = new G4EvaporationFactory();
44  theChannels = theChannelFactory->GetChannel();
45}
46
47G4Evaporation::G4Evaporation(const G4Evaporation &) : G4VEvaporation()
48{
49    throw G4HadronicException(__FILE__, __LINE__, "G4Evaporation::copy_constructor meant to not be accessable.");
50}
51
52
53G4Evaporation::~G4Evaporation()
54{
55  if (theChannels != 0) theChannels = 0;
56  if (theChannelFactory != 0) delete theChannelFactory;
57}
58
59const G4Evaporation & G4Evaporation::operator=(const G4Evaporation &)
60{
61    throw G4HadronicException(__FILE__, __LINE__, "G4Evaporation::operator= meant to not be accessable.");
62    return *this;
63}
64
65
66G4bool G4Evaporation::operator==(const G4Evaporation &) const
67{
68    return false;
69}
70
71G4bool G4Evaporation::operator!=(const G4Evaporation &) const
72{
73    return true;
74}
75
76
77void G4Evaporation::SetDefaultChannel()
78{
79  if (theChannelFactory != 0) delete theChannelFactory;
80  theChannelFactory = new G4EvaporationFactory();
81  theChannels = theChannelFactory->GetChannel();
82}
83
84void G4Evaporation::SetGEMChannel()
85{
86  if (theChannelFactory != 0) delete theChannelFactory;
87  theChannelFactory = new G4EvaporationGEMFactory();
88  theChannels = theChannelFactory->GetChannel();
89}
90
91
92G4FragmentVector * G4Evaporation::BreakItUp(const G4Fragment &theNucleus)
93{
94    G4FragmentVector * theResult = new G4FragmentVector;
95
96    // CHECK that Excitation Energy != 0
97    if (theNucleus.GetExcitationEnergy() <= 0.0) {
98        theResult->push_back(new G4Fragment(theNucleus));
99        return theResult;
100    }
101
102    // The residual nucleus (after evaporation of each fragment)
103    G4Fragment theResidualNucleus = theNucleus;
104
105    // Number of channels
106    G4int TotNumberOfChannels = theChannels->size(); 
107       
108
109    // Starts loop over evaporated particles
110    for (;;) 
111      {
112        // loop over evaporation channels
113        std::vector<G4VEvaporationChannel*>::iterator i;
114        for (i=theChannels->begin(); i != theChannels->end(); i++) 
115          {
116            (*i)->Initialize(theResidualNucleus);
117          }
118        // Can't use this form beacuse Initialize is a non const member function
119        //      for_each(theChannels->begin(),theChannels->end(),
120        //               bind2nd(mem_fun(&G4VEvaporationChannel::Initialize),theResidualNucleus));
121        // Work out total decay probability by summing over channels
122        G4double TotalProbability =  std::accumulate(theChannels->begin(),
123                                                       theChannels->end(),
124                                                       0.0,SumProbabilities());
125       
126        if (TotalProbability <= 0.0) 
127          {
128            // Will be no evaporation more
129            // write information about residual nucleus
130            theResult->push_back(new G4Fragment(theResidualNucleus));
131            break; 
132          } 
133        else 
134          {
135            // Selection of evaporation channel, fission or gamma
136            // G4double * EmissionProbChannel = new G4double(TotNumberOfChannels);
137            std::vector<G4double> EmissionProbChannel;
138            EmissionProbChannel.reserve(theChannels->size());
139           
140
141            // EmissionProbChannel[0] = theChannels->at(0)->GetEmissionProbability();
142
143
144            G4double first = theChannels->front()->GetEmissionProbability();
145
146            EmissionProbChannel.push_back(first >0 ? first : 0); // index 0
147
148
149            //      EmissionProbChannel.push_back(theChannels->front()->GetEmissionProbability()); // index 0
150           
151            for (i= (theChannels->begin()+1); i != theChannels->end(); i++) 
152              {
153                // EmissionProbChannel[i] = EmissionProbChannel[i-1] +
154                // theChannels->at(i)->GetEmissionProbability();
155                //              EmissionProbChannel.push_back(EmissionProbChannel.back() + (*i)->GetEmissionProbability());
156                first = (*i)->GetEmissionProbability();
157                EmissionProbChannel.push_back(first> 0? EmissionProbChannel.back() + first : EmissionProbChannel.back());
158              }
159
160            G4double shoot = G4UniformRand() * TotalProbability;
161            G4int j;
162            for (j=0; j < TotNumberOfChannels; j++) 
163              {
164                // if (shoot < EmissionProbChannel[i])
165                if (shoot < EmissionProbChannel[j]) 
166                  break;
167              }
168           
169            // delete [] EmissionProbChannel;
170            EmissionProbChannel.clear();
171                       
172            if( j >= TotNumberOfChannels ) 
173              {
174                throw G4HadronicException(__FILE__, __LINE__,  "G4Evaporation::BreakItUp: Can't define emission probability of the channels!" );
175              } 
176            else 
177              {
178                // Perform break-up
179                G4FragmentVector * theEvaporationResult = (*theChannels)[j]->BreakUp(theResidualNucleus);
180
181#ifdef debug
182                G4cout <<           "-----------------------------------------------------------\n"; 
183                G4cout << G4endl << " After the evaporation of a particle, testing conservation \n";
184                CheckConservation(theResidualNucleus,theEvaporationResult);
185                G4cout << G4endl
186                       <<           "------------------------------------------------------------\n"; 
187#endif 
188
189                // Check if chosen channel is fission (there are only two EXCITED fragments)
190                // or the channel could not evaporate anything
191                if ( theEvaporationResult->size() == 1 || 
192                     ((*(theEvaporationResult->begin()))->GetExcitationEnergy() > 0.0 && 
193                      (*(theEvaporationResult->end()-1))->GetExcitationEnergy() > 0.0) ) {
194                  // FISSION
195                  for (G4FragmentVector::iterator i = theEvaporationResult->begin();
196                       i != theEvaporationResult->end(); ++i) 
197                    {
198                      theResult->push_back(*(i));
199                    }
200                  delete theEvaporationResult;
201                  break;
202                } else {
203                  // EVAPORATION
204                  for (G4FragmentVector::iterator i = theEvaporationResult->begin();
205                       i != theEvaporationResult->end()-1; ++i) 
206                    {
207#ifdef PRECOMPOUND_TEST
208                      if ((*i)->GetA() == 0) (*i)->SetCreatorModel(G4String("G4PhotonEvaporation"));
209#endif
210                      theResult->push_back(*(i));
211                    }
212                  theResidualNucleus = *(theEvaporationResult->back());
213                  delete theEvaporationResult->back();
214                  delete theEvaporationResult;
215#ifdef PRECOMPOUND_TEST
216                  theResidualNucleus.SetCreatorModel(G4String("ResidualNucleus"));
217#endif
218                }
219              }
220          }
221      }
222   
223#ifdef debug
224    G4cout << "======== Evaporation Conservation Test ===========\n"
225           << "==================================================\n";
226    CheckConservation(theNucleus,theResult);
227    G4cout << "==================================================\n";
228#endif
229    return theResult;
230}
231
232
233
234#ifdef debug
235void G4Evaporation::CheckConservation(const G4Fragment & theInitialState,
236                                      G4FragmentVector * Result) const
237{
238    G4double ProductsEnergy =0;
239    G4ThreeVector ProductsMomentum;
240    G4int ProductsA = 0;
241    G4int ProductsZ = 0;
242    for (G4FragmentVector::iterator h = Result->begin(); h != Result->end(); h++) {
243        G4LorentzVector tmp = (*h)->GetMomentum();
244        ProductsEnergy += tmp.e();
245        ProductsMomentum += tmp.vect();
246        ProductsA += static_cast<G4int>((*h)->GetA());
247        ProductsZ += static_cast<G4int>((*h)->GetZ());
248    }
249
250    if (ProductsA != theInitialState.GetA()) {
251        G4cout << "!!!!!!!!!! Baryonic Number Conservation Violation !!!!!!!!!!" << G4endl;
252        G4cout << "G4Evaporation.cc: Barionic Number Conservation test for evaporation fragments" 
253               << G4endl; 
254        G4cout << "Initial A = " << theInitialState.GetA() 
255               << "   Fragments A = " << ProductsA << "   Diference --> " 
256               << theInitialState.GetA() - ProductsA << G4endl;
257    }
258    if (ProductsZ != theInitialState.GetZ()) {
259        G4cout << "!!!!!!!!!! Charge Conservation Violation !!!!!!!!!!" << G4endl;
260        G4cout << "G4Evaporation.cc: Charge Conservation test for evaporation fragments" 
261               << G4endl; 
262        G4cout << "Initial Z = " << theInitialState.GetZ() 
263               << "   Fragments Z = " << ProductsZ << "   Diference --> " 
264               << theInitialState.GetZ() - ProductsZ << G4endl;
265    }
266    if (std::abs(ProductsEnergy-theInitialState.GetMomentum().e()) > 1.0*keV) {
267        G4cout << "!!!!!!!!!! Energy Conservation Violation !!!!!!!!!!" << G4endl;
268        G4cout << "G4Evaporation.cc: Energy Conservation test for evaporation fragments" 
269               << G4endl; 
270        G4cout << "Initial E = " << theInitialState.GetMomentum().e()/MeV << " MeV"
271               << "   Fragments E = " << ProductsEnergy/MeV  << " MeV   Diference --> " 
272               << (theInitialState.GetMomentum().e() - ProductsEnergy)/MeV << " MeV" << G4endl;
273    } 
274    if (std::abs(ProductsMomentum.x()-theInitialState.GetMomentum().x()) > 1.0*keV || 
275        std::abs(ProductsMomentum.y()-theInitialState.GetMomentum().y()) > 1.0*keV ||
276        std::abs(ProductsMomentum.z()-theInitialState.GetMomentum().z()) > 1.0*keV) {
277        G4cout << "!!!!!!!!!! Momentum Conservation Violation !!!!!!!!!!" << G4endl;
278        G4cout << "G4Evaporation.cc: Momentum Conservation test for evaporation fragments" 
279               << G4endl; 
280        G4cout << "Initial P = " << theInitialState.GetMomentum().vect() << " MeV"
281               << "   Fragments P = " << ProductsMomentum  << " MeV   Diference --> " 
282               << theInitialState.GetMomentum().vect() - ProductsMomentum << " MeV" << G4endl;
283    }
284    return;
285}
286#endif
287
288
289
290
Note: See TracBrowser for help on using the repository browser.