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