| 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 | // $Id: G4GeneratorPrecompoundInterface.cc,v 1.10 2010/08/31 16:16:51 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: had-binary-V09-03-03 $
|
|---|
| 28 | //
|
|---|
| 29 | // -----------------------------------------------------------------------------
|
|---|
| 30 | // GEANT 4 class file
|
|---|
| 31 | //
|
|---|
| 32 | // History: first implementation
|
|---|
| 33 | // HPW, 10DEC 98, the decay part originally written by Gunter Folger
|
|---|
| 34 | // in his FTF-test-program.
|
|---|
| 35 | //
|
|---|
| 36 | //
|
|---|
| 37 | // -----------------------------------------------------------------------------
|
|---|
| 38 |
|
|---|
| 39 | #include "G4GeneratorPrecompoundInterface.hh"
|
|---|
| 40 | #include "G4DynamicParticleVector.hh"
|
|---|
| 41 | #include "G4KineticTrackVector.hh"
|
|---|
| 42 | #include "G4Proton.hh"
|
|---|
| 43 | #include "G4Neutron.hh"
|
|---|
| 44 | #include "G4V3DNucleus.hh"
|
|---|
| 45 | #include "G4Nucleon.hh"
|
|---|
| 46 | #include "G4FragmentVector.hh"
|
|---|
| 47 | #include "G4ReactionProduct.hh"
|
|---|
| 48 | #include "G4PreCompoundModel.hh"
|
|---|
| 49 | #include "G4ExcitationHandler.hh"
|
|---|
| 50 |
|
|---|
| 51 | G4GeneratorPrecompoundInterface::G4GeneratorPrecompoundInterface(G4VPreCompoundModel* p)
|
|---|
| 52 | : CaptureThreshold(80*MeV)
|
|---|
| 53 | {
|
|---|
| 54 | proton = G4Proton::Proton();
|
|---|
| 55 | neutron = G4Neutron::Neutron();
|
|---|
| 56 | if(p) { SetDeExcitation(p); }
|
|---|
| 57 | else { SetDeExcitation(new G4PreCompoundModel(new G4ExcitationHandler())); }
|
|---|
| 58 | }
|
|---|
| 59 |
|
|---|
| 60 | G4GeneratorPrecompoundInterface::~G4GeneratorPrecompoundInterface()
|
|---|
| 61 | {}
|
|---|
| 62 |
|
|---|
| 63 | G4ReactionProductVector* G4GeneratorPrecompoundInterface::
|
|---|
| 64 | Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus)
|
|---|
| 65 | {
|
|---|
| 66 | G4ReactionProductVector * theTotalResult = new G4ReactionProductVector;
|
|---|
| 67 |
|
|---|
| 68 | // decay the strong resonances
|
|---|
| 69 | G4KineticTrackVector *result1, *secondaries, *result;
|
|---|
| 70 | result1=theSecondaries;
|
|---|
| 71 | result=new G4KineticTrackVector();
|
|---|
| 72 | //G4cout << "### G4GeneratorPrecompoundInterface::Propagate "
|
|---|
| 73 | // << result1->size() << " tracks " << theDeExcitation << G4endl;
|
|---|
| 74 | for (unsigned int aResult=0; aResult < result1->size(); ++aResult)
|
|---|
| 75 | {
|
|---|
| 76 | G4ParticleDefinition * pdef;
|
|---|
| 77 | pdef=result1->operator[](aResult)->GetDefinition();
|
|---|
| 78 | secondaries=0;
|
|---|
| 79 | if ( pdef->IsShortLived() )
|
|---|
| 80 | {
|
|---|
| 81 | secondaries = result1->operator[](aResult)->Decay();
|
|---|
| 82 | }
|
|---|
| 83 | if ( 0 == secondaries )
|
|---|
| 84 | {
|
|---|
| 85 | result->push_back(result1->operator[](aResult));
|
|---|
| 86 | result1->operator[](aResult)=NULL; //protect for clearAndDestroy
|
|---|
| 87 | }
|
|---|
| 88 | else
|
|---|
| 89 | {
|
|---|
| 90 | unsigned int amax = secondaries->size();
|
|---|
| 91 | for (unsigned int aSecondary=0; aSecondary<amax; ++aSecondary)
|
|---|
| 92 | {
|
|---|
| 93 | result1->push_back(secondaries->operator[](aSecondary));
|
|---|
| 94 | }
|
|---|
| 95 | delete secondaries;
|
|---|
| 96 | }
|
|---|
| 97 | }
|
|---|
| 98 | //G4cout << "Delete tracks" << G4endl;
|
|---|
| 99 | std::for_each(result1->begin(), result1->end(), DeleteKineticTrack());
|
|---|
| 100 | delete result1;
|
|---|
| 101 |
|
|---|
| 102 | // prepare the fragment
|
|---|
| 103 | G4int anA=theNucleus->GetMassNumber();
|
|---|
| 104 | G4int aZ=theNucleus->GetCharge();
|
|---|
| 105 | G4int numberOfEx = 0;
|
|---|
| 106 | G4int numberOfCh = 0;
|
|---|
| 107 | G4int numberOfHoles = 0;
|
|---|
| 108 | G4double exEnergy = 0.0;
|
|---|
| 109 | G4double R = theNucleus->GetNuclearRadius();
|
|---|
| 110 | G4ThreeVector exciton3Momentum(0.,0.,0.);
|
|---|
| 111 |
|
|---|
| 112 | // loop over secondaries
|
|---|
| 113 | unsigned int amax = result->size();
|
|---|
| 114 | for(unsigned int list=0; list<amax; ++list)
|
|---|
| 115 | {
|
|---|
| 116 | G4KineticTrack *aTrack = result->operator[](list);
|
|---|
| 117 | G4ParticleDefinition* part = aTrack->GetDefinition();
|
|---|
| 118 | G4double e = aTrack->Get4Momentum().e();
|
|---|
| 119 | G4double mass = aTrack->Get4Momentum().mag();
|
|---|
| 120 | G4ThreeVector mom = aTrack->Get4Momentum().vect();
|
|---|
| 121 | if((part != proton && part != neutron) ||
|
|---|
| 122 | (e > mass + CaptureThreshold) ||
|
|---|
| 123 | (aTrack->GetPosition().mag() > R))
|
|---|
| 124 | {
|
|---|
| 125 | G4ReactionProduct * theNew = new G4ReactionProduct(part);
|
|---|
| 126 | theNew->SetMomentum(mom);
|
|---|
| 127 | theNew->SetTotalEnergy(e);
|
|---|
| 128 | theTotalResult->push_back(theNew);
|
|---|
| 129 | }
|
|---|
| 130 | else
|
|---|
| 131 | {
|
|---|
| 132 | // within the nucleus, neutron or proton
|
|---|
| 133 | // now calculate A, Z of the fragment, momentum, number of exciton states
|
|---|
| 134 | ++anA;
|
|---|
| 135 | ++numberOfEx;
|
|---|
| 136 | G4int Z = G4int(part->GetPDGCharge()/eplus + 0.1);
|
|---|
| 137 | aZ += Z;
|
|---|
| 138 | numberOfCh += Z;
|
|---|
| 139 | exciton3Momentum += mom;
|
|---|
| 140 | exEnergy += (e - mass);
|
|---|
| 141 | }
|
|---|
| 142 | }
|
|---|
| 143 |
|
|---|
| 144 | // loop over wounded nucleus
|
|---|
| 145 | G4Nucleon * theCurrentNucleon =
|
|---|
| 146 | theNucleus->StartLoop() ? theNucleus->GetNextNucleon() : 0;
|
|---|
| 147 | while(0 != theCurrentNucleon)
|
|---|
| 148 | {
|
|---|
| 149 | if(theCurrentNucleon->AreYouHit())
|
|---|
| 150 | {
|
|---|
| 151 | ++numberOfHoles;
|
|---|
| 152 | ++numberOfEx;
|
|---|
| 153 | --anA;
|
|---|
| 154 | aZ -= G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/eplus + 0.1);
|
|---|
| 155 | exciton3Momentum -= theCurrentNucleon->Get4Momentum().vect();
|
|---|
| 156 | exEnergy += theCurrentNucleon->GetBindingEnergy();
|
|---|
| 157 | }
|
|---|
| 158 | theCurrentNucleon = theNucleus->GetNextNucleon();
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | if(0!=anA && 0!=aZ)
|
|---|
| 162 | {
|
|---|
| 163 | G4double fMass = G4NucleiProperties::GetNuclearMass(anA, aZ);
|
|---|
| 164 | fMass += exEnergy;
|
|---|
| 165 |
|
|---|
| 166 | G4LorentzVector exciton4Momentum(exciton3Momentum,
|
|---|
| 167 | std::sqrt(exciton3Momentum.mag2() + fMass*fMass));
|
|---|
| 168 |
|
|---|
| 169 | G4Fragment anInitialState(anA, aZ, exciton4Momentum);
|
|---|
| 170 | anInitialState.SetNumberOfParticles(numberOfEx-numberOfHoles);
|
|---|
| 171 | anInitialState.SetNumberOfCharged(numberOfCh);
|
|---|
| 172 | anInitialState.SetNumberOfHoles(numberOfHoles);
|
|---|
| 173 | G4ReactionProductVector * aPreResult = theDeExcitation->DeExcite(anInitialState);
|
|---|
| 174 |
|
|---|
| 175 | // fill pre-compound part into the result, and return
|
|---|
| 176 | unsigned int amax = aPreResult->size();
|
|---|
| 177 | for(unsigned int ll=0; ll<amax; ++ll)
|
|---|
| 178 | {
|
|---|
| 179 | theTotalResult->push_back(aPreResult->operator[](ll));
|
|---|
| 180 | }
|
|---|
| 181 | delete aPreResult;
|
|---|
| 182 | }
|
|---|
| 183 |
|
|---|
| 184 | std::for_each(result->begin(), result->end(), DeleteKineticTrack());
|
|---|
| 185 | delete result;
|
|---|
| 186 | return theTotalResult;
|
|---|
| 187 | }
|
|---|
| 188 |
|
|---|
| 189 | G4HadFinalState* G4GeneratorPrecompoundInterface::
|
|---|
| 190 | ApplyYourself(const G4HadProjectile &, G4Nucleus & )
|
|---|
| 191 | {
|
|---|
| 192 | G4cout << "G4GeneratorPrecompoundInterface: ApplyYourself interface called stand-allone."
|
|---|
| 193 | << G4endl;
|
|---|
| 194 | G4cout << "This class is only a mediator between generator and precompound"<<G4endl;
|
|---|
| 195 | G4cout << "Please remove from your physics list."<<G4endl;
|
|---|
| 196 | throw G4HadronicException(__FILE__, __LINE__, "SEVERE: G4GeneratorPrecompoundInterface model interface called stand-allone.");
|
|---|
| 197 | return new G4HadFinalState;
|
|---|
| 198 | }
|
|---|