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