[1350] | 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 | // $Id: G4CascadeRecoilMaker.cc,v 1.7 2010/09/25 06:44:30 mkelsey Exp $ |
---|
| 26 | // Geant4 tag: $Name: hadr-casc-V09-03-85 $ |
---|
| 27 | // |
---|
| 28 | // Collects generated cascade data (using Collider::collide() interface) |
---|
| 29 | // and computes the nuclear recoil kinematics needed to balance the event. |
---|
| 30 | // |
---|
| 31 | // 20100909 M. Kelsey -- Inspired by G4CascadeCheckBalance |
---|
| 32 | // 20100909 M. Kelsey -- Move G4IntraNucleiCascader::goodCase() here, add |
---|
| 33 | // tolerance for "almost zero" excitation energy |
---|
| 34 | // 20100910 M. Kelsey -- Drop getRecoilFragment() in favor of user calling |
---|
| 35 | // makeRecoilFragment() with returned non-const pointer. Drop |
---|
| 36 | // handling of excitons. |
---|
| 37 | // 20100921 M. Kelsey -- Return G4InuclNuclei using "makeRecoilNuclei()". |
---|
| 38 | // Repurpose "makeRecoilFragment()" to return G4Fragment. |
---|
| 39 | // 20100924 M. Kelsey -- Remove unusable G4Fragment::SetExcitationEnergy(). |
---|
| 40 | // Add deltaM to compute mass difference, use excitationEnergy |
---|
| 41 | // to force G4Fragment four-vector to match. |
---|
| 42 | |
---|
| 43 | #include "G4CascadeRecoilMaker.hh" |
---|
| 44 | #include "globals.hh" |
---|
| 45 | #include "G4CascadParticle.hh" |
---|
| 46 | #include "G4CascadeCheckBalance.hh" |
---|
| 47 | #include "G4CollisionOutput.hh" |
---|
| 48 | #include "G4Fragment.hh" |
---|
| 49 | #include "G4InuclElementaryParticle.hh" |
---|
| 50 | #include "G4InuclNuclei.hh" |
---|
| 51 | #include "G4InuclParticle.hh" |
---|
| 52 | #include "G4InuclSpecialFunctions.hh" |
---|
| 53 | #include "G4LorentzVector.hh" |
---|
| 54 | #include <vector> |
---|
| 55 | |
---|
| 56 | using namespace G4InuclSpecialFunctions; |
---|
| 57 | |
---|
| 58 | |
---|
| 59 | // Constructor and destructor |
---|
| 60 | |
---|
| 61 | G4CascadeRecoilMaker::G4CascadeRecoilMaker(G4double tolerance) |
---|
| 62 | : G4VCascadeCollider("G4CascadeRecoilMaker"), |
---|
| 63 | excTolerance(tolerance), inputEkin(0.), |
---|
| 64 | recoilA(0), recoilZ(0), excitationEnergy(0.) { |
---|
| 65 | balance = new G4CascadeCheckBalance(tolerance, tolerance, theName); |
---|
| 66 | } |
---|
| 67 | |
---|
| 68 | G4CascadeRecoilMaker::~G4CascadeRecoilMaker() { |
---|
| 69 | delete balance; |
---|
| 70 | } |
---|
| 71 | |
---|
| 72 | |
---|
| 73 | // Standard Collider interface (non-const output "buffer") |
---|
| 74 | |
---|
| 75 | void G4CascadeRecoilMaker::collide(G4InuclParticle* bullet, |
---|
| 76 | G4InuclParticle* target, |
---|
| 77 | G4CollisionOutput& output) { |
---|
| 78 | if (verboseLevel > 1) |
---|
| 79 | G4cout << " >>> G4CascadeRecoilMaker::collide" << G4endl; |
---|
| 80 | |
---|
| 81 | // Available energy needed for "goodNucleus()" test at end |
---|
| 82 | inputEkin = bullet ? bullet->getKineticEnergy() : 0.; |
---|
| 83 | |
---|
| 84 | balance->setVerboseLevel(verboseLevel); |
---|
| 85 | balance->collide(bullet, target, output); |
---|
| 86 | fillRecoil(); |
---|
| 87 | } |
---|
| 88 | |
---|
| 89 | // This is for use with G4IntraNucleiCascader |
---|
| 90 | |
---|
| 91 | void G4CascadeRecoilMaker::collide(G4InuclParticle* bullet, |
---|
| 92 | G4InuclParticle* target, |
---|
| 93 | const std::vector<G4InuclElementaryParticle>& particles, |
---|
| 94 | const std::vector<G4CascadParticle>& cparticles) { |
---|
| 95 | if (verboseLevel > 1) |
---|
| 96 | G4cout << " >>> G4CascadeRecoilMaker::collide(<EP>,<CP>)" << G4endl; |
---|
| 97 | |
---|
| 98 | // Available energy needed for "goodNucleus()" test at end |
---|
| 99 | inputEkin = bullet ? bullet->getKineticEnergy() : 0.; |
---|
| 100 | |
---|
| 101 | balance->setVerboseLevel(verboseLevel); |
---|
| 102 | balance->collide(bullet, target, particles, cparticles); |
---|
| 103 | fillRecoil(); |
---|
| 104 | } |
---|
| 105 | |
---|
| 106 | |
---|
| 107 | // Used current event configuration to construct recoil nucleus |
---|
| 108 | // NOTE: CheckBalance uses "final-initial", we want "initial-final" |
---|
| 109 | |
---|
| 110 | void G4CascadeRecoilMaker::fillRecoil() { |
---|
| 111 | recoilZ = -(balance->deltaQ()); // Charge "non-conservation" |
---|
| 112 | recoilA = -(balance->deltaB()); // Baryon "non-conservation" |
---|
| 113 | recoilMomentum = -(balance->deltaLV()); |
---|
| 114 | |
---|
| 115 | theExcitons.clear(); // Discard previous exciton configuraiton |
---|
| 116 | |
---|
| 117 | // Bertini uses MeV for excitation energy |
---|
| 118 | if (!goodFragment()) excitationEnergy = 0.; |
---|
| 119 | else excitationEnergy = deltaM() * GeV; |
---|
| 120 | |
---|
| 121 | // Allow for very small negative mass difference, and round to zero |
---|
| 122 | if (std::abs(excitationEnergy) < excTolerance) excitationEnergy = 0.; |
---|
| 123 | |
---|
| 124 | if (verboseLevel > 2) { |
---|
| 125 | G4cout << " recoil px " << recoilMomentum.px() |
---|
| 126 | << " py " << recoilMomentum.py() << " pz " << recoilMomentum.pz() |
---|
| 127 | << " E " << recoilMomentum.e() << " baryon " << recoilA |
---|
| 128 | << " charge " << recoilZ |
---|
| 129 | << "\n recoil mass " << recoilMomentum.m() |
---|
| 130 | << " 'excitation' energy " << excitationEnergy << G4endl; |
---|
| 131 | } |
---|
| 132 | } |
---|
| 133 | |
---|
| 134 | |
---|
| 135 | // Construct physical nucleus from recoil parameters, if reasonable |
---|
| 136 | |
---|
| 137 | G4InuclNuclei* G4CascadeRecoilMaker::makeRecoilNuclei(G4int model) { |
---|
| 138 | if (verboseLevel > 1) |
---|
| 139 | G4cout << " >>> G4CascadeRecoilMaker::makeRecoilNuclei" << G4endl; |
---|
| 140 | |
---|
| 141 | if (!goodRecoil()) { |
---|
| 142 | if (verboseLevel > 2 && !wholeEvent()) |
---|
| 143 | G4cout << theName << ": event recoil is not a physical nucleus" << G4endl; |
---|
| 144 | |
---|
| 145 | return 0; // Null pointer means no fragment |
---|
| 146 | } |
---|
| 147 | |
---|
| 148 | theRecoilNuclei.fill(recoilMomentum, recoilA, recoilZ, |
---|
| 149 | excitationEnergy, model); |
---|
| 150 | theRecoilNuclei.setExitonConfiguration(theExcitons); |
---|
| 151 | |
---|
| 152 | return &theRecoilNuclei; |
---|
| 153 | } |
---|
| 154 | |
---|
| 155 | |
---|
| 156 | // Construct pre-compound nuclear fragment from recoil parameters |
---|
| 157 | |
---|
| 158 | G4Fragment* G4CascadeRecoilMaker::makeRecoilFragment() { |
---|
| 159 | if (verboseLevel > 1) |
---|
| 160 | G4cout << " >>> G4CascadeRecoilMaker::makeRecoilFragment" << G4endl; |
---|
| 161 | |
---|
| 162 | if (!goodRecoil()) { |
---|
| 163 | if (verboseLevel > 2 && !wholeEvent()) |
---|
| 164 | G4cout << theName << ": event recoil is not a physical nucleus" << G4endl; |
---|
| 165 | |
---|
| 166 | return 0; // Null pointer means no fragment |
---|
| 167 | } |
---|
| 168 | |
---|
| 169 | theRecoilFragment.SetZandA_asInt(recoilZ, recoilA); // Note convention! |
---|
| 170 | |
---|
| 171 | // User may have overridden excitation energy; force four-momentum to match |
---|
| 172 | G4double fragMass = |
---|
| 173 | G4InuclNuclei::getNucleiMass(recoilA,recoilZ) + excitationEnergy/GeV; |
---|
| 174 | |
---|
| 175 | G4LorentzVector fragMom; fragMom.setVectM(recoilMomentum.vect(), fragMass); |
---|
| 176 | theRecoilFragment.SetMomentum(fragMom*GeV); // Bertini uses GeV! |
---|
| 177 | |
---|
| 178 | // Note: exciton configuration has to be set piece by piece |
---|
| 179 | theRecoilFragment.SetNumberOfHoles(theExcitons.protonHoles |
---|
| 180 | + theExcitons.neutronHoles); |
---|
| 181 | |
---|
| 182 | theRecoilFragment.SetNumberOfParticles(theExcitons.protonQuasiParticles |
---|
| 183 | + theExcitons.neutronQuasiParticles); |
---|
| 184 | |
---|
| 185 | theRecoilFragment.SetNumberOfCharged(theExcitons.protonQuasiParticles); |
---|
| 186 | |
---|
| 187 | return &theRecoilFragment; |
---|
| 188 | } |
---|
| 189 | |
---|
| 190 | |
---|
| 191 | // Compute raw mass difference from recoil parameters |
---|
| 192 | |
---|
| 193 | G4double G4CascadeRecoilMaker::deltaM() const { |
---|
| 194 | G4double nucMass = G4InuclNuclei::getNucleiMass(recoilA,recoilZ); |
---|
| 195 | return (recoilMomentum.m() - nucMass); |
---|
| 196 | } |
---|
| 197 | |
---|
| 198 | |
---|
| 199 | // Data quality checks |
---|
| 200 | |
---|
| 201 | G4bool G4CascadeRecoilMaker::goodFragment() const { |
---|
| 202 | return (recoilA>0 && recoilZ>=0 && recoilA >= recoilZ); |
---|
| 203 | } |
---|
| 204 | |
---|
| 205 | G4bool G4CascadeRecoilMaker::goodRecoil() const { |
---|
| 206 | return (goodFragment() && excitationEnergy > -excTolerance); |
---|
| 207 | } |
---|
| 208 | |
---|
| 209 | G4bool G4CascadeRecoilMaker::wholeEvent() const { |
---|
| 210 | return (recoilA==0 && recoilZ==0 && |
---|
| 211 | recoilMomentum.rho() < excTolerance/GeV && |
---|
| 212 | std::abs(recoilMomentum.e()) < excTolerance/GeV); |
---|
| 213 | } |
---|
| 214 | |
---|
| 215 | // Determine whether desired nuclear fragment is constructable outcome |
---|
| 216 | |
---|
| 217 | G4bool G4CascadeRecoilMaker::goodNucleus() const { |
---|
| 218 | if (verboseLevel > 2) { |
---|
| 219 | G4cout << " >>> G4CascadeRecoilMaker::goodNucleus" << G4endl; |
---|
| 220 | } |
---|
| 221 | |
---|
| 222 | const G4double minExcitation = 0.1*keV; |
---|
| 223 | const G4double reasonableExcitation = 7.0; // Multiple of binding energy |
---|
| 224 | const G4double fractionalExcitation = 0.2; // Fraction of input to excite |
---|
| 225 | |
---|
| 226 | if (!goodRecoil()) return false; // Not a sensible nucleus |
---|
| 227 | |
---|
| 228 | if (excitationEnergy < -excTolerance) return false; // Negative mass-diff |
---|
| 229 | |
---|
| 230 | if (excitationEnergy <= minExcitation) return true; // Effectively zero |
---|
| 231 | |
---|
| 232 | // Maximum possible excitation energy determined by initial energy |
---|
| 233 | G4double dm = bindingEnergy(recoilA,recoilZ); |
---|
| 234 | G4double exc_max0z = fractionalExcitation * inputEkin*GeV; |
---|
| 235 | G4double exc_dm = reasonableExcitation * dm; |
---|
| 236 | G4double exc_max = (exc_max0z > exc_dm) ? exc_max0z : exc_dm; |
---|
| 237 | |
---|
| 238 | if (verboseLevel > 3) { |
---|
| 239 | G4cout << " eexs " << excitationEnergy << " max " << exc_max |
---|
| 240 | << " dm " << dm << G4endl; |
---|
| 241 | } |
---|
| 242 | |
---|
| 243 | return (excitationEnergy < exc_max); // Below maximum possible |
---|
| 244 | } |
---|