Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4PreCompoundFragmentVector.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // Hadronic Process: Nuclear Preequilibrium
 28 // by V. Lara 
 29 //
 30 // Modified:
 31 // 27.08.2010 V.Ivanchenko moved constructor and destructor to source, 
 32 //            simplify run time computations making inlined 
 33 //
 34 
 35 #include "G4PreCompoundFragmentVector.hh"
 36 
 37 G4PreCompoundFragmentVector::G4PreCompoundFragmentVector(pcfvector * avector) 
 38   : theChannels(nullptr), nChannels(0) 
 39 {
 40   SetVector(avector);
 41 }
 42 
 43 void G4PreCompoundFragmentVector::SetVector(pcfvector * avector)
 44 {
 45   if(avector != theChannels) {
 46     delete theChannels;
 47     theChannels = avector;
 48   }
 49   if(theChannels) {
 50     nChannels = (G4int)theChannels->size();
 51     probabilities.resize(nChannels, 0.0);
 52   } else {
 53     nChannels = 0;
 54     probabilities.clear();
 55   }
 56 }
 57 
 58 //for inverse cross section choice
 59 void G4PreCompoundFragmentVector::SetOPTxs(G4int opt)
 60 {    
 61   for (G4int i=0; i<nChannels; ++i) { 
 62     (*theChannels)[i]->SetOPTxs(opt);
 63   }
 64 }
 65 
 66 //for superimposed Coulomb Barrier for inverse  cross sections 
 67 void G4PreCompoundFragmentVector::UseSICB(G4bool use)
 68 {    
 69   for (G4int i=0; i< nChannels; ++i) { 
 70     (*theChannels)[i]->UseSICB(use);
 71   }
 72 }
 73 
 74 G4double G4PreCompoundFragmentVector::CalculateProbabilities(
 75          const G4Fragment & aFragment)
 76 {
 77   //G4cout << "## G4PreCompoundFragmentVector::CalculateProbabilities nCh= " 
 78   //   << nChannels << G4endl;
 79   G4double probtot = 0.0;  
 80   for (G4int i=0; i<nChannels; ++i) { 
 81     if ((*theChannels)[i]->Initialize(aFragment)) {
 82       G4double prob = (*theChannels)[i]->CalcEmissionProbability(aFragment);
 83       probtot += prob;
 84     }
 85     probabilities[i] = probtot;
 86     //G4cout<< "   probtot= " << probtot << " for "<< i << "-th channel" <<G4endl;
 87   }
 88   return probtot;
 89 }
 90 
 91 G4VPreCompoundFragment* G4PreCompoundFragmentVector::ChooseFragment()
 92 {
 93   //G4cout << "## G4PreCompoundFragmentVector::ChooseFragment nCh= " 
 94   //   << nChannels << G4endl;
 95   G4double x = probabilities[nChannels-1]*G4UniformRand();
 96   G4int i=0;
 97   for (; i<nChannels; ++i) { 
 98     if(x <= probabilities[i]) { break; }
 99   }
100   return (*theChannels)[i];  
101 }
102 
103