Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/lend/src/G4LENDFission.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 #include "G4LENDFission.hh"
 27 #include "G4SystemOfUnits.hh"
 28 #include "G4Nucleus.hh"
 29 #include "G4IonTable.hh"
 30   
 31 G4HadFinalState * G4LENDFission::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTarg )
 32 {
 33 
 34    G4double temp = aTrack.GetMaterial()->GetTemperature();
 35 
 36    //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this) 
 37    G4int iZ = aTarg.GetZ_asInt();
 38    G4int iA = aTarg.GetA_asInt();
 39    //G4int iM = aTarg.GetM_asInt();
 40    G4int iM = 0;
 41    if ( aTarg.GetIsotope() != NULL ) {
 42       iM = aTarg.GetIsotope()->Getm();
 43    }
 44 
 45    G4double ke = aTrack.GetKineticEnergy();
 46 
 47    G4HadFinalState* theResult = &theParticleChange;
 48    theResult->Clear();
 49 
 50    G4GIDI_target* aTarget = get_target_from_map( lend_manager->GetNucleusEncoding( iZ , iA , iM ) );
 51    if ( aTarget == NULL ) return returnUnchanged( aTrack , theResult );
 52    std::vector<G4GIDI_Product>* products = aTarget->getFissionFinalState( ke*MeV, temp, MyRNG, NULL );
 53    if ( products != NULL ) 
 54    {
 55       for ( G4int j = 0; j < int( products->size() ); j++ ) 
 56       {
 57          G4int jZ = (*products)[j].Z; 
 58          G4int jA = (*products)[j].A; 
 59          G4int jM = (*products)[j].m; 
 60 
 61          //G4cout << "Z = "    << (*products)[j].Z 
 62          //       << ", A = "  << (*products)[j].A 
 63          //       << ", EK = " << (*products)[j].kineticEnergy << " [MeV]" 
 64          //       << ", px = " << (*products)[j].px
 65          //       << ", py = " << (*products)[j].py
 66          //       << ", pz = " << (*products)[j].pz
 67          //       << ", birthTimeSec = " << (*products)[j].birthTimeSec << " [second]" 
 68          //       << G4endl;
 69 
 70          G4DynamicParticle* theSec = new G4DynamicParticle;
 71 
 72          if ( jZ > 0 )
 73          {
 74             theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ, jA , jM ) );
 75          } 
 76          else if ( jA == 1 && jZ == 0 )
 77          {
 78             theSec->SetDefinition( G4Neutron::Neutron() );
 79          } 
 80          else
 81          {
 82             theSec->SetDefinition( G4Gamma::Gamma() );
 83          } 
 84 
 85          theSec->SetMomentum( G4ThreeVector( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV ) );
 86          //G4cout << theSec->GetDefinition()->GetParticleName() << G4endl;
 87          theResult->AddSecondary( theSec, secID );
 88          //Set time for delayed neutrons
 89          //Current implementation is a little tricky, 
 90          if ( (*products)[j].birthTimeSec != 0 ) {
 91             G4double time = (*products)[j].birthTimeSec*second + aTrack.GetGlobalTime();
 92             theResult->GetSecondary(theResult->GetNumberOfSecondaries()-1)->SetTime(time);
 93          }
 94       } 
 95    }
 96    delete products;
 97 
 98    theResult->SetStatusChange( stopAndKill );
 99 
100    return theResult; 
101 
102 }
103 const std::pair<G4double, G4double> G4LENDFission::GetFatalEnergyCheckLevels() const
104 {
105         // max energy non-conservation is mass of heavy nucleus
106         //return std::pair<G4double, G4double>(5*perCent,250*GeV);
107         return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
108 }
109