[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 | // $Id: G4QCoherentChargeExchange.cc,v 1.5 2007/10/02 10:00:37 mkossov Exp $ |
---|
| 27 | // GEANT4 tag $Name: $ |
---|
| 28 | // |
---|
| 29 | // ---------------- G4QCoherentChargeExchange class ----------------- |
---|
| 30 | // by Mikhail Kossov, December 2003. |
---|
| 31 | // G4QCoherentChargeExchange class of the CHIPS Simulation Branch in GEANT4 |
---|
| 32 | // --------------------------------------------------------------- |
---|
| 33 | // **************************************************************************************** |
---|
| 34 | // ********** This CLASS is temporary moved from the photolepton_hadron directory ********* |
---|
| 35 | // **************************************************************************************** |
---|
| 36 | |
---|
| 37 | //#define debug |
---|
| 38 | //#define pdebug |
---|
| 39 | //#define tdebug |
---|
| 40 | //#define nandebug |
---|
| 41 | //#define ppdebug |
---|
| 42 | |
---|
| 43 | #include "G4QCoherentChargeExchange.hh" |
---|
| 44 | |
---|
| 45 | // Initialization of static vectors |
---|
| 46 | G4int G4QCoherentChargeExchange::nPartCWorld=152; // #of particles initialized in CHIPS |
---|
| 47 | std::vector<G4int> G4QCoherentChargeExchange::ElementZ; // Z of element(i) in theLastCalc |
---|
| 48 | std::vector<G4double> G4QCoherentChargeExchange::ElProbInMat; // SumProbOfElem in Material |
---|
| 49 | std::vector<std::vector<G4int>*> G4QCoherentChargeExchange::ElIsoN;// N of isotope(j), E(i) |
---|
| 50 | std::vector<std::vector<G4double>*>G4QCoherentChargeExchange::IsoProbInEl;//SumProbIsotE(i) |
---|
| 51 | |
---|
| 52 | // Constructor |
---|
| 53 | G4QCoherentChargeExchange::G4QCoherentChargeExchange(const G4String& processName) |
---|
| 54 | : G4VDiscreteProcess(processName) |
---|
| 55 | { |
---|
| 56 | #ifdef debug |
---|
| 57 | G4cout<<"G4QCohChargeEx::Constructor is called processName="<<processName<<G4endl; |
---|
| 58 | #endif |
---|
| 59 | if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl; |
---|
| 60 | |
---|
| 61 | //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max) |
---|
| 62 | } |
---|
| 63 | |
---|
| 64 | // Destructor |
---|
| 65 | G4QCoherentChargeExchange::~G4QCoherentChargeExchange() {} |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | G4LorentzVector G4QCoherentChargeExchange::GetEnegryMomentumConservation() |
---|
| 69 | {return EnMomConservation;} |
---|
| 70 | |
---|
| 71 | G4int G4QCoherentChargeExchange::GetNumberOfNeutronsInTarget() {return nOfNeutrons;} |
---|
| 72 | |
---|
| 73 | // output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)), |
---|
| 74 | // where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3) |
---|
| 75 | // ********** All CHIPS cross sections are calculated in the surface units ************ |
---|
| 76 | G4double G4QCoherentChargeExchange::GetMeanFreePath(const G4Track& aTrack,G4double Q, |
---|
| 77 | G4ForceCondition* Fc) |
---|
| 78 | { |
---|
| 79 | *Fc = NotForced; |
---|
| 80 | const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle(); |
---|
| 81 | G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition(); |
---|
| 82 | if( !IsApplicable(*incidentParticleDefinition)) |
---|
| 83 | G4cout<<"*W*G4QCohChargeEx::GetMeanFreePath called for notImplementedParticle"<<G4endl; |
---|
| 84 | // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material |
---|
| 85 | G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle |
---|
| 86 | #ifdef debug |
---|
| 87 | G4double KinEn = incidentParticle->GetKineticEnergy(); |
---|
| 88 | G4cout<<"G4QCohChEx::GetMeanFreePath: kinE="<<KinEn<<",Mom="<<Momentum<<G4endl; // Result |
---|
| 89 | #endif |
---|
| 90 | const G4Material* material = aTrack.GetMaterial(); // Get the current material |
---|
| 91 | const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume(); |
---|
| 92 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
| 93 | G4int nE=material->GetNumberOfElements(); |
---|
| 94 | #ifdef debug |
---|
| 95 | G4cout<<"G4QCohChargeExchange::GetMeanFreePath:"<<nE<<" Elem's in theMaterial"<<G4endl; |
---|
| 96 | #endif |
---|
| 97 | G4int pPDG=0; |
---|
| 98 | |
---|
| 99 | if (incidentParticleDefinition == G4Proton::Proton() ) pPDG=2212; |
---|
| 100 | else if(incidentParticleDefinition == G4Neutron::Neutron()) pPDG=2112; |
---|
| 101 | else G4cout<<"G4QCohChargeEx::GetMeanFreePath: only nA & pA are implemented"<<G4endl; |
---|
| 102 | |
---|
| 103 | G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton |
---|
| 104 | G4double sigma=0.; // Sums over elements for the material |
---|
| 105 | G4int IPIE=IsoProbInEl.size(); // How many old elements? |
---|
| 106 | if(IPIE) for(G4int ip=0; ip<IPIE; ++ip) // Clean up the SumProb's of Isotopes (SPI) |
---|
| 107 | { |
---|
| 108 | std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector |
---|
| 109 | SPI->clear(); |
---|
| 110 | delete SPI; |
---|
| 111 | std::vector<G4int>* IsN=ElIsoN[ip]; // Pointer to the N vector |
---|
| 112 | IsN->clear(); |
---|
| 113 | delete IsN; |
---|
| 114 | } |
---|
| 115 | ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE) |
---|
| 116 | ElementZ.clear(); // Clear the body vector for Z of Elements |
---|
| 117 | IsoProbInEl.clear(); // Clear the body vector for SPI |
---|
| 118 | ElIsoN.clear(); // Clear the body vector for N of Isotopes |
---|
| 119 | for(G4int i=0; i<nE; ++i) |
---|
| 120 | { |
---|
| 121 | G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element |
---|
| 122 | G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element |
---|
| 123 | ElementZ.push_back(Z); // Remember Z of the Element |
---|
| 124 | G4int isoSize=0; // The default for the isoVectorLength is 0 |
---|
| 125 | G4int indEl=0; // Index of non-natural element or 0(default) |
---|
| 126 | G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect |
---|
| 127 | if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector |
---|
| 128 | #ifdef debug |
---|
| 129 | G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath:isovectorLength="<<isoSize<<G4endl; |
---|
| 130 | #endif |
---|
| 131 | if(isoSize) // The Element has non-trivial abundance set |
---|
| 132 | { |
---|
| 133 | indEl=pElement->GetIndex()+1; // Index of the non-trivial element is an order |
---|
| 134 | #ifdef debug |
---|
| 135 | G4cout<<"G4QCCX::GetMFP: iE="<<indEl<<", def="<<Isotopes->IsDefined(Z,indEl)<<G4endl; |
---|
| 136 | #endif |
---|
| 137 | if(!Isotopes->IsDefined(Z,indEl)) // This index is not defined for this Z: define |
---|
| 138 | { |
---|
| 139 | std::vector<std::pair<G4int,G4double>*>* newAbund = |
---|
| 140 | new std::vector<std::pair<G4int,G4double>*>; |
---|
| 141 | G4double* abuVector=pElement->GetRelativeAbundanceVector(); |
---|
| 142 | for(G4int j=0; j<isoSize; j++) // Calculation of abundance vector for isotopes |
---|
| 143 | { |
---|
| 144 | G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z ! |
---|
| 145 | if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QCohChX::GetMeanFreePath: Z=" |
---|
| 146 | <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl; |
---|
| 147 | G4double abund=abuVector[j]; |
---|
| 148 | std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund); |
---|
| 149 | #ifdef debug |
---|
| 150 | G4cout<<"G4QCohChEx::GetMeanFreePath:pair#="<<j<<",N="<<N<<",ab="<<abund<<G4endl; |
---|
| 151 | #endif |
---|
| 152 | newAbund->push_back(pr); |
---|
| 153 | } |
---|
| 154 | #ifdef debug |
---|
| 155 | G4cout<<"G4QCohChEx::GetMeanFreePath: pairVectorLength="<<newAbund->size()<<G4endl; |
---|
| 156 | #endif |
---|
| 157 | indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd |
---|
| 158 | for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k]; // Cleaning temporary |
---|
| 159 | delete newAbund; // Was "new" in the beginning of the name space |
---|
| 160 | } |
---|
| 161 | } |
---|
| 162 | std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer |
---|
| 163 | std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector |
---|
| 164 | IsoProbInEl.push_back(SPI); |
---|
| 165 | std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector |
---|
| 166 | ElIsoN.push_back(IsN); |
---|
| 167 | G4int nIs=cs->size(); // A#Of Isotopes in the Element |
---|
| 168 | #ifdef debug |
---|
| 169 | G4cout<<"G4QCohChargEx::GMFP:=***=>,#isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl; |
---|
| 170 | #endif |
---|
| 171 | G4double susi=0.; // sum of CS over isotopes |
---|
| 172 | if(nIs) for(G4int j=0; j<nIs; j++) // Calculate CS for eachIsotope of El |
---|
| 173 | { |
---|
| 174 | std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice |
---|
| 175 | G4int N=curIs->first; // #of Neuterons in the isotope j of El i |
---|
| 176 | IsN->push_back(N); // Remember Min N for the Element |
---|
| 177 | #ifdef debug |
---|
| 178 | G4cout<<"G4QCCX::GMFP:true, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl; |
---|
| 179 | #endif |
---|
| 180 | G4bool ccsf=true; |
---|
| 181 | if(Q==-27.) ccsf=false; |
---|
| 182 | #ifdef debug |
---|
| 183 | G4cout<<"G4QCoherentChargeExchange::GMFP: GetCS #1 j="<<j<<G4endl; |
---|
| 184 | #endif |
---|
| 185 | G4double CSI=CalculateXSt(ccsf, true, Momentum, Z, N, pPDG);// CS(j,i) for theIsotope |
---|
| 186 | |
---|
| 187 | #ifdef debug |
---|
| 188 | G4cout<<"G4QCohChEx::GetMeanFreePath:jI="<<j<<",Zt="<<Z<<",Nt="<<N<<",Mom="<<Momentum |
---|
| 189 | <<", XSec="<<CSI/millibarn<<G4endl; |
---|
| 190 | #endif |
---|
| 191 | curIs->second = CSI; |
---|
| 192 | susi+=CSI; // Make a sum per isotopes |
---|
| 193 | SPI->push_back(susi); // Remember summed cross-section |
---|
| 194 | } // End of temporary initialization of the cross sections in the G4QIsotope singeltone |
---|
| 195 | sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV) |
---|
| 196 | #ifdef debug |
---|
| 197 | G4cout<<"G4QCohChEx::GMFP:<S>="<<Isotopes->GetMeanCrossSection(Z,indEl)<<",AddToSigma=" |
---|
| 198 | <<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl; |
---|
| 199 | #endif |
---|
| 200 | ElProbInMat.push_back(sigma); |
---|
| 201 | } // End of LOOP over Elements |
---|
| 202 | // Check that cross section is not zero and return the mean free path |
---|
| 203 | #ifdef debug |
---|
| 204 | G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl; |
---|
| 205 | #endif |
---|
| 206 | if(sigma > 0.) return 1./sigma; // Mean path [distance] |
---|
| 207 | return DBL_MAX; |
---|
| 208 | } |
---|
| 209 | |
---|
| 210 | G4bool G4QCoherentChargeExchange::IsApplicable(const G4ParticleDefinition& particle) |
---|
| 211 | { |
---|
| 212 | if (particle == *( G4Proton::Proton() )) return true; |
---|
| 213 | else if (particle == *( G4Neutron::Neutron() )) return true; |
---|
| 214 | //else if (particle == *( G4MuonMinus::MuonMinus() )) return true; |
---|
| 215 | //else if (particle == *( G4TauPlus::TauPlus() )) return true; |
---|
| 216 | //else if (particle == *( G4TauMinus::TauMinus() )) return true; |
---|
| 217 | //else if (particle == *( G4Electron::Electron() )) return true; |
---|
| 218 | //else if (particle == *( G4Positron::Positron() )) return true; |
---|
| 219 | //else if (particle == *( G4Gamma::Gamma() )) return true; |
---|
| 220 | //else if (particle == *( G4MuonPlus::MuonPlus() )) return true; |
---|
| 221 | //else if (particle == *(G4AntiNeutrinoMu::AntiNeutrinoMu())) return true; |
---|
| 222 | //else if (particle == *( G4NeutrinoMu::NeutrinoMu() )) return true; |
---|
| 223 | //else if (particle == *( G4PionMinus::PionMinus() )) return true; |
---|
| 224 | //else if (particle == *( G4PionPlus::PionPlus() )) return true; |
---|
| 225 | //else if (particle == *( G4KaonPlus::KaonPlus() )) return true; |
---|
| 226 | //else if (particle == *( G4KaonMinus::KaonMinus() )) return true; |
---|
| 227 | //else if (particle == *( G4KaonZeroLong::KaonZeroLong() )) return true; |
---|
| 228 | //else if (particle == *( G4KaonZeroShort::KaonZeroShort() )) return true; |
---|
| 229 | //else if (particle == *( G4Lambda::Lambda() )) return true; |
---|
| 230 | //else if (particle == *( G4SigmaPlus::SigmaPlus() )) return true; |
---|
| 231 | //else if (particle == *( G4SigmaMinus::SigmaMinus() )) return true; |
---|
| 232 | //else if (particle == *( G4SigmaZero::SigmaZero() )) return true; |
---|
| 233 | //else if (particle == *( G4XiMinus::XiMinus() )) return true; |
---|
| 234 | //else if (particle == *( G4XiZero::XiZero() )) return true; |
---|
| 235 | //else if (particle == *( G4OmegaMinus::OmegaMinus() )) return true; |
---|
| 236 | //else if (particle == *( G4AntiNeutron::AntiNeutron() )) return true; |
---|
| 237 | //else if (particle == *( G4AntiProton::AntiProton() )) return true; |
---|
| 238 | #ifdef debug |
---|
| 239 | G4cout<<"***>>G4QCoherChargeExch::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl; |
---|
| 240 | #endif |
---|
| 241 | return false; |
---|
| 242 | } |
---|
| 243 | |
---|
| 244 | G4VParticleChange* G4QCoherentChargeExchange::PostStepDoIt(const G4Track& track, |
---|
| 245 | const G4Step& step) |
---|
| 246 | { |
---|
| 247 | static const G4double mProt= G4QPDGCode(2212).GetMass(); // CHIPS pMass in MeV |
---|
| 248 | static const G4double mNeut= G4QPDGCode(2212).GetMass(); // CHIPS pMass in MeV |
---|
| 249 | // |
---|
| 250 | //------------------------------------------------------------------------------------- |
---|
| 251 | static G4bool CWinit = true; // CHIPS Warld needs to be initted |
---|
| 252 | if(CWinit) |
---|
| 253 | { |
---|
| 254 | CWinit=false; |
---|
| 255 | G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max) |
---|
| 256 | } |
---|
| 257 | //------------------------------------------------------------------------------------- |
---|
| 258 | const G4DynamicParticle* projHadron = track.GetDynamicParticle(); |
---|
| 259 | const G4ParticleDefinition* particle=projHadron->GetDefinition(); |
---|
| 260 | #ifdef debug |
---|
| 261 | G4cout<<"G4QCohChargeExchange::PostStepDoIt: Before the GetMeanFreePath is called In4M=" |
---|
| 262 | <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type=" |
---|
| 263 | <<particle->GetParticleType()<<", Subtp="<<particle->GetParticleSubType()<<G4endl; |
---|
| 264 | #endif |
---|
| 265 | G4ForceCondition cond=NotForced; |
---|
| 266 | GetMeanFreePath(track, -27., &cond); // @@ ?? jus to update parameters? |
---|
| 267 | #ifdef debug |
---|
| 268 | G4cout<<"G4QCohChargeExchange::PostStepDoIt: After GetMeanFreePath is called"<<G4endl; |
---|
| 269 | #endif |
---|
| 270 | G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV! |
---|
| 271 | G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV |
---|
| 272 | G4double Momentum = proj4M.rho(); // @@ Just for the test purposes |
---|
| 273 | if(std::fabs(Momentum-momentum)>.000001) |
---|
| 274 | G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt:P(IU)="<<Momentum<<"#"<<momentum<<G4endl; |
---|
| 275 | #ifdef pdebug |
---|
| 276 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum |
---|
| 277 | <<",pM="<<pM<<",proj4M="<<proj4M<<proj4M.m()<<G4endl; |
---|
| 278 | #endif |
---|
| 279 | if (!IsApplicable(*particle)) // Check applicability |
---|
| 280 | { |
---|
| 281 | G4cerr<<"G4QCoherentChargeExchange::PostStepDoIt: Only NA is implemented."<<G4endl; |
---|
| 282 | return 0; |
---|
| 283 | } |
---|
| 284 | const G4Material* material = track.GetMaterial(); // Get the current material |
---|
| 285 | G4int Z=0; |
---|
| 286 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
| 287 | G4int nE=material->GetNumberOfElements(); |
---|
| 288 | #ifdef debug |
---|
| 289 | G4cout<<"G4QCohChargeExchange::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl; |
---|
| 290 | #endif |
---|
| 291 | G4int projPDG=0; // PDG Code prototype for the captured hadron |
---|
| 292 | // Not all these particles are implemented yet (see Is Applicable) |
---|
| 293 | if (particle == G4Proton::Proton() ) projPDG= 2212; |
---|
| 294 | else if (particle == G4Neutron::Neutron() ) projPDG= 2112; |
---|
| 295 | //else if (particle == G4PionMinus::PionMinus() ) projPDG= -211; |
---|
| 296 | //else if (particle == G4PionPlus::PionPlus() ) projPDG= 211; |
---|
| 297 | //else if (particle == G4KaonPlus::KaonPlus() ) projPDG= 2112; |
---|
| 298 | //else if (particle == G4KaonMinus::KaonMinus() ) projPDG= -321; |
---|
| 299 | //else if (particle == G4KaonZeroLong::KaonZeroLong() ) projPDG= 130; |
---|
| 300 | //else if (particle == G4KaonZeroShort::KaonZeroShort() ) projPDG= 310; |
---|
| 301 | //else if (particle == G4MuonPlus::MuonPlus() ) projPDG= -13; |
---|
| 302 | //else if (particle == G4MuonMinus::MuonMinus() ) projPDG= 13; |
---|
| 303 | //else if (particle == G4NeutrinoMu::NeutrinoMu() ) projPDG= 14; |
---|
| 304 | //else if (particle == G4AntiNeutrinoMu::AntiNeutrinoMu() ) projPDG= -14; |
---|
| 305 | //else if (particle == G4Electron::Electron() ) projPDG= 11; |
---|
| 306 | //else if (particle == G4Positron::Positron() ) projPDG= -11; |
---|
| 307 | //else if (particle == G4NeutrinoE::NeutrinoE() ) projPDG= 12; |
---|
| 308 | //else if (particle == G4AntiNeutrinoE::AntiNeutrinoE() ) projPDG= -12; |
---|
| 309 | //else if (particle == G4Gamma::Gamma() ) projPDG= 22; |
---|
| 310 | //else if (particle == G4TauPlus::TauPlus() ) projPDG= -15; |
---|
| 311 | //else if (particle == G4TauMinus::TauMinus() ) projPDG= 15; |
---|
| 312 | //else if (particle == G4NeutrinoTau::NeutrinoTau() ) projPDG= 16; |
---|
| 313 | //else if (particle == G4AntiNeutrinoTau::AntiNeutrinoTau()) projPDG= -16; |
---|
| 314 | //else if (particle == G4Lambda::Lambda() ) projPDG= 3122; |
---|
| 315 | //else if (particle == G4SigmaPlus::SigmaPlus() ) projPDG= 3222; |
---|
| 316 | //else if (particle == G4SigmaMinus::SigmaMinus() ) projPDG= 3112; |
---|
| 317 | //else if (particle == G4SigmaZero::SigmaZero() ) projPDG= 3212; |
---|
| 318 | //else if (particle == G4XiMinus::XiMinus() ) projPDG= 3312; |
---|
| 319 | //else if (particle == G4XiZero::XiZero() ) projPDG= 3322; |
---|
| 320 | //else if (particle == G4OmegaMinus::OmegaMinus() ) projPDG= 3334; |
---|
| 321 | //else if (particle == G4AntiNeutron::AntiNeutron() ) projPDG=-2112; |
---|
| 322 | //else if (particle == G4AntiProton::AntiProton() ) projPDG=-2212; |
---|
| 323 | #ifdef debug |
---|
| 324 | G4int prPDG=particle->GetPDGEncoding(); |
---|
| 325 | G4cout<<"G4QCohChrgExchange::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl; |
---|
| 326 | #endif |
---|
| 327 | if(!projPDG) |
---|
| 328 | { |
---|
| 329 | G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:UndefinedProjHadron"<<G4endl; |
---|
| 330 | return 0; |
---|
| 331 | } |
---|
| 332 | //G4double pM2=proj4M.m2(); // in MeV^2 |
---|
| 333 | //G4double pM=std::sqrt(pM2); // in MeV |
---|
| 334 | G4double pM=mNeut; |
---|
| 335 | G4int fPDG=2112; |
---|
| 336 | if(projPDG==2112) |
---|
| 337 | { |
---|
| 338 | pM=mProt; |
---|
| 339 | fPDG=2212; |
---|
| 340 | } |
---|
| 341 | G4double pM2=pM*pM; |
---|
| 342 | // Element treatment |
---|
| 343 | G4int EPIM=ElProbInMat.size(); |
---|
| 344 | #ifdef debug |
---|
| 345 | G4cout<<"G4QCohChEx::PostStDoIt:m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl; |
---|
| 346 | #endif |
---|
| 347 | G4int i=0; |
---|
| 348 | if(EPIM>1) |
---|
| 349 | { |
---|
| 350 | G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand(); |
---|
| 351 | for(i=0; i<nE; ++i) |
---|
| 352 | { |
---|
| 353 | #ifdef debug |
---|
| 354 | G4cout<<"G4QCohChEx::PostStepDoIt:EPM["<<i<<"]="<<ElProbInMat[i]<<",r="<<rnd<<G4endl; |
---|
| 355 | #endif |
---|
| 356 | if (rnd<ElProbInMat[i]) break; |
---|
| 357 | } |
---|
| 358 | if(i>=nE) i=nE-1; // Top limit for the Element |
---|
| 359 | } |
---|
| 360 | G4Element* pElement=(*theElementVector)[i]; |
---|
| 361 | Z=static_cast<G4int>(pElement->GetZ()); |
---|
| 362 | #ifdef debug |
---|
| 363 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl; |
---|
| 364 | #endif |
---|
| 365 | if(Z<=0) |
---|
| 366 | { |
---|
| 367 | G4cerr<<"-Warning-G4QCoherentChargeExchange::PostStepDoIt: Element with Z="<<Z<<G4endl; |
---|
| 368 | if(Z<0) return 0; |
---|
| 369 | } |
---|
| 370 | std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes |
---|
| 371 | std::vector<G4int>* IsN = ElIsoN[i]; // Vector of "#of neutrons" in the isotope El[i] |
---|
| 372 | G4int nofIsot=SPI->size(); // #of isotopes in the element i |
---|
| 373 | #ifdef debug |
---|
| 374 | G4cout<<"G4QCohChargeExchange::PosStDoIt:nI="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl; |
---|
| 375 | #endif |
---|
| 376 | G4int j=0; |
---|
| 377 | if(nofIsot>1) |
---|
| 378 | { |
---|
| 379 | G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element |
---|
| 380 | for(j=0; j<nofIsot; ++j) |
---|
| 381 | { |
---|
| 382 | #ifdef debug |
---|
| 383 | G4cout<<"G4QCohChargEx::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl; |
---|
| 384 | #endif |
---|
| 385 | if(rndI < (*SPI)[j]) break; |
---|
| 386 | } |
---|
| 387 | if(j>=nofIsot) j=nofIsot-1; // Top limit for the isotope |
---|
| 388 | } |
---|
| 389 | G4int N =(*IsN)[j]; ; // Randomized number of neutrons |
---|
| 390 | #ifdef debug |
---|
| 391 | G4cout<<"G4QCohChargeEx::PostStepDoIt: j="<<i<<", N(isotope)="<<N<<", MeV="<<MeV<<G4endl; |
---|
| 392 | #endif |
---|
| 393 | if(N<0) |
---|
| 394 | { |
---|
| 395 | G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt: Isotope with Z="<<Z<<", 0>N="<<N<<G4endl; |
---|
| 396 | return 0; |
---|
| 397 | } |
---|
| 398 | nOfNeutrons=N; // Remember it for the energy-momentum check |
---|
| 399 | #ifdef debug |
---|
| 400 | G4cout<<"G4QCohChargeExchange::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl; |
---|
| 401 | #endif |
---|
| 402 | if(N<0) |
---|
| 403 | { |
---|
| 404 | G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:Element with N="<<N<< G4endl; |
---|
| 405 | return 0; |
---|
| 406 | } |
---|
| 407 | aParticleChange.Initialize(track); |
---|
| 408 | #ifdef debug |
---|
| 409 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: track is initialized"<<G4endl; |
---|
| 410 | #endif |
---|
| 411 | G4double weight = track.GetWeight(); |
---|
| 412 | G4double localtime = track.GetGlobalTime(); |
---|
| 413 | G4ThreeVector position = track.GetPosition(); |
---|
| 414 | #ifdef debug |
---|
| 415 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: before Touchable extraction"<<G4endl; |
---|
| 416 | #endif |
---|
| 417 | G4TouchableHandle trTouchable = track.GetTouchableHandle(); |
---|
| 418 | #ifdef debug |
---|
| 419 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Touchable is extracted"<<G4endl; |
---|
| 420 | #endif |
---|
| 421 | // |
---|
| 422 | G4int targPDG=90000000+Z*1000+N; // CHIPS PDG Code of the target nucleus |
---|
| 423 | if(projPDG==2212) targPDG+=999; // convert to final nucleus code |
---|
| 424 | else if(projPDG==2112) targPDG-=999; |
---|
| 425 | G4QPDGCode targQPDG(targPDG); // @@ use G4Ion and get rid of CHIPS World |
---|
| 426 | G4double tM=targQPDG.GetMass(); // CHIPS final nucleus mass in MeV |
---|
| 427 | G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?) |
---|
| 428 | G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector |
---|
| 429 | G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction |
---|
| 430 | #ifdef debug |
---|
| 431 | G4cout<<"G4QCohChrgEx::PostStepDoIt: tM="<<tM<<", p4M="<<proj4M<<", t4M="<<tot4M<<G4endl; |
---|
| 432 | #endif |
---|
| 433 | EnMomConservation=tot4M; // Total 4-mom of reaction for E/M conservation |
---|
| 434 | // @@ Probably this is not necessary any more |
---|
| 435 | #ifdef debug |
---|
| 436 | G4cout<<"G4QCCX::PSDI:false, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<projPDG<<G4endl; |
---|
| 437 | #endif |
---|
| 438 | G4double xSec=CalculateXSt(false, true, Momentum, Z, N, projPDG); // Recalc. CrossSection |
---|
| 439 | #ifdef debug |
---|
| 440 | G4cout<<"G4QCoChEx::PSDI:PDG="<<projPDG<<",P="<<Momentum<<",CS="<<xSec/millibarn<<G4endl; |
---|
| 441 | #endif |
---|
| 442 | #ifdef nandebug |
---|
| 443 | if(xSec>0. || xSec<0. || xSec==0); |
---|
| 444 | else G4cout<<"*Warning*G4QCohChargeExchange::PSDI: xSec="<<xSec/millibarn<<G4endl; |
---|
| 445 | #endif |
---|
| 446 | // @@ check a possibility to separate p, n, or alpha (!) |
---|
| 447 | if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing |
---|
| 448 | { |
---|
| 449 | #ifdef pdebug |
---|
| 450 | G4cerr<<"*Warning*G4QCoherentChargeExchange::PSDoIt:*Zero cross-section* PDG="<<projPDG |
---|
| 451 | <<",tPDG="<<targPDG<<",P="<<Momentum<<G4endl; |
---|
| 452 | #endif |
---|
| 453 | //Do Nothing Action insead of the reaction |
---|
| 454 | aParticleChange.ProposeEnergy(kinEnergy); |
---|
| 455 | aParticleChange.ProposeLocalEnergyDeposit(0.); |
---|
| 456 | aParticleChange.ProposeMomentumDirection(dir) ; |
---|
| 457 | return G4VDiscreteProcess::PostStepDoIt(track,step); |
---|
| 458 | } |
---|
| 459 | G4double mint=CalculateXSt(false, false, Momentum, Z, N, projPDG);// randomize t in MeV^2 |
---|
| 460 | #ifdef pdebug |
---|
| 461 | G4cout<<"G4QCohChEx::PoStDoIt:pPDG="<<projPDG<<",tPDG="<<targPDG<<",P="<<Momentum<<",CS=" |
---|
| 462 | <<xSec<<",-t="<<mint<<G4endl; |
---|
| 463 | #endif |
---|
| 464 | #ifdef nandebug |
---|
| 465 | if(mint>-.0000001); |
---|
| 466 | else G4cout<<"*Warning*G4QCoherentChargeExchange::PostStDoIt:-t="<<mint<<G4endl; |
---|
| 467 | #endif |
---|
| 468 | G4double maxt=CalculateXSt(true, false, Momentum, Z, N, projPDG); |
---|
| 469 | if(maxt<=0.) maxt=1.e-27; |
---|
| 470 | G4double cost=1.-mint/maxt; // cos(theta) in CMS (@@ we neglect mass diff for ChEx) |
---|
| 471 | // |
---|
| 472 | #ifdef ppdebug |
---|
| 473 | G4cout<<"G4QCoherentChargeExchange::PoStDoIt:t="<<mint<<",dpcm2="<<maxt |
---|
| 474 | <<",Ek="<<kinEnergy<<",tM="<<tM<<",pM="<<pM<<",cost="<<cost<<G4endl; |
---|
| 475 | #endif |
---|
| 476 | if(cost>1. || cost<-1. || !(cost>-1. || cost<1.)) |
---|
| 477 | { |
---|
| 478 | if(cost>1.000001 || cost<-1.000001 || !(cost>-1. || cost<1.)) |
---|
| 479 | { |
---|
| 480 | G4double tM2=tM*tM; // Squared target mass |
---|
| 481 | G4double pEn=pM+kinEnergy; // tot projectile Energy in MeV |
---|
| 482 | G4double sM=(tM+tM)*pEn+tM2+pM2; // Mondelstam s |
---|
| 483 | G4double twop2cm=(tM2+tM2)*(pEn*pEn-pM2)/sM;// Max_t/2 (2*p^2_cm) |
---|
| 484 | G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:cos="<<cost<<",t="<<mint<<",T="<<kinEnergy |
---|
| 485 | <<",tM="<<tM<<",tmax="<<2*kinEnergy*tM<<",p="<<projPDG<<",t="<<targPDG<<G4endl; |
---|
| 486 | G4cout<<"..G4QCohChEx::PoStDoI: dpcm2="<<twop2cm<<"="<<maxt<<G4endl; |
---|
| 487 | } |
---|
| 488 | if (cost>1.) cost=1.; |
---|
| 489 | else if(cost<-1.) cost=-1.; |
---|
| 490 | } |
---|
| 491 | G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM); // 4mom of the recoil target |
---|
| 492 | G4LorentzVector scat4M=G4LorentzVector(0.,0.,0.,pM); // 4mom of the recoil target |
---|
| 493 | G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01); |
---|
| 494 | if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost)) |
---|
| 495 | { |
---|
| 496 | G4cerr<<"G4QCohChEx::PSDI:t4M="<<tot4M<<",pM="<<pM<<",tM="<<tM<<",cost="<<cost<<G4endl; |
---|
| 497 | //G4Exception("G4QCoherentChargeExchange::PostStepDoIt:","009",FatalException,"Decay"); |
---|
| 498 | } |
---|
| 499 | #ifdef debug |
---|
| 500 | G4cout<<"G4QCohChEx::PoStDoIt:s4M="<<scat4M<<"+r4M="<<reco4M<<"="<<scat4M+reco4M<<G4endl; |
---|
| 501 | G4cout<<"G4QCohChEx::PoStDoIt: scatE="<<scat4M.e()-pM<<", recoE="<<reco4M.e()-tM<<",d4M=" |
---|
| 502 | <<tot4M-scat4M-reco4M<<G4endl; |
---|
| 503 | #endif |
---|
| 504 | // Kill scattered hadron |
---|
| 505 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
| 506 | // Definition of the scattered nucleon |
---|
| 507 | G4DynamicParticle* theSec = new G4DynamicParticle; // A secondary for the recoil hadron |
---|
| 508 | G4ParticleDefinition* theDefinition=G4Proton::Proton(); |
---|
| 509 | if(projPDG==2212) theDefinition=G4Neutron::Neutron(); |
---|
| 510 | theSec->SetDefinition(theDefinition); |
---|
| 511 | EnMomConservation-=scat4M; |
---|
| 512 | theSec->Set4Momentum(scat4M); |
---|
| 513 | G4Track* aNewTrack = new G4Track(theSec, localtime, position ); |
---|
| 514 | aNewTrack->SetWeight(weight); // weighted |
---|
| 515 | aNewTrack->SetTouchableHandle(trTouchable); |
---|
| 516 | aParticleChange.AddSecondary( aNewTrack ); |
---|
| 517 | // Filling the recoil nucleus |
---|
| 518 | theSec = new G4DynamicParticle; // A secondary for the recoil hadron |
---|
| 519 | G4int aA = Z+N; |
---|
| 520 | #ifdef pdebug |
---|
| 521 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Ion Z="<<Z<<", A="<<aA<<G4endl; |
---|
| 522 | #endif |
---|
| 523 | theDefinition=G4ParticleTable::GetParticleTable()->FindIon(Z,aA,0,Z); |
---|
| 524 | if(!theDefinition)G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:drop PDG="<<targPDG<<G4endl; |
---|
| 525 | #ifdef pdebug |
---|
| 526 | G4cout<<"G4QCohChEx::PostStepDoIt:RecoilName="<<theDefinition->GetParticleName()<<G4endl; |
---|
| 527 | #endif |
---|
| 528 | theSec->SetDefinition(theDefinition); |
---|
| 529 | EnMomConservation-=reco4M; |
---|
| 530 | #ifdef tdebug |
---|
| 531 | G4cout<<"G4QCohChEx::PostSDoIt:"<<targPDG<<reco4M<<reco4M.m()<<EnMomConservation<<G4endl; |
---|
| 532 | #endif |
---|
| 533 | theSec->Set4Momentum(reco4M); |
---|
| 534 | #ifdef debug |
---|
| 535 | G4ThreeVector curD=theSec->GetMomentumDirection(); |
---|
| 536 | G4double curM=theSec->GetMass(); |
---|
| 537 | G4double curE=theSec->GetKineticEnergy()+curM; |
---|
| 538 | G4cout<<"G4QCohChEx::PostStpDoIt:p="<<curD<<curD.mag()<<",e="<<curE<<",m="<<curM<<G4endl; |
---|
| 539 | #endif |
---|
| 540 | // Make a recoil nucleus |
---|
| 541 | aNewTrack = new G4Track(theSec, localtime, position ); |
---|
| 542 | aNewTrack->SetWeight(weight); // weighted |
---|
| 543 | aNewTrack->SetTouchableHandle(trTouchable); |
---|
| 544 | aParticleChange.AddSecondary( aNewTrack ); |
---|
| 545 | #ifdef debug |
---|
| 546 | G4cout<<"G4QCoherentChargeExchange::PostStepDoIt:*** PostStepDoIt is done ***"<<G4endl; |
---|
| 547 | #endif |
---|
| 548 | return G4VDiscreteProcess::PostStepDoIt(track, step); |
---|
| 549 | } |
---|
| 550 | |
---|
| 551 | G4double G4QCoherentChargeExchange::CalculateXSt(G4bool oxs, G4bool xst, G4double p, |
---|
| 552 | G4int Z, G4int N, G4int pPDG) |
---|
| 553 | { |
---|
| 554 | static G4bool init=false; |
---|
| 555 | static G4bool first=true; |
---|
| 556 | static G4VQCrossSection* CSmanager; |
---|
| 557 | G4QuasiFreeRatios* qfMan=G4QuasiFreeRatios::GetPointer(); |
---|
| 558 | if(first) // Connection with a singletone |
---|
| 559 | { |
---|
| 560 | CSmanager=G4QElasticCrossSection::GetPointer(); |
---|
| 561 | first=false; |
---|
| 562 | } |
---|
| 563 | G4double res=0.; |
---|
| 564 | if(oxs && xst) // Only the Cross-Section can be returened |
---|
| 565 | { |
---|
| 566 | res=CSmanager->GetCrossSection(true, p, Z, N, pPDG); // XS for isotope |
---|
| 567 | res*=qfMan->ChExElCoef(p*MeV, Z, N, pPDG); |
---|
| 568 | } |
---|
| 569 | else if(!oxs && xst) // Calculate Cross-Section & prepare differential |
---|
| 570 | { |
---|
| 571 | res=CSmanager->GetCrossSection(false, p, Z, N, pPDG);// XS for isotope + init t-distr. |
---|
| 572 | res*=qfMan->ChExElCoef(p*MeV, Z, N, pPDG); |
---|
| 573 | // The XS for the nucleus must be calculated the last |
---|
| 574 | init=true; |
---|
| 575 | } |
---|
| 576 | else if(init) // Return t-value for scattering (=G4QElastic) |
---|
| 577 | { |
---|
| 578 | if(oxs) res=CSmanager->GetHMaxT(); // Calculate the max_t value |
---|
| 579 | else res=CSmanager->GetExchangeT(Z, N, pPDG); // fanctionally randomized -t in MeV^2 |
---|
| 580 | } |
---|
| 581 | else G4cout<<"*Warning*G4QCohChrgExchange::CalculateXSt: NotInitiatedScattering"<<G4endl; |
---|
| 582 | return res; |
---|
| 583 | } |
---|