[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 | // |
---|
[962] | 26 | // $Id: G4QLowEnergy.cc,v 1.7 2008/10/02 21:10:07 dennis Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
[819] | 28 | // |
---|
| 29 | // ---------------- G4QLowEnergy class ----------------- |
---|
| 30 | // by Mikhail Kossov, Aug 2007. |
---|
| 31 | // G4QLowEnergy class of the CHIPS Simulation Branch in GEANT4 |
---|
| 32 | // --------------------------------------------------------------- |
---|
| 33 | // **************************************************************************************** |
---|
| 34 | // ********** This CLASS is temporary moved from the "chips/interface" directory ********* |
---|
| 35 | // **************************************************************************************** |
---|
| 36 | |
---|
| 37 | //#define debug |
---|
| 38 | //#define pdebug |
---|
| 39 | //#define tdebug |
---|
| 40 | //#define nandebug |
---|
| 41 | //#define ppdebug |
---|
| 42 | |
---|
| 43 | #include "G4QLowEnergy.hh" |
---|
| 44 | |
---|
| 45 | // Initialization of static vectors |
---|
| 46 | G4int G4QLowEnergy::nPartCWorld=152; // #of particles initialized in CHIPS |
---|
| 47 | std::vector<G4int> G4QLowEnergy::ElementZ; // Z of element(i) in theLastCalc |
---|
| 48 | std::vector<G4double> G4QLowEnergy::ElProbInMat; // SumProbOfElem in Material |
---|
| 49 | std::vector<std::vector<G4int>*> G4QLowEnergy::ElIsoN;// N of isotope(j), E(i) |
---|
| 50 | std::vector<std::vector<G4double>*>G4QLowEnergy::IsoProbInEl;//SumProbIsotE(i) |
---|
| 51 | |
---|
| 52 | // Constructor |
---|
| 53 | G4QLowEnergy::G4QLowEnergy(const G4String& processName): |
---|
[962] | 54 | G4VDiscreteProcess(processName, fHadronic), evaporate(true) |
---|
[819] | 55 | { |
---|
| 56 | #ifdef debug |
---|
| 57 | G4cout<<"G4QLowEnergy::Constructor is called processName="<<processName<<G4endl; |
---|
| 58 | #endif |
---|
| 59 | if (verboseLevel>0) G4cout<<GetProcessName()<<" process is created "<<G4endl; |
---|
[962] | 60 | SetProcessSubType(fHadronInelastic); |
---|
[819] | 61 | G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max) |
---|
| 62 | } |
---|
| 63 | |
---|
| 64 | // Destructor |
---|
| 65 | G4QLowEnergy::~G4QLowEnergy() {} |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | G4LorentzVector G4QLowEnergy::GetEnegryMomentumConservation(){return EnMomConservation;} |
---|
| 69 | |
---|
| 70 | G4int G4QLowEnergy::GetNumberOfNeutronsInTarget() {return nOfNeutrons;} |
---|
| 71 | |
---|
| 72 | // output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)), |
---|
| 73 | // where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3) |
---|
| 74 | // ********** All CHIPS cross sections are calculated in the surface units ************ |
---|
| 75 | G4double G4QLowEnergy::GetMeanFreePath(const G4Track&Track, G4double, G4ForceCondition*F) |
---|
| 76 | { |
---|
| 77 | *F = NotForced; |
---|
| 78 | const G4DynamicParticle* incidentParticle = Track.GetDynamicParticle(); |
---|
| 79 | G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition(); |
---|
| 80 | if( !IsApplicable(*incidentParticleDefinition)) |
---|
| 81 | G4cout<<"-Warning-G4QLowEnergy::GetMeanFreePath for notImplemented Particle"<<G4endl; |
---|
| 82 | // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material |
---|
| 83 | G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle |
---|
| 84 | #ifdef debug |
---|
| 85 | G4double KinEn = incidentParticle->GetKineticEnergy(); |
---|
| 86 | G4cout<<"G4QLowEnergy::GetMeanFreePath:Prpj, kinE="<<KinEn<<", Mom="<<Momentum<<G4endl; |
---|
| 87 | #endif |
---|
| 88 | const G4Material* material = Track.GetMaterial(); // Get the current material |
---|
| 89 | const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume(); |
---|
| 90 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
| 91 | G4int nE=material->GetNumberOfElements(); |
---|
| 92 | #ifdef debug |
---|
| 93 | G4cout<<"G4QLowEnergy::GetMeanFreePath:"<<nE<<" Elems"<<G4endl; |
---|
| 94 | #endif |
---|
| 95 | G4VQCrossSection* CSmanager=G4QIonIonCrossSection::GetPointer(); |
---|
| 96 | G4int pPDG=0; |
---|
| 97 | // @@ At present it is made only for n & p, but can be extended if inXS are available |
---|
| 98 | if ( incidentParticleDefinition == G4Deuteron::Deuteron() ) pPDG = 100001002; |
---|
| 99 | else if ( incidentParticleDefinition == G4Alpha::Alpha() ) pPDG = 100002004; |
---|
| 100 | else if ( incidentParticleDefinition == G4Triton::Triton() ) pPDG = 100001003; |
---|
| 101 | else if ( incidentParticleDefinition == G4He3::He3() ) pPDG = 100002003; |
---|
| 102 | else if ( incidentParticleDefinition == G4GenericIon::GenericIon() ) |
---|
| 103 | { |
---|
| 104 | pPDG=incidentParticleDefinition->GetPDGEncoding(); |
---|
| 105 | #ifdef debug |
---|
| 106 | G4int B=incidentParticleDefinition->GetBaryonNumber(); |
---|
| 107 | G4int C=incidentParticleDefinition->GetPDGCharge(); |
---|
| 108 | prPDG=100000000+1000*C+B; |
---|
| 109 | G4cout<<"G4QIonIonElastic::GetMeanFreePath: PDG="<<prPDG<<"="<<pPDG<<G4endl; |
---|
| 110 | #endif |
---|
| 111 | } |
---|
| 112 | else G4cout<<"-Warning-G4QLowEnergy::GetMeanFreePath: only AA are implemented"<<G4endl; |
---|
| 113 | Momentum/=incidentParticleDefinition->GetBaryonNumber(); // Divide Mom by projectile A |
---|
| 114 | G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton |
---|
| 115 | G4double sigma=0.; // Sums over elements for the material |
---|
| 116 | G4int IPIE=IsoProbInEl.size(); // How many old elements? |
---|
| 117 | if(IPIE) for(G4int ip=0; ip<IPIE; ++ip) // Clean up the SumProb's of Isotopes (SPI) |
---|
| 118 | { |
---|
| 119 | std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector |
---|
| 120 | SPI->clear(); |
---|
| 121 | delete SPI; |
---|
| 122 | std::vector<G4int>* IsN=ElIsoN[ip]; // Pointer to the N vector |
---|
| 123 | IsN->clear(); |
---|
| 124 | delete IsN; |
---|
| 125 | } |
---|
| 126 | ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE) |
---|
| 127 | ElementZ.clear(); // Clear the body vector for Z of Elements |
---|
| 128 | IsoProbInEl.clear(); // Clear the body vector for SPI |
---|
| 129 | ElIsoN.clear(); // Clear the body vector for N of Isotopes |
---|
| 130 | for(G4int i=0; i<nE; ++i) |
---|
| 131 | { |
---|
| 132 | G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element |
---|
| 133 | G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element |
---|
| 134 | ElementZ.push_back(Z); // Remember Z of the Element |
---|
| 135 | G4int isoSize=0; // The default for the isoVectorLength is 0 |
---|
| 136 | G4int indEl=0; // Index of non-natural element or 0(default) |
---|
| 137 | G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect |
---|
| 138 | if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector |
---|
| 139 | #ifdef debug |
---|
| 140 | G4cout<<"G4QLowEnergy::GetMeanFreePath: isovector Length="<<isoSize<<G4endl; |
---|
| 141 | #endif |
---|
| 142 | if(isoSize) // The Element has non-trivial abundance set |
---|
| 143 | { |
---|
| 144 | indEl=pElement->GetIndex()+1; // Index of the non-trivial element is an order |
---|
| 145 | #ifdef debug |
---|
| 146 | G4cout<<"G4QLowEn::GetMFP:iE="<<indEl<<",def="<<Isotopes->IsDefined(Z,indEl)<<G4endl; |
---|
| 147 | #endif |
---|
| 148 | if(!Isotopes->IsDefined(Z,indEl)) // This index is not defined for this Z: define |
---|
| 149 | { |
---|
| 150 | std::vector<std::pair<G4int,G4double>*>* newAbund = |
---|
| 151 | new std::vector<std::pair<G4int,G4double>*>; |
---|
| 152 | G4double* abuVector=pElement->GetRelativeAbundanceVector(); |
---|
| 153 | for(G4int j=0; j<isoSize; j++) // Calculation of abundance vector for isotopes |
---|
| 154 | { |
---|
| 155 | G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z ! |
---|
| 156 | if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QDiffract::GetMeanFreePath: Z=" |
---|
| 157 | <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl; |
---|
| 158 | G4double abund=abuVector[j]; |
---|
| 159 | std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund); |
---|
| 160 | #ifdef debug |
---|
| 161 | G4cout<<"G4QLowEnergy::GetMeanFreePath:pair#"<<j<<",N="<<N<<",a="<<abund<<G4endl; |
---|
| 162 | #endif |
---|
| 163 | newAbund->push_back(pr); |
---|
| 164 | } |
---|
| 165 | #ifdef debug |
---|
| 166 | G4cout<<"G4QLowEnergy::GetMeanFreePath: pairVectLength="<<newAbund->size()<<G4endl; |
---|
| 167 | #endif |
---|
| 168 | indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd |
---|
| 169 | for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k]; // Cleaning temporary |
---|
| 170 | delete newAbund; // Was "new" in the beginning of the name space |
---|
| 171 | } |
---|
| 172 | } |
---|
| 173 | std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer |
---|
| 174 | std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector |
---|
| 175 | IsoProbInEl.push_back(SPI); |
---|
| 176 | std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector |
---|
| 177 | ElIsoN.push_back(IsN); |
---|
| 178 | G4int nIs=cs->size(); // A#Of Isotopes in the Element |
---|
| 179 | #ifdef debug |
---|
| 180 | G4cout<<"G4QLowEnergy::GetMFP:***=>,#isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl; |
---|
| 181 | #endif |
---|
| 182 | G4double susi=0.; // sum of CS over isotopes |
---|
| 183 | if(nIs) for(G4int j=0; j<nIs; j++) // Calculate CS for eachIsotope of El |
---|
| 184 | { |
---|
| 185 | std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice |
---|
| 186 | G4int N=curIs->first; // #of Neuterons in the isotope j of El i |
---|
| 187 | IsN->push_back(N); // Remember Min N for the Element |
---|
| 188 | #ifdef debug |
---|
| 189 | G4cout<<"G4QLowE::GMFP:true,P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl; |
---|
| 190 | #endif |
---|
| 191 | G4bool ccsf=true; // Extract inelastic Ion-Ion cross-section |
---|
| 192 | #ifdef debug |
---|
| 193 | G4cout<<"G4QLowEnergy::GMFP: GetCS #1 j="<<j<<G4endl; |
---|
| 194 | #endif |
---|
| 195 | G4double CSI=CSmanager->GetCrossSection(ccsf,Momentum,Z,N,pPDG);//CS(j,i) for isotope |
---|
| 196 | #ifdef debug |
---|
| 197 | G4cout<<"G4QLowEnergy::GetMeanFreePath: jI="<<j<<", Zt="<<Z<<", Nt="<<N<<", Mom=" |
---|
| 198 | <<Momentum<<", XSec="<<CSI/millibarn<<G4endl; |
---|
| 199 | #endif |
---|
| 200 | curIs->second = CSI; |
---|
| 201 | susi+=CSI; // Make a sum per isotopes |
---|
| 202 | SPI->push_back(susi); // Remember summed cross-section |
---|
| 203 | } // End of temporary initialization of the cross sections in the G4QIsotope singeltone |
---|
| 204 | sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV) |
---|
| 205 | #ifdef debug |
---|
| 206 | G4cout<<"G4QLowEnergy::GetMeanFreePath:<XS>="<<Isotopes->GetMeanCrossSection(Z,indEl) |
---|
| 207 | <<",AddSigm="<<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl; |
---|
| 208 | #endif |
---|
| 209 | ElProbInMat.push_back(sigma); |
---|
| 210 | } // End of LOOP over Elements |
---|
| 211 | // Check that cross section is not zero and return the mean free path |
---|
| 212 | #ifdef debug |
---|
| 213 | G4cout<<"G4QLowEnergy::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl; |
---|
| 214 | #endif |
---|
| 215 | if(sigma > 0.) return 1./sigma; // Mean path [distance] |
---|
| 216 | return DBL_MAX; |
---|
| 217 | } |
---|
| 218 | |
---|
| 219 | G4bool G4QLowEnergy::IsApplicable(const G4ParticleDefinition& particle) |
---|
| 220 | { |
---|
| 221 | if (particle == *( G4Proton::Proton() )) return true; |
---|
| 222 | else if (particle == *( G4Neutron::Neutron() )) return true; |
---|
| 223 | else if (particle == *( G4Deuteron::Deuteron() )) return true; |
---|
| 224 | else if (particle == *( G4Alpha::Alpha() )) return true; |
---|
| 225 | else if (particle == *( G4Triton::Triton() )) return true; |
---|
| 226 | else if (particle == *( G4He3::He3() )) return true; |
---|
| 227 | else if (particle == *( G4GenericIon::GenericIon() )) return true; |
---|
| 228 | #ifdef debug |
---|
| 229 | G4cout<<"***>>G4QLowEnergy::IsApplicable: projPDG="<<particle.GetPDGEncoding()<<G4endl; |
---|
| 230 | #endif |
---|
| 231 | return false; |
---|
| 232 | } |
---|
| 233 | |
---|
| 234 | G4VParticleChange* G4QLowEnergy::PostStepDoIt(const G4Track& track, const G4Step& step) |
---|
| 235 | { |
---|
| 236 | static const G4double mProt= G4QPDGCode(2212).GetMass()/MeV; // CHIPS proton Mass in MeV |
---|
| 237 | static const G4double mNeut= G4QPDGCode(2112).GetMass()/MeV; // CHIPS neutron Mass in MeV |
---|
| 238 | static const G4double mLamb= G4QPDGCode(3122).GetMass()/MeV; // CHIPS Lambda Mass in MeV |
---|
| 239 | static const G4double mDeut= G4QPDGCode(2112).GetNuclMass(1,1,0)/MeV; |
---|
| 240 | static const G4double mTrit= G4QPDGCode(2112).GetNuclMass(1,2,0)/MeV; |
---|
| 241 | static const G4double mHel3= G4QPDGCode(2112).GetNuclMass(2,1,0)/MeV; |
---|
| 242 | static const G4double mAlph= G4QPDGCode(2112).GetNuclMass(2,2,0)/MeV; |
---|
| 243 | static const G4ThreeVector zeroMom(0.,0.,0.); |
---|
| 244 | static G4ParticleDefinition* aGamma = G4Gamma::Gamma(); |
---|
| 245 | static G4ParticleDefinition* aProton = G4Proton::Proton(); |
---|
| 246 | static G4ParticleDefinition* aNeutron = G4Neutron::Neutron(); |
---|
| 247 | static G4ParticleDefinition* aLambda = G4Lambda::Lambda(); |
---|
| 248 | static G4ParticleDefinition* aDeuteron = G4Deuteron::Deuteron(); |
---|
| 249 | static G4ParticleDefinition* aTriton = G4Triton::Triton(); |
---|
| 250 | static G4ParticleDefinition* aHe3 = G4He3::He3(); |
---|
| 251 | static G4ParticleDefinition* anAlpha = G4Alpha::Alpha(); |
---|
| 252 | static const G4int nCh=26; // #of combinations |
---|
| 253 | static G4QNucleus Nuc; // A fake nucleus to call Evaporation |
---|
| 254 | // |
---|
| 255 | //------------------------------------------------------------------------------------- |
---|
| 256 | static G4bool CWinit = true; // CHIPS Warld needs to be initted |
---|
| 257 | if(CWinit) |
---|
| 258 | { |
---|
| 259 | CWinit=false; |
---|
| 260 | G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max) |
---|
| 261 | } |
---|
| 262 | //------------------------------------------------------------------------------------- |
---|
| 263 | const G4DynamicParticle* projHadron = track.GetDynamicParticle(); |
---|
| 264 | const G4ParticleDefinition* particle=projHadron->GetDefinition(); |
---|
| 265 | #ifdef debug |
---|
| 266 | G4cout<<"G4QLowEnergy::PostStepDoIt: Before the GetMeanFreePath is called In4M=" |
---|
| 267 | <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type=" |
---|
| 268 | <<particle->GetParticleType()<<",SubType="<<particle->GetParticleSubType()<<G4endl; |
---|
| 269 | #endif |
---|
| 270 | //G4ForceCondition cond=NotForced; |
---|
| 271 | //GetMeanFreePath(track, -27., &cond); // @@ ?? jus to update parameters? |
---|
| 272 | #ifdef debug |
---|
| 273 | G4cout<<"G4QLowEnergy::PostStepDoIt: After GetMeanFreePath is called"<<G4endl; |
---|
| 274 | #endif |
---|
| 275 | G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV! |
---|
| 276 | G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV |
---|
| 277 | G4double Momentum = proj4M.rho(); // @@ Just for the test purposes |
---|
| 278 | if(std::fabs(Momentum-momentum)>.000001) |
---|
| 279 | G4cerr<<"-Warning-G4QLowEnergy::PostStepDoIt:P_IU="<<Momentum<<"#"<<momentum<<G4endl; |
---|
| 280 | #ifdef pdebug |
---|
| 281 | G4cout<<"G4QLowEnergy::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum<<",proj4M/m=" |
---|
| 282 | <<proj4M<<proj4M.m()<<G4endl; |
---|
| 283 | #endif |
---|
| 284 | if (!IsApplicable(*particle)) // Check applicability |
---|
| 285 | { |
---|
| 286 | G4cerr<<"G4QLowEnergy::PostStepDoIt: Only NA is implemented."<<G4endl; |
---|
| 287 | return 0; |
---|
| 288 | } |
---|
| 289 | const G4Material* material = track.GetMaterial(); // Get the current material |
---|
| 290 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
| 291 | G4int nE=material->GetNumberOfElements(); |
---|
| 292 | #ifdef debug |
---|
| 293 | G4cout<<"G4QLowEnergy::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl; |
---|
| 294 | #endif |
---|
| 295 | G4int projPDG=0; // PDG Code prototype for the captured hadron |
---|
| 296 | // Not all these particles are implemented yet (see Is Applicable) |
---|
| 297 | if (particle == G4Proton::Proton() ) projPDG= 2212; |
---|
| 298 | else if (particle == G4Neutron::Neutron() ) projPDG= 2112; |
---|
| 299 | else if (particle == G4Deuteron::Deuteron() ) projPDG= 100001002; |
---|
| 300 | else if (particle == G4Alpha::Alpha() ) projPDG= 100002004; |
---|
| 301 | else if (particle == G4Triton::Triton() ) projPDG= 100001003; |
---|
| 302 | else if (particle == G4He3::He3() ) projPDG= 100002003; |
---|
| 303 | else if (particle == G4GenericIon::GenericIon() ) |
---|
| 304 | { |
---|
| 305 | projPDG=particle->GetPDGEncoding(); |
---|
| 306 | #ifdef debug |
---|
| 307 | G4int B=particle->GetBaryonNumber(); |
---|
| 308 | G4int C=particle->GetPDGCharge(); |
---|
| 309 | prPDG=100000000+1000*C+B; |
---|
| 310 | G4cout<<"G4QLowEnergy::PostStepDoIt: PDG="<<prPDG<<"="<<projPDG<<G4endl; |
---|
| 311 | #endif |
---|
| 312 | } |
---|
| 313 | else G4cout<<"-Warning-G4QLowEnergy::PostStepDoIt:Unknown projectile Ion"<<G4endl; |
---|
| 314 | #ifdef debug |
---|
| 315 | G4int prPDG=particle->GetPDGEncoding(); |
---|
| 316 | G4cout<<"G4QLowEnergy::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl; |
---|
| 317 | #endif |
---|
| 318 | if(!projPDG) |
---|
| 319 | { |
---|
| 320 | G4cerr<<"-Warning-G4QLowEnergy::PostStepDoIt:UndefProjHadron(PDG=0) ->ret 0"<<G4endl; |
---|
| 321 | return 0; |
---|
| 322 | } |
---|
| 323 | // Element treatment |
---|
| 324 | G4int EPIM=ElProbInMat.size(); |
---|
| 325 | #ifdef debug |
---|
| 326 | G4cout<<"G4QLowEn::PostStDoIt: m="<<EPIM<<", n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl; |
---|
| 327 | #endif |
---|
| 328 | G4int i=0; |
---|
| 329 | if(EPIM>1) |
---|
| 330 | { |
---|
| 331 | G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand(); |
---|
| 332 | for(i=0; i<nE; ++i) |
---|
| 333 | { |
---|
| 334 | #ifdef debug |
---|
| 335 | G4cout<<"G4QLowEn::PostStepDoIt: EPM["<<i<<"]="<<ElProbInMat[i]<<", r="<<rnd<<G4endl; |
---|
| 336 | #endif |
---|
| 337 | if (rnd<ElProbInMat[i]) break; |
---|
| 338 | } |
---|
| 339 | if(i>=nE) i=nE-1; // Top limit for the Element |
---|
| 340 | } |
---|
| 341 | G4Element* pElement=(*theElementVector)[i]; |
---|
| 342 | G4int tZ=static_cast<G4int>(pElement->GetZ()); |
---|
| 343 | #ifdef debug |
---|
| 344 | G4cout<<"G4QLowEnergy::PostStepDoIt: i="<<i<<", Z(element)="<<tZ<<G4endl; |
---|
| 345 | #endif |
---|
| 346 | if(tZ<=0) |
---|
| 347 | { |
---|
| 348 | G4cerr<<"-Warning-G4QLowEnergy::PostStepDoIt: Element with Z="<<tZ<<G4endl; |
---|
| 349 | if(tZ<0) return 0; |
---|
| 350 | } |
---|
| 351 | std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes |
---|
| 352 | std::vector<G4int>* IsN = ElIsoN[i]; // Vector of "#of neutrons" in the isotope El[i] |
---|
| 353 | G4int nofIsot=SPI->size(); // #of isotopes in the element i |
---|
| 354 | #ifdef debug |
---|
| 355 | G4cout<<"G4QLowEnergy::PostStepDoIt: nI="<<nofIsot<<", T="<<(*SPI)[nofIsot-1]<<G4endl; |
---|
| 356 | #endif |
---|
| 357 | G4int j=0; |
---|
| 358 | if(nofIsot>1) |
---|
| 359 | { |
---|
| 360 | G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element |
---|
| 361 | for(j=0; j<nofIsot; ++j) |
---|
| 362 | { |
---|
| 363 | #ifdef debug |
---|
| 364 | G4cout<<"G4QLowEnergy::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<",r="<<rndI<<G4endl; |
---|
| 365 | #endif |
---|
| 366 | if(rndI < (*SPI)[j]) break; |
---|
| 367 | } |
---|
| 368 | if(j>=nofIsot) j=nofIsot-1; // Top limit for the isotope |
---|
| 369 | } |
---|
| 370 | G4int tN =(*IsN)[j]; ; // Randomized number of neutrons |
---|
| 371 | #ifdef debug |
---|
| 372 | G4cout<<"G4QLowEnergy::PostStepDoIt: j="<<i<<", N(isotope)="<<tN<<", MeV="<<MeV<<G4endl; |
---|
| 373 | #endif |
---|
| 374 | if(tN<0) |
---|
| 375 | { |
---|
| 376 | G4cerr<<"-Warning-G4QLowEnergy::PostStepDoIt: Isotope Z="<<tZ<<" has 0>N="<<tN<<G4endl; |
---|
| 377 | return 0; |
---|
| 378 | } |
---|
| 379 | nOfNeutrons=tN; // Remember it for the energy-momentum check |
---|
| 380 | #ifdef debug |
---|
| 381 | G4cout<<"G4QLowEnergy::PostStepDoIt: N="<<tN<<" for element with Z="<<tZ<<G4endl; |
---|
| 382 | #endif |
---|
| 383 | if(tN<0) |
---|
| 384 | { |
---|
| 385 | G4cerr<<"*Warning*G4QLowEnergy::PostStepDoIt:Element with N="<<tN<< G4endl; |
---|
| 386 | return 0; |
---|
| 387 | } |
---|
| 388 | aParticleChange.Initialize(track); |
---|
| 389 | #ifdef debug |
---|
| 390 | G4cout<<"G4QLowEnergy::PostStepDoIt: track is initialized"<<G4endl; |
---|
| 391 | #endif |
---|
| 392 | G4double weight = track.GetWeight(); |
---|
| 393 | G4double localtime = track.GetGlobalTime(); |
---|
| 394 | G4ThreeVector position = track.GetPosition(); |
---|
| 395 | #ifdef debug |
---|
| 396 | G4cout<<"G4QLowEnergy::PostStepDoIt: before Touchable extraction"<<G4endl; |
---|
| 397 | #endif |
---|
| 398 | G4TouchableHandle trTouchable = track.GetTouchableHandle(); |
---|
| 399 | #ifdef debug |
---|
| 400 | G4cout<<"G4QLowEnergy::PostStepDoIt: Touchable is extracted"<<G4endl; |
---|
| 401 | #endif |
---|
| 402 | G4QPDGCode targQPDG(90000000+tZ*1000+tN); // @@ use G4Ion and get rid of CHIPS World |
---|
| 403 | G4double tM=targQPDG.GetMass(); // CHIPS target nucleus mass in MeV |
---|
| 404 | G4int pL=particle->GetQuarkContent(3)-particle->GetAntiQuarkContent(3); // Strangeness |
---|
| 405 | G4int pZ=static_cast<G4int>(particle->GetPDGCharge()); // Charge of the projectile |
---|
| 406 | G4int pN=particle->GetBaryonNumber()-pZ-pL;// #of neutrons in projectile |
---|
| 407 | G4double pM=targQPDG.GetNuclMass(pZ,pN,0); // CHIPS projectile nucleus mass in MeV |
---|
| 408 | G4double cosp=-14*Momentum*(tM-pM)/tM/pM; // Asymmetry power for angular distribution |
---|
| 409 | #ifdef debug |
---|
| 410 | G4cout<<"G4QLowEnergy::PoStDoIt: Proj("<<pZ<<","<<pN<<","<<pL<<")p="<<pM<<",Targ=("<<tZ |
---|
| 411 | <<","<<tN<<"), cosp="<<cosp<<G4endl; |
---|
| 412 | #endif |
---|
| 413 | G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?) |
---|
| 414 | G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector |
---|
| 415 | G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction |
---|
| 416 | #ifdef debug |
---|
| 417 | G4cout<<"G4QLowEnergy::PostStepDoIt: tM="<<tM<<",p4M="<<proj4M<<",t4M="<<tot4M<<G4endl; |
---|
| 418 | #endif |
---|
| 419 | EnMomConservation=tot4M; // Total 4-mom of reaction for E/M conservation |
---|
| 420 | // @@ Probably this is not necessary any more |
---|
| 421 | #ifdef debug |
---|
| 422 | G4cout<<"G4QLE::PSDI:false,P="<<Momentum<<",Z="<<pZ<<",N="<<pN<<",PDG="<<projPDG<<G4endl; |
---|
| 423 | #endif |
---|
| 424 | G4double xSec=CalculateXS(Momentum, tZ, tN, projPDG); // Recalculate CrossSection |
---|
| 425 | #ifdef debug |
---|
| 426 | G4cout<<"G4QLowEn::PSDI:PDG="<<projPDG<<",P="<<Momentum<<",XS="<<xSec/millibarn<<G4endl; |
---|
| 427 | #endif |
---|
| 428 | #ifdef nandebug |
---|
| 429 | if(xSec>0. || xSec<0. || xSec==0); |
---|
| 430 | else G4cout<<"-Warning-G4QLowEnergy::PostSDI: *NAN* xSec="<<xSec/millibarn<<G4endl; |
---|
| 431 | #endif |
---|
| 432 | // @@ check a possibility to separate p, n, or alpha (!) |
---|
| 433 | if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing |
---|
| 434 | { |
---|
| 435 | #ifdef pdebug |
---|
| 436 | G4cerr<<"-Warning-G4QLowEnergy::PSDoIt:*Zero cross-section* PDG="<<projPDG |
---|
| 437 | <<",Z="<<tZ<<",tN="<<tN<<",P="<<Momentum<<G4endl; |
---|
| 438 | #endif |
---|
| 439 | //Do Nothing Action insead of the reaction |
---|
| 440 | aParticleChange.ProposeEnergy(kinEnergy); |
---|
| 441 | aParticleChange.ProposeLocalEnergyDeposit(0.); |
---|
| 442 | aParticleChange.ProposeMomentumDirection(dir) ; |
---|
| 443 | return G4VDiscreteProcess::PostStepDoIt(track,step); |
---|
| 444 | } |
---|
| 445 | // Kill interacting hadron |
---|
| 446 | aParticleChange.ProposeEnergy(0.); |
---|
| 447 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
| 448 | #ifdef debug |
---|
| 449 | G4cout<<"G4QLowEn::PSDI: Projectile track is killed"<<G4endl; |
---|
| 450 | #endif |
---|
| 451 | // algorithm implementation --- STARTS HERE --- All calculations are in IU -------- |
---|
| 452 | G4double totM=tot4M.m(); // total CMS mass of the reaction |
---|
| 453 | G4int totN=tN+pN; |
---|
| 454 | G4int totZ=tZ+pZ; |
---|
| 455 | // @@ Here mass[i] can be calculated if mass=0 |
---|
| 456 | G4double mass[nCh]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0., |
---|
| 457 | 0.,0.,0.,0.,0.,0.}; |
---|
| 458 | mass[0] = targQPDG.GetNuclMass(totZ,totN,0); // gamma+gamma |
---|
| 459 | #ifdef debug |
---|
| 460 | G4cout<<"G4QLowEn::PSDI: Nuclear Mass="<<mass[0]<<G4endl; |
---|
| 461 | #endif |
---|
| 462 | if (totN>0 && totZ>0) |
---|
| 463 | { |
---|
| 464 | mass[1] = targQPDG.GetNuclMass(totZ,totN-1,0); // gamma+neutron |
---|
| 465 | mass[2] = targQPDG.GetNuclMass(totZ-1,totN,0); // gamma+proton |
---|
| 466 | } |
---|
| 467 | if ( (totZ > 0 && totN > 1) || (totN > 0 && totZ > 1) ) |
---|
| 468 | mass[3] = targQPDG.GetNuclMass(totZ-1,totN-1,0); //g+d |
---|
| 469 | if ( (totZ > 0 && totN > 2) || (totN > 1 && totZ > 1) ) |
---|
| 470 | mass[4] = targQPDG.GetNuclMass(totZ-1,totN-2,0); //g+t |
---|
| 471 | if ( (totZ > 2 && totN > 0) || (totN > 1 && totZ > 1) ) |
---|
| 472 | mass[5] = targQPDG.GetNuclMass(totZ-2,totN-1,0); //g+3 |
---|
| 473 | if ( (totZ > 2 && totN > 1) || (totN > 2 && totZ > 1) ) |
---|
| 474 | mass[6] = targQPDG.GetNuclMass(totZ-2,totN-2,0); //g+a |
---|
| 475 | if ( (totZ > 0 && totN > 1) || totN > 2 ) |
---|
| 476 | mass[7] = targQPDG.GetNuclMass(totZ,totN-2,0); // n+n |
---|
| 477 | mass[ 8] = mass[3]; // neutron+proton (the same as a deuteron) |
---|
| 478 | if ( (totZ > 1 && totN > 0) || totZ > 2 ) |
---|
| 479 | mass[9] = targQPDG.GetNuclMass(totZ-2,totN,0); // p+p |
---|
| 480 | mass[10] = mass[5]; // proton+deuteron (the same as He3) |
---|
| 481 | mass[11] = mass[4]; // neutron+deuteron (the same as t) |
---|
| 482 | mass[12] = mass[6]; // deuteron+deuteron (the same as alpha) |
---|
| 483 | mass[13] = mass[6]; // proton+tritium (the same as alpha) |
---|
| 484 | if ( totN > 3 || (totN > 2 && totZ > 0) ) |
---|
| 485 | mass[14] = targQPDG.GetNuclMass(totZ-1,totN-3,0); // n+t |
---|
| 486 | if ( totZ > 3 || (totZ > 2 && totN > 0) ) |
---|
| 487 | mass[15] = targQPDG.GetNuclMass(totZ-3,totN-1,0); // He3+p |
---|
| 488 | mass[16] = mass[6]; // neutron+He3 (the same as alpha) |
---|
| 489 | if ( (totZ > 3 && totN > 1) || (totZ > 2 && totN > 2) ) |
---|
| 490 | mass[17] = targQPDG.GetNuclMass(totZ-3,totN-2,0); // pa |
---|
| 491 | if ( (totZ > 1 && totN > 3) || (totZ > 2 && totN > 2) ) |
---|
| 492 | mass[18] = targQPDG.GetNuclMass(totZ-2,totN-3,0); // na |
---|
| 493 | if(pL>0) |
---|
| 494 | { |
---|
| 495 | G4int pL1=pL-1; |
---|
| 496 | if(totN>0||totZ>0) mass[19] = targQPDG.GetNuclMass(totZ ,totN ,pL1);// Lambda+gamma |
---|
| 497 | if( (totN > 0 && totZ > 0) || totZ > 1 ) |
---|
| 498 | mass[20]=targQPDG.GetNuclMass(totZ-1,totN ,pL1);//Lp |
---|
| 499 | if( (totN > 0 && totZ > 0) || totN > 0 ) |
---|
| 500 | mass[21]=targQPDG.GetNuclMass(totZ ,totN-1,pL1);//Ln |
---|
| 501 | if( (totN > 1 && totZ > 0) || (totN > 0 && totZ > 1) ) |
---|
| 502 | mass[22]=targQPDG.GetNuclMass(totZ-1,totN-1,pL1);//Ld |
---|
| 503 | if( (totN > 2 && totZ > 0) || (totN > 1 && totZ > 1) ) |
---|
| 504 | mass[23]=targQPDG.GetNuclMass(totZ-1,totN-2,pL1);//Lt |
---|
| 505 | if( (totN > 0 && totZ > 2) || (totN > 1 && totZ > 1) ) |
---|
| 506 | mass[24]=targQPDG.GetNuclMass(totZ-2,totN-1,pL1);//L3 |
---|
| 507 | if( (totN > 1 && totZ > 2) || (totN > 2 && totZ > 1) ) |
---|
| 508 | mass[25]=targQPDG.GetNuclMass(totZ-2,totN-2,pL1);//La |
---|
| 509 | } |
---|
| 510 | #ifdef debug |
---|
| 511 | G4cout<<"G4QLowEn::PSDI: Residual masses are calculated"<<G4endl; |
---|
| 512 | #endif |
---|
| 513 | G4double tA=tZ+tN; |
---|
| 514 | G4double pA=pZ+pN; |
---|
| 515 | G4double prZ=pZ/pA+tZ/tA; |
---|
| 516 | G4double prN=pN/pA+tN/tA; |
---|
| 517 | G4double prD=prN*prZ; |
---|
| 518 | G4double prA=prD*prD; |
---|
| 519 | G4double prH=prD*prZ; |
---|
| 520 | G4double prT=prD*prN; |
---|
| 521 | G4double fhe3=6.*std::sqrt(tA); |
---|
| 522 | G4double prL=0.; |
---|
| 523 | if(pL>0) prL=pL/pA; |
---|
| 524 | G4double qval[nCh]={0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0., |
---|
| 525 | 0.,0.,0.,0.,0.,0.}; |
---|
| 526 | qval[ 0] = totM - mass[ 0]; |
---|
| 527 | qval[ 1] = (totM - mass[ 1] - mNeut)*prN; |
---|
| 528 | qval[ 2] = (totM - mass[ 2] - mProt)*prZ; |
---|
| 529 | qval[ 3] = (totM - mass[ 3] - mDeut)*prD/3.; |
---|
| 530 | qval[ 4] = (totM - mass[ 4] - mTrit)*prT/6.; |
---|
| 531 | qval[ 5] = (totM - mass[ 5] - mHel3)*prH/fhe3; |
---|
| 532 | qval[ 6] = (totM - mass[ 6] - mAlph)*prA/9.; |
---|
| 533 | qval[ 7] = (totM - mass[ 7] - mNeut - mNeut)*prN*prN; |
---|
| 534 | qval[ 8] = (totM - mass[ 8] - mNeut - mProt)*prD; |
---|
| 535 | qval[ 9] = (totM - mass[ 9] - mProt - mProt)*prZ*prZ; |
---|
| 536 | qval[10] = (totM - mass[10] - mProt - mDeut)*prH/3.; |
---|
| 537 | qval[11] = (totM - mass[11] - mNeut - mDeut)*prT/3.; |
---|
| 538 | qval[12] = (totM - mass[12] - mDeut - mDeut)*prA/3./3.; |
---|
| 539 | qval[13] = (totM - mass[13] - mProt - mTrit)*prA/6.; |
---|
| 540 | qval[14] = (totM - mass[14] - mNeut - mTrit)*prT*prN/6.; |
---|
| 541 | qval[15] = (totM - mass[15] - mProt - mHel3)*prH*prZ/fhe3; |
---|
| 542 | qval[16] = (totM - mass[16] - mNeut - mHel3)*prA/fhe3; |
---|
| 543 | qval[17] = (totM - mass[17] - mProt - mAlph)*prZ*prA/9.; |
---|
| 544 | qval[18] = (totM - mass[18] - mNeut - mAlph)*prN*prA/9.; |
---|
| 545 | if(pZ>0) |
---|
| 546 | { |
---|
| 547 | qval[19] = (totM - mass[19] - mLamb)*prL; |
---|
| 548 | qval[20] = (totM - mass[20] - mProt - mLamb)*prL*prZ; |
---|
| 549 | qval[21] = (totM - mass[21] - mNeut - mLamb)*prL*prN; |
---|
| 550 | qval[22] = (totM - mass[22] - mDeut - mLamb)*prL*prD/2.; |
---|
| 551 | qval[23] = (totM - mass[23] - mTrit - mLamb)*prL*prT/3.; |
---|
| 552 | qval[24] = (totM - mass[24] - mHel3 - mLamb)*prL*prH/fhe3; |
---|
| 553 | qval[25] = (totM - mass[25] - mAlph - mLamb)*prL*prA/4; |
---|
| 554 | } |
---|
| 555 | #ifdef debug |
---|
| 556 | G4cout<<"G4QLowEn::PSDI: Q-values are calculated, tgA="<<tA<<"prA="<<pA<<G4endl; |
---|
| 557 | #endif |
---|
| 558 | |
---|
| 559 | if( !pZ && pN==1) |
---|
| 560 | { |
---|
| 561 | if(G4UniformRand()>(tA-1.)*(tA-1.)/52900.) qval[0] = 0.0; |
---|
| 562 | if(G4UniformRand()>kinEnergy/7.925*tA) qval[2]=qval[3]=qval[4]=qval[5]=qval[9]=0.; |
---|
| 563 | } |
---|
| 564 | else qval[0] = 0.0; |
---|
| 565 | |
---|
| 566 | G4double qv = 0.0; // Total sum of probabilities (q-values) |
---|
| 567 | for(G4int i=0; i<nCh; ++i ) |
---|
| 568 | { |
---|
| 569 | if( mass[i] < 500.*MeV ) qval[i] = 0.0; // Close A/Z impossible channels |
---|
| 570 | if( qval[i] < 0.0 ) qval[i] = 0.0; // Close the splitting impossible channels |
---|
| 571 | qv += qval[i]; |
---|
| 572 | } |
---|
| 573 | // Select the channel |
---|
| 574 | G4double qv1 = 0.0; |
---|
| 575 | G4double ran = G4UniformRand(); |
---|
| 576 | G4int index = 0; |
---|
| 577 | for( index=0; index<nCh; ++index ) |
---|
| 578 | { |
---|
| 579 | if( qval[index] > 0.0 ) |
---|
| 580 | { |
---|
| 581 | qv1 += qval[index]/qv; |
---|
| 582 | if( ran <= qv1 ) break; |
---|
| 583 | } |
---|
| 584 | } |
---|
| 585 | #ifdef debug |
---|
| 586 | G4cout<<"G4QLowEn::PSDI: index="<<index<<" < "<<nCh<<G4endl; |
---|
| 587 | #endif |
---|
| 588 | if(index == nCh) |
---|
| 589 | { |
---|
| 590 | G4cout<<"***G4QLowEnergy::PoStDI:Decay is impossible,totM="<<totM<<",GSM="<<tM<<G4endl; |
---|
| 591 | throw G4QException("G4QLowEnergy::PostStepDoIt: Can't decay the Compound"); |
---|
| 592 | } |
---|
| 593 | // @@ Convert it to G4QHadrons |
---|
| 594 | G4DynamicParticle* ResSec = new G4DynamicParticle; |
---|
| 595 | G4DynamicParticle* FstSec = new G4DynamicParticle; |
---|
| 596 | G4DynamicParticle* SecSec = new G4DynamicParticle; |
---|
| 597 | #ifdef debug |
---|
| 598 | G4cout<<"G4QLowEn::PSDI: Dynamic particles are created pL="<<pL<<G4endl; |
---|
| 599 | #endif |
---|
| 600 | |
---|
| 601 | G4LorentzVector res4Mom(zeroMom,mass[index]*MeV); // The recoil nucleus prototype |
---|
| 602 | G4double mF=0.; |
---|
| 603 | G4double mS=0.; |
---|
| 604 | G4int rA=totZ+totN; |
---|
| 605 | G4int rZ=totZ; |
---|
| 606 | G4int rL=pL; |
---|
| 607 | G4int complete=3; |
---|
| 608 | G4ParticleDefinition* theDefinition; // Prototype for qfNucleon |
---|
| 609 | switch( index ) |
---|
| 610 | { |
---|
| 611 | case 0: |
---|
| 612 | if(!evaporate || rA<2) |
---|
| 613 | { |
---|
| 614 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 615 | if(!theDefinition) |
---|
| 616 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 617 | ResSec->SetDefinition( theDefinition ); |
---|
| 618 | FstSec->SetDefinition( aGamma ); |
---|
| 619 | SecSec->SetDefinition( aGamma ); |
---|
| 620 | } |
---|
| 621 | else |
---|
| 622 | { |
---|
| 623 | delete ResSec; |
---|
| 624 | delete FstSec; |
---|
| 625 | delete SecSec; |
---|
| 626 | complete=0; |
---|
| 627 | } |
---|
| 628 | break; |
---|
| 629 | case 1: |
---|
| 630 | rA-=1; |
---|
| 631 | if(!evaporate || rA<2) |
---|
| 632 | { |
---|
| 633 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 634 | if(!theDefinition) |
---|
| 635 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 636 | ResSec->SetDefinition( theDefinition ); |
---|
| 637 | SecSec->SetDefinition( aGamma ); |
---|
| 638 | } |
---|
| 639 | else |
---|
| 640 | { |
---|
| 641 | delete ResSec; |
---|
| 642 | delete SecSec; |
---|
| 643 | complete=1; |
---|
| 644 | } |
---|
| 645 | FstSec->SetDefinition( aNeutron ); |
---|
| 646 | mF=mNeut; // First hadron 4-momentum |
---|
| 647 | break; |
---|
| 648 | case 2: |
---|
| 649 | rA-=1; |
---|
| 650 | rZ-=1; |
---|
| 651 | if(!evaporate || rA<2) |
---|
| 652 | { |
---|
| 653 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 654 | if(!theDefinition) |
---|
| 655 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 656 | ResSec->SetDefinition( theDefinition ); |
---|
| 657 | SecSec->SetDefinition( aGamma ); |
---|
| 658 | } |
---|
| 659 | else |
---|
| 660 | { |
---|
| 661 | delete ResSec; |
---|
| 662 | delete SecSec; |
---|
| 663 | complete=1; |
---|
| 664 | } |
---|
| 665 | FstSec->SetDefinition( aProton ); |
---|
| 666 | mF=mProt; // First hadron 4-momentum |
---|
| 667 | break; |
---|
| 668 | case 3: |
---|
| 669 | rA-=2; |
---|
| 670 | rZ-=1; |
---|
| 671 | if(!evaporate || rA<2) |
---|
| 672 | { |
---|
| 673 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 674 | if(!theDefinition) |
---|
| 675 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 676 | ResSec->SetDefinition( theDefinition ); |
---|
| 677 | SecSec->SetDefinition( aGamma ); |
---|
| 678 | } |
---|
| 679 | else |
---|
| 680 | { |
---|
| 681 | delete ResSec; |
---|
| 682 | delete SecSec; |
---|
| 683 | complete=1; |
---|
| 684 | } |
---|
| 685 | FstSec->SetDefinition( aDeuteron ); |
---|
| 686 | mF=mDeut; // First hadron 4-momentum |
---|
| 687 | break; |
---|
| 688 | case 4: |
---|
| 689 | rA-=3; |
---|
| 690 | rZ-=1; |
---|
| 691 | if(!evaporate || rA<2) |
---|
| 692 | { |
---|
| 693 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 694 | if(!theDefinition) |
---|
| 695 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 696 | ResSec->SetDefinition( theDefinition ); |
---|
| 697 | SecSec->SetDefinition( aGamma ); |
---|
| 698 | } |
---|
| 699 | else |
---|
| 700 | { |
---|
| 701 | delete ResSec; |
---|
| 702 | delete SecSec; |
---|
| 703 | complete=1; |
---|
| 704 | } |
---|
| 705 | FstSec->SetDefinition( aTriton ); |
---|
| 706 | mF=mTrit; // First hadron 4-momentum |
---|
| 707 | break; |
---|
| 708 | case 5: |
---|
| 709 | rA-=3; |
---|
| 710 | rZ-=2; |
---|
| 711 | if(!evaporate || rA<2) |
---|
| 712 | { |
---|
| 713 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 714 | if(!theDefinition) |
---|
| 715 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 716 | ResSec->SetDefinition( theDefinition ); |
---|
| 717 | SecSec->SetDefinition( aGamma ); |
---|
| 718 | } |
---|
| 719 | else |
---|
| 720 | { |
---|
| 721 | delete ResSec; |
---|
| 722 | delete SecSec; |
---|
| 723 | complete=1; |
---|
| 724 | } |
---|
| 725 | FstSec->SetDefinition( aHe3); |
---|
| 726 | mF=mHel3; // First hadron 4-momentum |
---|
| 727 | break; |
---|
| 728 | case 6: |
---|
| 729 | rA-=4; |
---|
| 730 | rZ-=2; |
---|
| 731 | if(!evaporate || rA<2) |
---|
| 732 | { |
---|
| 733 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 734 | if(!theDefinition) |
---|
| 735 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 736 | ResSec->SetDefinition( theDefinition ); |
---|
| 737 | SecSec->SetDefinition( aGamma ); |
---|
| 738 | } |
---|
| 739 | else |
---|
| 740 | { |
---|
| 741 | delete ResSec; |
---|
| 742 | delete SecSec; |
---|
| 743 | complete=1; |
---|
| 744 | } |
---|
| 745 | FstSec->SetDefinition( anAlpha ); |
---|
| 746 | mF=mAlph; // First hadron 4-momentum |
---|
| 747 | break; |
---|
| 748 | case 7: |
---|
| 749 | rA-=2; |
---|
| 750 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 751 | if(!theDefinition) |
---|
| 752 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 753 | ResSec->SetDefinition( theDefinition ); |
---|
| 754 | FstSec->SetDefinition( aNeutron ); |
---|
| 755 | SecSec->SetDefinition( aNeutron ); |
---|
| 756 | mF=mNeut; // First hadron 4-momentum |
---|
| 757 | mS=mNeut; // Second hadron 4-momentum |
---|
| 758 | break; |
---|
| 759 | case 8: |
---|
| 760 | rZ-=1; |
---|
| 761 | rA-=2; |
---|
| 762 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 763 | if(!theDefinition) |
---|
| 764 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 765 | ResSec->SetDefinition( theDefinition ); |
---|
| 766 | FstSec->SetDefinition( aNeutron ); |
---|
| 767 | SecSec->SetDefinition( aProton ); |
---|
| 768 | mF=mNeut; // First hadron 4-momentum |
---|
| 769 | mS=mProt; // Second hadron 4-momentum |
---|
| 770 | break; |
---|
| 771 | case 9: |
---|
| 772 | rZ-=2; |
---|
| 773 | rA-=2; |
---|
| 774 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 775 | if(!theDefinition) |
---|
| 776 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 777 | ResSec->SetDefinition( theDefinition ); |
---|
| 778 | FstSec->SetDefinition( aProton ); |
---|
| 779 | SecSec->SetDefinition( aProton ); |
---|
| 780 | mF=mProt; // First hadron 4-momentum |
---|
| 781 | mS=mProt; // Second hadron 4-momentum |
---|
| 782 | break; |
---|
| 783 | case 10: |
---|
| 784 | rZ-=2; |
---|
| 785 | rA-=3; |
---|
| 786 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 787 | if(!theDefinition) |
---|
| 788 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 789 | ResSec->SetDefinition( theDefinition ); |
---|
| 790 | FstSec->SetDefinition( aProton ); |
---|
| 791 | SecSec->SetDefinition( aDeuteron ); |
---|
| 792 | mF=mProt; // First hadron 4-momentum |
---|
| 793 | mS=mDeut; // Second hadron 4-momentum |
---|
| 794 | break; |
---|
| 795 | case 11: |
---|
| 796 | rZ-=1; |
---|
| 797 | rA-=3; |
---|
| 798 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 799 | if(!theDefinition) |
---|
| 800 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 801 | ResSec->SetDefinition( theDefinition ); |
---|
| 802 | FstSec->SetDefinition( aNeutron ); |
---|
| 803 | SecSec->SetDefinition( aDeuteron ); |
---|
| 804 | mF=mNeut; // First hadron 4-momentum |
---|
| 805 | mS=mDeut; // Second hadron 4-momentum |
---|
| 806 | break; |
---|
| 807 | case 12: |
---|
| 808 | rZ-=2; |
---|
| 809 | rA-=4; |
---|
| 810 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 811 | if(!theDefinition) |
---|
| 812 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 813 | ResSec->SetDefinition( theDefinition ); |
---|
| 814 | FstSec->SetDefinition( aDeuteron ); |
---|
| 815 | SecSec->SetDefinition( aDeuteron ); |
---|
| 816 | mF=mDeut; // First hadron 4-momentum |
---|
| 817 | mS=mDeut; // Second hadron 4-momentum |
---|
| 818 | break; |
---|
| 819 | case 13: |
---|
| 820 | rZ-=2; |
---|
| 821 | rA-=4; |
---|
| 822 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 823 | if(!theDefinition) |
---|
| 824 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 825 | ResSec->SetDefinition( theDefinition ); |
---|
| 826 | FstSec->SetDefinition( aProton ); |
---|
| 827 | SecSec->SetDefinition( aTriton ); |
---|
| 828 | mF=mProt; // First hadron 4-momentum |
---|
| 829 | mS=mTrit; // Second hadron 4-momentum |
---|
| 830 | break; |
---|
| 831 | case 14: |
---|
| 832 | rZ-=1; |
---|
| 833 | rA-=4; |
---|
| 834 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 835 | if(!theDefinition) |
---|
| 836 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 837 | ResSec->SetDefinition( theDefinition ); |
---|
| 838 | FstSec->SetDefinition( aNeutron ); |
---|
| 839 | SecSec->SetDefinition( aTriton ); |
---|
| 840 | mF=mNeut; // First hadron 4-momentum |
---|
| 841 | mS=mTrit; // Second hadron 4-momentum |
---|
| 842 | break; |
---|
| 843 | case 15: |
---|
| 844 | rZ-=3; |
---|
| 845 | rA-=4; |
---|
| 846 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 847 | if(!theDefinition) |
---|
| 848 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 849 | ResSec->SetDefinition( theDefinition ); |
---|
| 850 | FstSec->SetDefinition( aProton); |
---|
| 851 | SecSec->SetDefinition( aHe3 ); |
---|
| 852 | mF=mProt; // First hadron 4-momentum |
---|
| 853 | mS=mHel3; // Second hadron 4-momentum |
---|
| 854 | break; |
---|
| 855 | case 16: |
---|
| 856 | rZ-=2; |
---|
| 857 | rA-=4; |
---|
| 858 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 859 | if(!theDefinition) |
---|
| 860 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 861 | ResSec->SetDefinition( theDefinition ); |
---|
| 862 | FstSec->SetDefinition( aNeutron ); |
---|
| 863 | SecSec->SetDefinition( aHe3 ); |
---|
| 864 | mF=mNeut; // First hadron 4-momentum |
---|
| 865 | mS=mHel3; // Second hadron 4-momentum |
---|
| 866 | break; |
---|
| 867 | case 17: |
---|
| 868 | rZ-=3; |
---|
| 869 | rA-=5; |
---|
| 870 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 871 | if(!theDefinition) |
---|
| 872 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 873 | ResSec->SetDefinition( theDefinition ); |
---|
| 874 | FstSec->SetDefinition( aProton ); |
---|
| 875 | SecSec->SetDefinition( anAlpha ); |
---|
| 876 | mF=mProt; // First hadron 4-momentum |
---|
| 877 | mS=mAlph; // Second hadron 4-momentum |
---|
| 878 | break; |
---|
| 879 | case 18: |
---|
| 880 | rZ-=2; |
---|
| 881 | rA-=5; |
---|
| 882 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 883 | if(!theDefinition) |
---|
| 884 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 885 | ResSec->SetDefinition( theDefinition ); |
---|
| 886 | FstSec->SetDefinition( aNeutron ); |
---|
| 887 | SecSec->SetDefinition( anAlpha ); |
---|
| 888 | mF=mNeut; // First hadron 4-momentum |
---|
| 889 | mS=mAlph; // Second hadron 4-momentum |
---|
| 890 | break; |
---|
| 891 | case 19: |
---|
| 892 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 893 | if(!theDefinition) |
---|
| 894 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 895 | ResSec->SetDefinition( theDefinition ); |
---|
| 896 | FstSec->SetDefinition( aLambda ); |
---|
| 897 | SecSec->SetDefinition( aGamma ); |
---|
| 898 | mF=mLamb; // First hadron 4-momentum |
---|
| 899 | break; |
---|
| 900 | case 20: |
---|
| 901 | rZ-=1; |
---|
| 902 | rA-=1; |
---|
| 903 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 904 | if(!theDefinition) |
---|
| 905 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 906 | ResSec->SetDefinition( theDefinition ); |
---|
| 907 | FstSec->SetDefinition( aProton ); |
---|
| 908 | SecSec->SetDefinition( aLambda ); |
---|
| 909 | mF=mProt; // First hadron 4-momentum |
---|
| 910 | mS=mLamb; // Second hadron 4-momentum |
---|
| 911 | break; |
---|
| 912 | case 21: |
---|
| 913 | rA-=1; |
---|
| 914 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 915 | if(!theDefinition) |
---|
| 916 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 917 | ResSec->SetDefinition( theDefinition ); |
---|
| 918 | FstSec->SetDefinition( aNeutron ); |
---|
| 919 | SecSec->SetDefinition( aLambda ); |
---|
| 920 | mF=mNeut; // First hadron 4-momentum |
---|
| 921 | mS=mLamb; // Second hadron 4-momentum |
---|
| 922 | break; |
---|
| 923 | case 22: |
---|
| 924 | rZ-=1; |
---|
| 925 | rA-=2; |
---|
| 926 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 927 | if(!theDefinition) |
---|
| 928 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 929 | ResSec->SetDefinition( theDefinition ); |
---|
| 930 | FstSec->SetDefinition( aDeuteron ); |
---|
| 931 | SecSec->SetDefinition( aLambda ); |
---|
| 932 | mF=mDeut; // First hadron 4-momentum |
---|
| 933 | mS=mLamb; // Second hadron 4-momentum |
---|
| 934 | break; |
---|
| 935 | case 23: |
---|
| 936 | rZ-=1; |
---|
| 937 | rA-=3; |
---|
| 938 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 939 | if(!theDefinition) |
---|
| 940 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 941 | ResSec->SetDefinition( theDefinition ); |
---|
| 942 | FstSec->SetDefinition( aTriton ); |
---|
| 943 | SecSec->SetDefinition( aLambda ); |
---|
| 944 | mF=mTrit; // First hadron 4-momentum |
---|
| 945 | mS=mLamb; // Second hadron 4-momentum |
---|
| 946 | break; |
---|
| 947 | case 24: |
---|
| 948 | rZ-=2; |
---|
| 949 | rA-=3; |
---|
| 950 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 951 | if(!theDefinition) |
---|
| 952 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 953 | ResSec->SetDefinition( theDefinition ); |
---|
| 954 | FstSec->SetDefinition( aHe3 ); |
---|
| 955 | SecSec->SetDefinition( aLambda ); |
---|
| 956 | mF=mHel3; // First hadron 4-momentum |
---|
| 957 | mS=mLamb; // Second hadron 4-momentum |
---|
| 958 | break; |
---|
| 959 | case 25: |
---|
| 960 | rZ-=2; |
---|
| 961 | rA-=4; |
---|
| 962 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(rZ,rA,rL,rZ); |
---|
| 963 | if(!theDefinition) |
---|
| 964 | G4cerr<<"-Warning-G4LE::PSDI: notDefined, Z="<<rZ<<", A="<<rA<<", L="<<rL<<G4endl; |
---|
| 965 | ResSec->SetDefinition( theDefinition ); |
---|
| 966 | FstSec->SetDefinition( anAlpha ); |
---|
| 967 | SecSec->SetDefinition( aLambda ); |
---|
| 968 | mF=mAlph; // First hadron 4-momentum |
---|
| 969 | mS=mLamb; // Second hadron 4-momentum |
---|
| 970 | break; |
---|
| 971 | } |
---|
| 972 | #ifdef debug |
---|
| 973 | G4cout<<"G4QLowEn::PSDI:F="<<mF<<",S="<<mS<<",com="<<complete<<",ev="<<evaporate<<G4endl; |
---|
| 974 | #endif |
---|
| 975 | G4LorentzVector fst4Mom(zeroMom,mF); // Prototype of the first hadron 4-momentum |
---|
| 976 | G4LorentzVector snd4Mom(zeroMom,mS); // Prototype of the second hadron 4-momentum |
---|
| 977 | G4LorentzVector dir4Mom=tot4M; // Prototype of the resN decay direction 4-momentum |
---|
| 978 | dir4Mom.setE(tot4M.e()/2.); // Get half energy and total 3-momentum |
---|
| 979 | // @@ Can be repeated to take into account the Coulomb Barrier |
---|
| 980 | if(!G4QHadron(tot4M).CopDecayIn3(fst4Mom,snd4Mom,res4Mom,dir4Mom,cosp)) |
---|
| 981 | { // |
---|
| 982 | G4cerr<<"**G4LowEnergy::PoStDoIt:i="<<index<<",tM="<<totM<<"->M1="<<res4Mom.m()<<"+M2=" |
---|
| 983 | <<fst4Mom.m()<<"+M3="<<snd4Mom.m()<<"=="<<res4Mom.m()+fst4Mom.m()+snd4Mom.m()<<G4endl; |
---|
| 984 | throw G4QException("G4QLowEnergy::PostStepDoIt: Can't decay the Compound"); |
---|
| 985 | } // |
---|
| 986 | #ifdef debug |
---|
| 987 | G4cout<<"G4QLowEn::PSDI:r4M="<<res4Mom<<",f4M="<<fst4Mom<<",s4M="<<snd4Mom<<G4endl; |
---|
| 988 | #endif |
---|
| 989 | G4Track* aNewTrack = 0; |
---|
| 990 | if(complete) |
---|
| 991 | { |
---|
| 992 | FstSec->Set4Momentum(fst4Mom); |
---|
| 993 | aNewTrack = new G4Track(FstSec, localtime, position ); |
---|
| 994 | aNewTrack->SetWeight(weight); // weighted |
---|
| 995 | aNewTrack->SetTouchableHandle(trTouchable); |
---|
| 996 | aParticleChange.AddSecondary( aNewTrack ); |
---|
| 997 | EnMomConservation-=fst4Mom; |
---|
| 998 | #ifdef debug |
---|
| 999 | G4cout<<"G4QLowEn::PSDI: ***Filled*** 1stH4M="<<fst4Mom<<G4endl; |
---|
| 1000 | #endif |
---|
| 1001 | if(complete>2) // Final solution |
---|
| 1002 | { |
---|
| 1003 | ResSec->Set4Momentum(res4Mom); |
---|
| 1004 | aNewTrack = new G4Track(ResSec, localtime, position ); |
---|
| 1005 | aNewTrack->SetWeight(weight); // weighted |
---|
| 1006 | aNewTrack->SetTouchableHandle(trTouchable); |
---|
| 1007 | aParticleChange.AddSecondary( aNewTrack ); |
---|
| 1008 | EnMomConservation-=res4Mom; |
---|
| 1009 | #ifdef debug |
---|
| 1010 | G4cout<<"G4QLowEn::PSDI: ***Filled*** rA4M="<<res4Mom<<",rZ="<<rZ<<",rA="<<rA<<",rL=" |
---|
| 1011 | <<rL<<G4endl; |
---|
| 1012 | #endif |
---|
| 1013 | SecSec->Set4Momentum(snd4Mom); |
---|
| 1014 | aNewTrack = new G4Track(SecSec, localtime, position ); |
---|
| 1015 | aNewTrack->SetWeight(weight); // weighted |
---|
| 1016 | aNewTrack->SetTouchableHandle(trTouchable); |
---|
| 1017 | aParticleChange.AddSecondary( aNewTrack ); |
---|
| 1018 | EnMomConservation-=snd4Mom; |
---|
| 1019 | #ifdef debug |
---|
| 1020 | G4cout<<"G4QLowEn::PSDI: ***Filled*** 2ndH4M="<<snd4Mom<<G4endl; |
---|
| 1021 | #endif |
---|
| 1022 | } |
---|
| 1023 | else res4Mom+=snd4Mom; |
---|
| 1024 | } |
---|
| 1025 | else res4Mom=tot4M; |
---|
| 1026 | if(complete<3) // Evaporation of the residual must be done |
---|
| 1027 | { |
---|
| 1028 | G4QHadron* rHadron = new G4QHadron(90000000+999*rZ+rA,res4Mom); // Input hadron-nucleus |
---|
| 1029 | G4QHadronVector* evaHV = new G4QHadronVector; // Output vector of hadrons (delete!) |
---|
| 1030 | Nuc.EvaporateNucleus(rHadron, evaHV); |
---|
| 1031 | G4int nOut=evaHV->size(); |
---|
| 1032 | for(G4int i=0; i<nOut; i++) |
---|
| 1033 | { |
---|
| 1034 | G4QHadron* curH = (*evaHV)[i]; |
---|
| 1035 | G4int hPDG=curH->GetPDGCode(); |
---|
| 1036 | G4LorentzVector h4Mom=curH->Get4Momentum(); |
---|
| 1037 | EnMomConservation-=h4Mom; |
---|
| 1038 | #ifdef debug |
---|
| 1039 | G4cout<<"G4QLowEn::PSDI: ***FillingCand#"<<i<<"*** evaH="<<hPDG<<h4Mom<<G4endl; |
---|
| 1040 | #endif |
---|
| 1041 | if (hPDG==90000001 || hPDG==2112) theDefinition = aNeutron; |
---|
| 1042 | else if(hPDG==90001000 || hPDG==2212) theDefinition = aProton; |
---|
| 1043 | else if(hPDG==91000000 || hPDG==3122) theDefinition = aLambda; |
---|
| 1044 | else if(hPDG== 22 ) theDefinition = aGamma; |
---|
| 1045 | else |
---|
| 1046 | { |
---|
| 1047 | G4int hZ=curH->GetCharge(); |
---|
| 1048 | G4int hA=curH->GetBaryonNumber(); |
---|
| 1049 | G4int hS=curH->GetStrangeness(); |
---|
| 1050 | theDefinition = G4ParticleTable::GetParticleTable()->FindIon(hZ,hA,hS,0); // ion |
---|
| 1051 | } |
---|
| 1052 | if(theDefinition) |
---|
| 1053 | { |
---|
| 1054 | G4DynamicParticle* theEQH = new G4DynamicParticle(theDefinition,h4Mom); |
---|
| 1055 | G4Track* evaQH = new G4Track(theEQH, localtime, position ); |
---|
| 1056 | evaQH->SetWeight(weight); // weighted |
---|
| 1057 | evaQH->SetTouchableHandle(trTouchable); |
---|
| 1058 | aParticleChange.AddSecondary( evaQH ); |
---|
| 1059 | } |
---|
| 1060 | else G4cerr<<"-Warning-G4QLowEnergy::PostStepDoIt: Bad secondary PDG="<<hPDG<<G4endl; |
---|
| 1061 | } |
---|
| 1062 | } |
---|
| 1063 | // algorithm implementation --- STOPS HERE |
---|
| 1064 | #ifdef debug |
---|
| 1065 | G4cout<<"G4QLowEnergy::PostStepDoIt:*** PostStepDoIt is done ***"<<G4endl; |
---|
| 1066 | #endif |
---|
| 1067 | return G4VDiscreteProcess::PostStepDoIt(track, step); |
---|
| 1068 | } |
---|
| 1069 | |
---|
| 1070 | G4double G4QLowEnergy::CalculateXS(G4double p, G4int Z, G4int N, G4int PDG) |
---|
| 1071 | { |
---|
| 1072 | static G4bool first=true; |
---|
| 1073 | static G4VQCrossSection* CSmanager; |
---|
| 1074 | if(first) // Connection with a singletone |
---|
| 1075 | { |
---|
| 1076 | CSmanager=G4QIonIonCrossSection::GetPointer(); |
---|
| 1077 | first=false; |
---|
| 1078 | } |
---|
| 1079 | #ifdef debug |
---|
| 1080 | G4cout<<"G4QLowE::CXS: *DONE* p="<<p<<",Z="<<Z<<",N="<<N<<",PDG="<<PDG<<G4endl; |
---|
| 1081 | #endif |
---|
| 1082 | return CSmanager->GetCrossSection(true, p, Z, N, PDG); |
---|
| 1083 | } |
---|