Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4LowEGammaNuclearModel.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 // Physics model class G4LowEGammaNuclearModel 
 27 // Created:  15 May 2019
 28 // Author  V.Ivanchenko
 29 //  
 30 //
 31 
 32 #include "G4LowEGammaNuclearModel.hh"
 33 #include "G4SystemOfUnits.hh"
 34 #include "G4ParticleDefinition.hh"
 35 #include "G4Fragment.hh"
 36 #include "G4NucleiProperties.hh"
 37 #include "G4DynamicParticle.hh"
 38 #include "G4HadSecondary.hh"
 39 #include "G4ReactionProduct.hh"
 40 #include "G4ReactionProductVector.hh"
 41 #include "G4HadronicParameters.hh"
 42 #include "G4HadronicInteractionRegistry.hh"
 43 #include "G4PreCompoundModel.hh"
 44 #include "G4ThreeVector.hh"
 45 #include "G4PhysicsModelCatalog.hh"
 46 
 47 G4LowEGammaNuclearModel::G4LowEGammaNuclearModel() 
 48   : G4HadronicInteraction("GammaNPreco"),lab4mom(0.,0.,0.,0.), secID(-1)
 49 {
 50   secID = G4PhysicsModelCatalog::GetModelID("model_" + GetModelName());
 51   SetMinEnergy( 0.0*CLHEP::GeV );
 52   SetMaxEnergy( G4HadronicParameters::Instance()->GetMaxEnergy() );
 53 
 54   // reuse existing pre-compound model
 55   G4HadronicInteraction* p =
 56     G4HadronicInteractionRegistry::Instance()->FindModel("PRECO");
 57   fPreco = static_cast<G4PreCompoundModel*>(p);
 58   if(!fPreco) { fPreco = new G4PreCompoundModel(); }
 59 }
 60 
 61 G4LowEGammaNuclearModel::~G4LowEGammaNuclearModel()
 62 {}
 63 
 64 void G4LowEGammaNuclearModel::InitialiseModel()
 65 {}
 66 
 67 G4HadFinalState* G4LowEGammaNuclearModel::ApplyYourself(
 68      const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
 69 {
 70   theParticleChange.Clear();
 71 
 72   G4int A = theNucleus.GetA_asInt();
 73   G4int Z = theNucleus.GetZ_asInt();
 74 
 75   // Create initial state
 76   lab4mom.set(0.,0.,0.,G4NucleiProperties::GetNuclearMass(A, Z));
 77   lab4mom += aTrack.Get4Momentum();
 78 
 79   G4Fragment frag(A, Z, lab4mom);
 80 
 81   frag.SetCreatorModelID(secID);
 82 
 83   if (verboseLevel > 1) {
 84     G4cout << "G4LowEGammaNuclearModel::ApplyYourself initial G4Fragmet:" 
 85      << G4endl;
 86     G4cout << frag << G4endl;
 87   }
 88   G4ReactionProductVector* res = fPreco->DeExcite(frag);
 89 
 90   // secondaries produced
 91   if(res) {
 92 
 93     theParticleChange.SetStatusChange(stopAndKill);
 94     std::size_t nsec = res->size();
 95     if (verboseLevel > 1) {
 96       G4cout << "G4LowEGammaNuclearModel: " << nsec << " secondaries" << G4endl;
 97     }
 98     for(std::size_t i=0; i<nsec; ++i) {
 99       if((*res)[i]) {
100   G4double ekin = (*res)[i]->GetKineticEnergy();
101   G4ThreeVector dir(0.,0.,1.);
102   if(ekin > 0.0) { dir = (*res)[i]->GetMomentum().unit(); }
103   G4HadSecondary* news = new G4HadSecondary(
104           new G4DynamicParticle((*res)[i]->GetDefinition(), dir, ekin));
105   news->SetTime((*res)[i]->GetTOF());
106   news->SetCreatorModelID(secID);
107   theParticleChange.AddSecondary(*news);
108   if (verboseLevel > 1) {
109     G4cout << i << ". " << (*res)[i]->GetDefinition()->GetParticleName()
110      << " Ekin(MeV)= " << ekin/MeV
111      << " dir: " << dir << G4endl;
112   }
113   delete (*res)[i];
114         delete news;
115       }
116     } 
117     delete res;
118   } 
119   return &theParticleChange;
120 }
121 
122