source: trunk/source/processes/hadronic/models/leading_particle/include/G4Mars5GeV.hh@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 9.5 KB
Line 
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//
27// $Id: G4Mars5GeV.hh,v 1.5 2006/06/29 20:43:22 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34//
35// ------------------------------------------------------------
36// First Implemention 17 Nov. 1998 M.Asai, H.Kurahige
37//
38// History:
39// modified as hadronic model 28 Oct 2001 N.Kanaya
40// remodified by T. Koi 18 Nov 2004
41// ------------------------------------------------------------
42// Class Description
43// This is a Event Biasing mechanism based on MARS code
44// This model is applicable to
45// proton/neutron/pi+-/K+-/gamma/anti_proton
46// with energy < 5.0GeV
47//*
48// Original code is MARS13 written by Nikolai Mokhov (FNAL)
49//**************************************************************
50//* MARS13: 9. hA EVENT GENERATOR:
51//* Copyright Nikolai Mokhov (Fermilab)
52//*
53//* LAST CHANGE: 14-NOV-1998
54//**************************************************************
55//* Copyright Nikolai Mokhov (Fermilab)
56//*
57//* MARS13(98)
58//*
59//* INCLUSIVE HADRON(photon)-NUCLEUS VERTEX AT E < 5 GEV !!!
60//* THREE WEIGHTED HADRONS IN FINAL STATE: !!!
61//* IP+A -> N/P(CASC)+ PI+/PI-(K+/K-) + PI0
62//
63
64#ifndef G4Mars5GeV_h
65#define G4Mars5GeV_h 1
66#include "G4Element.hh"
67#include "G4ElementVector.hh"
68#include "G4ElementTable.hh"
69#include "G4PhysicsTable.hh"
70#include "G4PhysicsVector.hh"
71#include "G4LPhysicsFreeVector.hh"
72#include "G4LightMedia.hh"
73#include "G4Step.hh"
74#include "G4TrackStatus.hh"
75#include "G4InelasticInteraction.hh"
76#include "G4ParticleTable.hh"
77#include "G4ParticleDefinition.hh"
78
79class G4Mars5GeV : public G4InelasticInteraction
80{
81public: // with description
82 // constructors
83 G4Mars5GeV();
84
85 //destructor
86 ~G4Mars5GeV() { ; };
87
88 // virtual methods derived from G4VEvtBiasMechanism
89 G4HadFinalState* ApplyYourself( const G4HadProjectile &aTrack,
90 G4Nucleus &targetNucleus );
91
92 private:
93 void Treem5();
94 // This is the mothod which invoke MARS
95 // secondary information will be filled up
96
97 G4bool CoulombBarrier(G4int pType, G4double pE);
98 // Check if coulomb barrier exists
99
100 void CreateNucleon(G4int ib, G4int pType, G4double pE);
101 void CreatePion(G4int ib, G4int pType, G4double pE);
102 void CreatePionZero(G4int ib, G4int pType, G4double pE);
103 // Create secondary particles and add them into the list
104
105 //void AddSecondary();
106 void AddSecondaryToMarsList();
107 // Add a secondary particle into the list
108
109 G4double SelBS(G4int pType, G4double aNucl, G4double zNucl);
110 // Calculate weight of secondary
111
112 G4double D2N2(G4int pType, G4double incidentE,
113 G4double prodE, G4double tin,
114 G4int reacType, G4int proType,
115 G4double ai, G4double z);
116 // Calculate Hadron Inclusive Yield
117
118 G4double Rkaon(G4int ib, G4int jp, G4double eRaw);
119 // Calculate energy dependent K/pi ratio
120
121 void Trans(G4ThreeVector* d1, G4ThreeVector* d2);
122 // Direction cosine transformation using selec2(cs,ss,ch,sh)
123
124 G4double GetMaxWeight() const;
125 G4double GetMinWeight() const;
126 void SetMaxWeight(G4double);
127 void SetMinWeight(G4double);
128
129 public:
130 enum {FastVectorSize = 16};
131 typedef G4FastVector<G4DynamicParticle ,FastVectorSize> G4MarsSecondaryVector;
132
133 private:
134 // weight max/min
135 G4double maxWeight;
136 G4double minWeight;
137
138 // information of secondary
139 G4int numberOfSecondaries;
140 G4double weightOfSecondaries[FastVectorSize];
141 G4MarsSecondaryVector secondaries;
142
143 private:
144 G4bool IsApplicable(const G4HadProjectile& , G4Nucleus& ) { return true; };
145 G4bool IsApplicable(G4int marsEncoding) const;
146
147 private:
148 // Particle Change
149
150 // Particle Table
151 G4ParticleTable* theParticleTable;
152
153 // particle encoding for MARS
154 enum { MarsUndefined =0,
155 MarsP, MarsN, MarsPIplus, MarsPIminus, MarsKplus, MarsKminus,
156 MarsMUplus, MarsMUminus, MarsGAM, MarsEminus, MarsEplus, MarsAP,
157 MarsPI0, MarsD, MarsT, MarsHe3, MarsHe4 };
158
159 G4int GetMarsEncoding(const G4ParticleDefinition* )const;
160 const G4String& GetParticleName(G4int marsEncoding) const;
161 G4ParticleDefinition* GetParticleDefinition(G4int marsEncoding) const;
162 G4double ProtonMass;
163
164 private:
165 // incident information
166 G4double incidentWeight;
167 const G4HadProjectile* incidentParticle;
168 G4int incidentMarsEncoding;
169
170 void GetTargetNuclei(const G4Material*); //fill up fANucle and fZnucl
171 G4double fANucl, fZNucl; // target nucleus
172
173 private:
174 // these class is to define common blocks in original code
175
176 class Selec1
177 {
178 public:
179 G4double Einc, EN, V, V10;
180 G4int Treac, Tprod;
181 } selec1;
182 class Selec2
183 {
184 public:
185 G4double Cs, Ss, Ch, Sh;
186 } selec2;
187 class Selec3
188 {
189 public:
190 G4double Eth, Emax, Sqs, X, Pt2, Pt, P;
191 } selec3;
192};
193
194#include "G4ParticleDefinition.hh"
195inline
196 G4double G4Mars5GeV::GetMaxWeight() const
197{
198 return maxWeight;
199}
200inline
201 G4double G4Mars5GeV::GetMinWeight() const
202{
203 return minWeight;
204}
205
206inline
207 void G4Mars5GeV::SetMaxWeight(G4double value)
208{
209 maxWeight = value;
210}
211inline
212 void G4Mars5GeV::SetMinWeight(G4double value)
213{
214 minWeight = value;
215}
216
217inline
218 G4int G4Mars5GeV::GetMarsEncoding(const G4ParticleDefinition* particle) const{
219 const G4String& name = particle->GetParticleName();
220 G4int encoding = MarsUndefined;
221 if (name == "proton") {
222 encoding = MarsP;
223 } else if (name == "neutron") {
224 encoding = MarsN;
225 } else if (name == "pi+") {
226 encoding = MarsPIplus;
227 } else if (name == "pi-") {
228 encoding = MarsPIminus;
229 } else if (name == "kaon+") {
230 encoding = MarsKplus;
231 } else if (name == "kaon-") {
232 encoding = MarsKminus;
233 } else if (name == "mu+") {
234 encoding = MarsMUplus;
235 } else if (name == "mu-") {
236 encoding = MarsMUminus;
237 } else if (name == "gamma") {
238 encoding = MarsGAM;
239 } else if (name == "e+") {
240 encoding = MarsEplus;
241 } else if (name == "e-") {
242 encoding = MarsEminus;
243 } else if (name == "anti_proton") {
244 encoding = MarsAP;
245 } else if (name == "pi0") {
246 encoding = MarsPI0;
247 } else if (name == "deuteron") {
248 encoding = MarsD;
249 } else if (name == "triton") {
250 encoding = MarsT;
251 } else if (name == "He3") {
252 encoding = MarsHe3;
253 } else if (name == "alpha") {
254 encoding = MarsHe4;
255 }
256 return encoding;
257}
258
259inline
260 const G4String& G4Mars5GeV::GetParticleName(G4int encoding) const
261{
262 static G4String name;
263 name = "None";
264 switch (encoding)
265 {
266 case MarsP:
267 name = "proton";
268 break;
269 case MarsN:
270 name = "neutron";
271 break;
272 case MarsPIplus:
273 name = "pi+";
274 break;
275 case MarsPIminus:
276 name = "pi-";
277 break;
278 case MarsKplus:
279 name = "kaon+";
280 break;
281 case MarsKminus:
282 name = "kaon-";
283 break;
284 case MarsMUplus:
285 name = "mu+";
286 break;
287 case MarsMUminus:
288 name = "mu-";
289 break;
290 case MarsGAM:
291 name = "gamma";
292 break;
293 case MarsEplus:
294 name = "e+";
295 break;
296 case MarsEminus:
297 name = "e-";
298 break;
299 case MarsAP:
300 name = "anti_proton";
301 break;
302 case MarsPI0:
303 name = "pi0";
304 break;
305 case MarsD:
306 name = "deuteron";
307 break;
308 case MarsT:
309 name = "triton";
310 break;
311 case MarsHe3:
312 name = "He3";
313 break;
314 case MarsHe4:
315 name = "alpha";
316 break;
317 default:
318 break;
319 }
320 return name;
321}
322
323inline
324G4ParticleDefinition* G4Mars5GeV::GetParticleDefinition(G4int encoding) const
325{
326 G4String name = GetParticleName(encoding);
327 G4ParticleDefinition* particle = 0;
328 if (name != "None") {
329 particle = theParticleTable->FindParticle(name);
330 }
331 return particle;
332}
333
334inline
335 G4bool G4Mars5GeV::IsApplicable(G4int marsEncoding) const
336{
337 return ( ((marsEncoding!=MarsUndefined) && (marsEncoding<=MarsKminus) )||
338 (marsEncoding==MarsGAM) ||
339 (marsEncoding==MarsAP) );
340}
341
342#endif
Note: See TracBrowser for help on using the repository browser.