[819] | 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 | // neutron_hp -- source file |
---|
| 27 | // J.P. Wellisch, Nov-1996 |
---|
| 28 | // A prototype of the low energy neutron transport model. |
---|
| 29 | // |
---|
| 30 | #include "G4NeutronHPFSFissionFS.hh" |
---|
| 31 | #include "G4ReactionProduct.hh" |
---|
| 32 | #include "G4Nucleus.hh" |
---|
| 33 | #include "G4Proton.hh" |
---|
| 34 | #include "G4Deuteron.hh" |
---|
| 35 | #include "G4Triton.hh" |
---|
| 36 | #include "G4Alpha.hh" |
---|
| 37 | #include "G4ThreeVector.hh" |
---|
| 38 | #include "G4Poisson.hh" |
---|
| 39 | #include "G4LorentzVector.hh" |
---|
| 40 | #include "G4NeutronHPDataUsed.hh" |
---|
| 41 | |
---|
| 42 | void G4NeutronHPFSFissionFS::Init (G4double A, G4double Z, G4String & dirName, G4String & ) |
---|
| 43 | { |
---|
| 44 | G4String tString = "/FS/"; |
---|
| 45 | G4bool dbool; |
---|
| 46 | G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), dirName, tString, dbool); |
---|
| 47 | G4String filename = aFile.GetName(); |
---|
| 48 | if(!dbool) |
---|
| 49 | { |
---|
| 50 | hasAnyData = false; |
---|
| 51 | hasFSData = false; |
---|
| 52 | hasXsec = false; |
---|
| 53 | return; |
---|
| 54 | } |
---|
| 55 | std::ifstream theData(filename, std::ios::in); |
---|
| 56 | |
---|
| 57 | // here it comes |
---|
| 58 | G4int infoType, dataType; |
---|
| 59 | hasFSData = false; |
---|
| 60 | while (theData >> infoType) |
---|
| 61 | { |
---|
| 62 | hasFSData = true; |
---|
| 63 | theData >> dataType; |
---|
| 64 | switch(infoType) |
---|
| 65 | { |
---|
| 66 | case 1: |
---|
| 67 | if(dataType==4) theNeutronAngularDis.Init(theData); |
---|
| 68 | if(dataType==5) thePromptNeutronEnDis.Init(theData); |
---|
| 69 | if(dataType==12) theFinalStatePhotons.InitMean(theData); |
---|
| 70 | if(dataType==14) theFinalStatePhotons.InitAngular(theData); |
---|
| 71 | if(dataType==15) theFinalStatePhotons.InitEnergies(theData); |
---|
| 72 | break; |
---|
| 73 | case 2: |
---|
| 74 | if(dataType==1) theFinalStateNeutrons.InitMean(theData); |
---|
| 75 | break; |
---|
| 76 | case 3: |
---|
| 77 | if(dataType==1) theFinalStateNeutrons.InitDelayed(theData); |
---|
| 78 | if(dataType==5) theDelayedNeutronEnDis.Init(theData); |
---|
| 79 | break; |
---|
| 80 | case 4: |
---|
| 81 | if(dataType==1) theFinalStateNeutrons.InitPrompt(theData); |
---|
| 82 | break; |
---|
| 83 | case 5: |
---|
| 84 | if(dataType==1) theEnergyRelease.Init(theData); |
---|
| 85 | break; |
---|
| 86 | default: |
---|
| 87 | G4cout << "G4NeutronHPFSFissionFS::Init: unknown data type"<<dataType<<G4endl; |
---|
| 88 | throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPFSFissionFS::Init: unknown data type"); |
---|
| 89 | break; |
---|
| 90 | } |
---|
| 91 | } |
---|
| 92 | targetMass = theFinalStateNeutrons.GetTargetMass(); |
---|
| 93 | theData.close(); |
---|
| 94 | } |
---|
| 95 | |
---|
| 96 | |
---|
| 97 | G4DynamicParticleVector * G4NeutronHPFSFissionFS::ApplyYourself(G4int nPrompt, |
---|
| 98 | G4int nDelayed, G4double * theDecayConst) |
---|
| 99 | { |
---|
| 100 | G4int i; |
---|
| 101 | G4DynamicParticleVector * aResult = new G4DynamicParticleVector; |
---|
| 102 | G4ReactionProduct boosted; |
---|
| 103 | boosted.Lorentz(theNeutron, theTarget); |
---|
| 104 | G4double eKinetic = boosted.GetKineticEnergy(); |
---|
| 105 | |
---|
| 106 | // Build neutrons |
---|
| 107 | G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed]; |
---|
| 108 | for(i=0; i<nPrompt+nDelayed; i++) |
---|
| 109 | { |
---|
| 110 | theNeutrons[i].SetDefinition(G4Neutron::Neutron()); |
---|
| 111 | } |
---|
| 112 | |
---|
| 113 | // sample energies |
---|
| 114 | G4int it, dummy; |
---|
| 115 | G4double tempE; |
---|
| 116 | for(i=0; i<nPrompt; i++) |
---|
| 117 | { |
---|
| 118 | tempE = thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab |
---|
| 119 | theNeutrons[i].SetKineticEnergy(tempE); |
---|
| 120 | } |
---|
| 121 | for(i=nPrompt; i<nPrompt+nDelayed; i++) |
---|
| 122 | { |
---|
| 123 | theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito |
---|
| 124 | if(it==0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy)); |
---|
| 125 | theDecayConst[i-nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned |
---|
| 126 | } |
---|
| 127 | |
---|
| 128 | // sample neutron angular distribution |
---|
| 129 | for(i=0; i<nPrompt+nDelayed; i++) |
---|
| 130 | { |
---|
| 131 | theNeutronAngularDis.SampleAndUpdate(theNeutrons[i]); // angular comes back in lab automatically |
---|
| 132 | } |
---|
| 133 | |
---|
| 134 | // already in lab. Add neutrons to dynamic particle vector |
---|
| 135 | for(i=0; i<nPrompt+nDelayed; i++) |
---|
| 136 | { |
---|
| 137 | G4DynamicParticle * it = new G4DynamicParticle; |
---|
| 138 | it->SetDefinition(theNeutrons[i].GetDefinition()); |
---|
| 139 | it->SetMomentum(theNeutrons[i].GetMomentum()); |
---|
| 140 | aResult->push_back(it); |
---|
| 141 | } |
---|
| 142 | delete [] theNeutrons; |
---|
| 143 | // return the result |
---|
| 144 | return aResult; |
---|
| 145 | } |
---|
| 146 | |
---|
| 147 | void G4NeutronHPFSFissionFS::SampleNeutronMult(G4int&all, G4int&Prompt, G4int&delayed, G4double eKinetic, G4int off) |
---|
| 148 | { |
---|
| 149 | G4double promptNeutronMulti = 0; |
---|
| 150 | promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic); |
---|
| 151 | G4double delayedNeutronMulti = 0; |
---|
| 152 | delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic); |
---|
| 153 | |
---|
| 154 | if(delayedNeutronMulti==0&&promptNeutronMulti==0) |
---|
| 155 | { |
---|
| 156 | Prompt = 0; |
---|
| 157 | delayed = 0; |
---|
| 158 | G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic); |
---|
| 159 | all = G4Poisson(totalNeutronMulti-off); |
---|
| 160 | all += off; |
---|
| 161 | } |
---|
| 162 | else |
---|
| 163 | { |
---|
| 164 | Prompt = G4Poisson(promptNeutronMulti-off); |
---|
| 165 | Prompt += off; |
---|
| 166 | delayed = G4Poisson(delayedNeutronMulti); |
---|
| 167 | all = Prompt+delayed; |
---|
| 168 | } |
---|
| 169 | } |
---|
| 170 | |
---|
| 171 | G4DynamicParticleVector * G4NeutronHPFSFissionFS::GetPhotons() |
---|
| 172 | { |
---|
| 173 | // sample photons |
---|
| 174 | G4ReactionProductVector * temp; |
---|
| 175 | G4ReactionProduct boosted; |
---|
| 176 | // the photon distributions are in the Nucleus rest frame. |
---|
| 177 | boosted.Lorentz(theNeutron, theTarget); |
---|
| 178 | G4double anEnergy = boosted.GetKineticEnergy(); |
---|
| 179 | temp = theFinalStatePhotons.GetPhotons(anEnergy); |
---|
| 180 | if(temp == 0) { return 0; } |
---|
| 181 | |
---|
| 182 | // lorentz transform, and add photons to final state |
---|
| 183 | unsigned int i; |
---|
| 184 | G4DynamicParticleVector * result = new G4DynamicParticleVector; |
---|
| 185 | for(i=0; i<temp->size(); i++) |
---|
| 186 | { |
---|
| 187 | // back to lab |
---|
| 188 | temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1.*theTarget); |
---|
| 189 | G4DynamicParticle * theOne = new G4DynamicParticle; |
---|
| 190 | theOne->SetDefinition(temp->operator[](i)->GetDefinition()); |
---|
| 191 | theOne->SetMomentum(temp->operator[](i)->GetMomentum()); |
---|
| 192 | result->push_back(theOne); |
---|
| 193 | delete temp->operator[](i); |
---|
| 194 | } |
---|
| 195 | delete temp; |
---|
| 196 | return result; |
---|
| 197 | } |
---|