FairRoot/PandaRoot
Test30Physics.cxx
Go to the documentation of this file.
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 // -------------------------------------------------------------
28 // GEANT 4 class for test30
29 //
30 // History: based on object model of
31 // ---------- Test30Physics -------
32 // by Vladimir Ivanchenko, 12 March 2002
33 //
34 // Modified:
35 // 11.10.2007 Added INCL cascade and RPG parameterized model (V.Ivanchenko)
36 //
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39 
40 #include "Test30Physics.hh"
41 
42 #include "G4SystemOfUnits.hh"
43 #include "G4UnitsTable.hh"
44 #include "Test30VSecondaryGenerator.hh"
45 #include "G4PhysicsTable.hh"
46 #include "G4ProcessManager.hh"
47 #include "G4ProcessVector.hh"
48 #include "G4ParticleTypes.hh"
49 
50 #include "G4ProcessManager.hh"
51 #include "G4ParticleDefinition.hh"
52 #include "G4DecayPhysics.hh"
53 
54 #include "G4PreCompoundModel.hh"
55 #include "G4ExcitationHandler.hh"
56 #include "G4BinaryCascade.hh"
57 #include "G4BinaryLightIonReaction.hh"
58 #include "G4CascadeInterface.hh"
59 #include "G4WilsonAbrasionModel.hh"
60 #include "G4QMDReaction.hh"
61 #include "G4INCLXXInterface.hh"
62 
63 #include "G4TheoFSGenerator.hh"
64 #include "G4FTFModel.hh"
65 #include "G4ExcitedStringDecay.hh"
66 #include "G4LundStringFragmentation.hh"
67 
68 #include "G4QGSModel.hh"
69 #include "G4QGSParticipants.hh"
70 #include "G4QGSMFragmentation.hh"
71 #include "G4QuasiElasticChannel.hh"
72 
73 #include "G4GeneratorPrecompoundInterface.hh"
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
77 Test30Physics::Test30Physics()
78 {
79  Initialise();
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
84 Test30Physics::~Test30Physics()
85 {}
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
89 void Test30Physics::Initialise()
90 {
91  G4DecayPhysics dp;
92  dp.ConstructParticle();
93  theProcess = 0;
94  theDeExcitation = 0;
95  thePreCompound = 0;
96  theQuasiElastic = 0;
97  hkmod = 0;
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
102 G4VProcess* Test30Physics::GetProcess(const G4String& gen_name,
103  const G4ParticleDefinition* part,
104  G4Material* mat)
105 {
106  G4cout << "Test30Physics entry" << G4endl;
107  if(theProcess) { G4cout<<"AAHHH theProcess is there and will be destroyed! " << theProcess<<G4endl; delete theProcess; }
108  theProcess = 0;
109 
110  G4String part_name = part->GetParticleName();
111  G4ProcessManager* man = new G4ProcessManager(part);
112  if(!man) { G4cout<<"Test30Physics::GetProcess(): ERROR G4ProcessManager is not there!"<< G4endl; return 0; }
113 
114  theProcess = new Test30HadronProduction();
115 
116  G4cout << "Process is created ("<<theProcess<<"); gen= " << gen_name << G4endl;
117 
118  if(!theDeExcitation) {
119  theDeExcitation = new G4ExcitationHandler();
120  G4cout << "New deexcitation" << G4endl;
121  }
122  if(!thePreCompound) {
123  thePreCompound = new G4PreCompoundModel(theDeExcitation);
124  G4cout << "New pre-compound" << G4endl;
125  }
126 
127  // Physics list for the given run
128  Test30VSecondaryGenerator* sg = 0;
129 
130  // Choose generator
131  if(gen_name == "preco") {
132  sg = new Test30VSecondaryGenerator(thePreCompound,mat);
133  theProcess->SetSecondaryGenerator(sg);
134  man->AddDiscreteProcess(theProcess);
135 
136  } else if(gen_name == "binary") {
137  G4BinaryCascade* hkm = new G4BinaryCascade();
138  hkm->SetDeExcitation(thePreCompound);
139  sg = new Test30VSecondaryGenerator(hkm, mat);
140  theProcess->SetSecondaryGenerator(sg);
141  man->AddDiscreteProcess(theProcess);
142 
143  } else if(gen_name == "binary_ion") {
144  G4BinaryLightIonReaction* hkm = new G4BinaryLightIonReaction();
145  //hkm->SetDeExcitation(theDeExcitation);
146  // hkm->SetPrecompound(thePreCompound);
147  sg = new Test30VSecondaryGenerator(hkm, mat);
148  theProcess->SetSecondaryGenerator(sg);
149  man->AddDiscreteProcess(theProcess);
150 
151  } else if(gen_name == "ftfp") {
152 
153  G4TheoFSGenerator* theModel = new G4TheoFSGenerator;
154  G4FTFModel* theStringModel = new G4FTFModel();
155  G4GeneratorPrecompoundInterface* theCascade =
156  new G4GeneratorPrecompoundInterface;
157  theCascade->SetDeExcitation(thePreCompound);
158  G4ExcitedStringDecay* theStringDecay =
159  new G4ExcitedStringDecay(new G4LundStringFragmentation());
160  theStringModel->SetFragmentationModel(theStringDecay);
161 
162  theModel->SetTransport(theCascade);
163  theModel->SetHighEnergyGenerator(theStringModel);
164  theModel->SetMinEnergy(GeV);
165 
166  sg = new Test30VSecondaryGenerator(theModel, mat);
167  theProcess->SetSecondaryGenerator(sg);
168  man->AddDiscreteProcess(theProcess);
169 
170  } else if(gen_name == "ftfb") {
171 
172  G4TheoFSGenerator* theModel = new G4TheoFSGenerator();
173  G4FTFModel* theStringModel= new G4FTFModel();
174  G4BinaryCascade* theCascade = new G4BinaryCascade();
175  theCascade->SetDeExcitation(thePreCompound);
176  G4ExcitedStringDecay * theStringDecay =
177  new G4ExcitedStringDecay(new G4LundStringFragmentation());
178  theModel->SetHighEnergyGenerator(theStringModel);
179  theStringModel->SetFragmentationModel(theStringDecay);
180 
181  theModel->SetTransport(theCascade);
182  theModel->SetHighEnergyGenerator(theStringModel);
183  theModel->SetMinEnergy(GeV);
184 
185  sg = new Test30VSecondaryGenerator(theModel, mat);
186  theProcess->SetSecondaryGenerator(sg);
187  man->AddDiscreteProcess(theProcess);
188 
189  } else if(gen_name == "qgsp") {
190  G4TheoFSGenerator* theModel = new G4TheoFSGenerator();
191  G4GeneratorPrecompoundInterface* theCascade =
192  new G4GeneratorPrecompoundInterface();
193  theCascade->SetDeExcitation(thePreCompound);
194  G4QGSModel< G4QGSParticipants > * theStringModel =
195  new G4QGSModel< G4QGSParticipants >;
196  G4ExcitedStringDecay* theQGStringDecay =
197  new G4ExcitedStringDecay(new G4QGSMFragmentation());
198  theStringModel->SetFragmentationModel(theQGStringDecay);
199  theModel->SetTransport(theCascade);
200  theModel->SetHighEnergyGenerator(theStringModel);
201  G4QuasiElasticChannel* quasiElastic = new G4QuasiElasticChannel;
202  theModel->SetQuasiElasticChannel(quasiElastic);
203  theModel->SetMinEnergy(0.5*GeV);
204  theModel->SetMaxEnergy(100*TeV);
205  sg = new Test30VSecondaryGenerator(theModel, mat);
206  theProcess->SetSecondaryGenerator(sg);
207  man->AddDiscreteProcess(theProcess);
208 
209  } else if(gen_name == "qgsb") {
210  G4TheoFSGenerator* theModel = new G4TheoFSGenerator();
211  G4BinaryCascade* theCascade = new G4BinaryCascade();
212  theCascade->SetDeExcitation(thePreCompound);
213  G4QGSModel< G4QGSParticipants > * theStringModel =
214  new G4QGSModel< G4QGSParticipants >;
215  G4ExcitedStringDecay* theQGStringDecay =
216  new G4ExcitedStringDecay(new G4QGSMFragmentation());
217  theStringModel->SetFragmentationModel(theQGStringDecay);
218  theModel->SetTransport(theCascade);
219  theModel->SetHighEnergyGenerator(theStringModel);
220  G4QuasiElasticChannel* quasiElastic = new G4QuasiElasticChannel;
221  theModel->SetQuasiElasticChannel(quasiElastic);
222  theModel->SetMinEnergy(0.5*GeV);
223  theModel->SetMaxEnergy(100*TeV);
224  sg = new Test30VSecondaryGenerator(theModel, mat);
225  theProcess->SetSecondaryGenerator(sg);
226  man->AddDiscreteProcess(theProcess);
227 
228  } else if(gen_name == "bertini") {
229  G4CascadeInterface* hkm = new G4CascadeInterface();
230  sg = new Test30VSecondaryGenerator(hkm, mat);
231  theProcess->SetSecondaryGenerator(sg);
232  man->AddDiscreteProcess(theProcess);
233 
234  } else if(gen_name == "bertini_preco") {
235  G4CascadeInterface* hkm = new G4CascadeInterface();
236  hkm->usePreCompoundDeexcitation();
237  sg = new Test30VSecondaryGenerator(hkm, mat);
238  theProcess->SetSecondaryGenerator(sg);
239  man->AddDiscreteProcess(theProcess);
240 
241  } else if(gen_name == "incl") {
242  G4INCLXXInterface* hkm = new G4INCLXXInterface();
243  sg = new Test30VSecondaryGenerator(hkm, mat);
244  theProcess->SetSecondaryGenerator(sg);
245  man->AddDiscreteProcess(theProcess);
246 
247  } else if(gen_name == "incl_ion") {
248  G4INCLXXInterface* hkm = new G4INCLXXInterface();
249  sg = new Test30VSecondaryGenerator(hkm, mat);
250  theProcess->SetSecondaryGenerator(sg);
251  man->AddDiscreteProcess(theProcess);
252 
253  } else if(gen_name == "abrasion") {
254  G4WilsonAbrasionModel* wam = new G4WilsonAbrasionModel();
255  size_t ne = mat->GetNumberOfElements();
256  const G4ElementVector* ev = mat->GetElementVector();
257  for(size_t i=0; i<ne; i++) {
258  G4Element* elm = (*ev)[i];
259  wam->ActivateFor(elm);
260  }
261  sg = new Test30VSecondaryGenerator(wam, mat);
262  theProcess->SetSecondaryGenerator(sg);
263  man->AddDiscreteProcess(theProcess);
264 
265  } else if(gen_name == "qmd") {
266  G4QMDReaction* qmd = new G4QMDReaction();
267  sg = new Test30VSecondaryGenerator(qmd, mat);
268  theProcess->SetSecondaryGenerator(sg);
269  man->AddDiscreteProcess(theProcess);
270 
271  } else {
272  G4cout << "WARNING: <" << gen_name << "> generator is unknown" << G4endl;
273  }
274 
275  G4cout << "Secondary generator <"
276  << gen_name << "> is initialized"
277  << G4endl;
278 
279  TestPointers();
280  return theProcess;
281 
282 }
283 
284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
285 
286 
287 
288 
289 
290 
291 
Int_t i
Definition: run_full.C:25
int ev
int ne
Definition: toy_core.C:124