Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GammaParticipants.cc
Go to the documentation of this file.
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 "globals.hh"
27 #include "G4GammaParticipants.hh"
28 #include "G4LorentzVector.hh"
29 #include "G4V3DNucleus.hh"
30 #include <utility>
31 
32 // Class G4GammaParticipants
33 
34 // J.P. Wellisch, April 2002
35 // new participants class for gamma nuclear, with this design more can come with
36 // cross-section based, and quasi-eiconal model based modelling
37 //
38 // 20110805 M. Kelsey -- Follow change to G4V3DNucleus::GetNucleons()
39 
40 G4VSplitableHadron* G4GammaParticipants::SelectInteractions(const G4ReactionProduct &thePrimary)
41 {
42  // Check reaction threshold - goes to CheckThreshold
43  G4VSplitableHadron* aProjectile = new G4QGSMSplitableHadron(thePrimary, TRUE); // @@@ check the TRUE
44 
45  const std::vector<G4Nucleon>& theTargetNuc = theNucleus->GetNucleons();
46  G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
47  if((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
48  {
49  throw G4HadronicException(__FILE__, __LINE__,
50  "G4GammaParticipants::SelectInteractions: primary nan energy.");
51  }
52  G4double S = (aPrimaryMomentum + theTargetNuc[0].Get4Momentum()).mag2();
53  G4double ThresholdMass = thePrimary.GetMass() + theTargetNuc[0].GetDefinition()->GetPDGMass();
54  ModelMode = SOFT;
55  if (sqr(ThresholdMass + ThresholdParameter) > S)
56  {
58  //throw G4HadronicException(__FILE__, __LINE__,
59  // "Initial energy is too low. The 4-vectors of the input are inconsistant with the particle masses.");
60  }
61  if (sqr(ThresholdMass + QGSMThreshold) > S) // thus only diffractive in cascade!
62  {
64  }
65 
66  // first find the collisions HPW
67  std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
68  theInteractions.clear();
69  G4int totalCuts = 0;
70 
71  #ifdef debug_G4GammaParticipants
72  G4double eK = thePrimary.GetKineticEnergy()/GeV;
73  G4int nucleonCount = theTargetNuc.size(); // debug
74  #endif
75 
76  G4int theCurrent = static_cast<G4int> (theTargetNuc.size()*G4UniformRand());
77  const G4Nucleon& pNucleon = theTargetNuc[theCurrent];
78  G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(pNucleon);
79  theTargets.push_back(aTarget);
80  const_cast<G4Nucleon&>(pNucleon).Hit(aTarget);
81  if ( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) )
82  {
83  // diffractive interaction occurs
85  {
86  theSingleDiffExcitation.ExciteParticipants(aProjectile, aTarget);
87  } else {
88  theDiffExcitaton.ExciteParticipants(aProjectile, aTarget);
89  }
90  G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
91  aInteraction->SetTarget(aTarget);
92  theInteractions.push_back(aInteraction);
93  aInteraction->SetNumberOfDiffractiveCollisions(1);
94  totalCuts += 1;
95  } else {
96  // nondiffractive soft interaction occurs
97  aTarget->IncrementCollisionCount(1);
98  aProjectile->IncrementCollisionCount(1);
99  G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
100  aInteraction->SetTarget(aTarget);
101  aInteraction->SetNumberOfSoftCollisions(1);
102  theInteractions.push_back(aInteraction);
103  totalCuts += 1;
104  }
105  return aProjectile;
106 }
107 
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
double S(double temp)
const G4double QGSMThreshold
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
#define TRUE
Definition: globals.hh:55
void IncrementCollisionCount(G4int aCount)
G4QGSDiffractiveExcitation theDiffExcitaton
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
G4V3DNucleus * theNucleus
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void SetNumberOfDiffractiveCollisions(int)
std::vector< G4VSplitableHadron * > theTargets
static constexpr double GeV
Definition: G4SIunits.hh:217
G4ThreeVector GetMomentum() const
virtual const std::vector< G4Nucleon > & GetNucleons()=0
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
G4SingleDiffractiveExcitation theSingleDiffExcitation
std::vector< G4InteractionContent * > theInteractions
G4double GetMass() const
const G4double ThresholdParameter
void SetTarget(G4VSplitableHadron *aTarget)