source: trunk/examples/advanced/raredecay_calorimetry/src/PhotInPhysicsList.cc@ 1036

Last change on this file since 1036 was 807, checked in by garnier, 17 years ago

update

File size: 28.1 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: PhotInPhysicsList.cc,v 1.4 2006/06/29 16:25:19 gunter Exp $
28// GEANT4 tag $Name: $
29//
30
31#define debug
32
33#include "PhotInPhysicsList.hh"
34
35PhotInPhysicsList::PhotInPhysicsList(): G4VUserPhysicsList() { SetVerboseLevel(1); }
36
37PhotInPhysicsList::~PhotInPhysicsList() {}
38
39void PhotInPhysicsList::ConstructParticle()
40{
41 // In this method, static member functions for particles should be called
42 // for all particles which user is going to use in the simulation. If not
43 // defined particle appear in the simulation, it can cause a WORNING which
44 // means that in the simulation appeard unexpected particles. Then add them.
45
46 //// @@ Word "Definition" can be skipped. - Old fashion(M.K.)
47 //
48 //// pseudo-particles
49 //G4Geantino::GeantinoDefinition();
50 //G4ChargedGeantino::ChargedGeantinoDefinition();
51 //
52 //// gammas
53 //G4Gamma::GammaDefinition();
54 //
55 //// leptons (without tau and it's neutrino)
56 //G4Electron::ElectronDefinition();
57 //G4Positron::PositronDefinition();
58 //G4MuonPlus::MuonPlusDefinition();
59 //G4MuonMinus::MuonMinusDefinition();
60 //
61 //G4NeutrinoE::NeutrinoEDefinition();
62 //G4AntiNeutrinoE::AntiNeutrinoEDefinition();
63 //G4NeutrinoMu::NeutrinoMuDefinition();
64 //G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
65 //
66 //// mesons
67 //G4PionPlus::PionPlusDefinition();
68 //G4PionMinus::PionMinusDefinition();
69 //G4PionZero::PionZeroDefinition();
70 //G4Eta::EtaDefinition();
71 //G4EtaPrime::EtaPrimeDefinition();
72 //G4KaonPlus::KaonPlusDefinition();
73 //G4KaonMinus::KaonMinusDefinition();
74 //G4KaonZero::KaonZeroDefinition();
75 //G4AntiKaonZero::AntiKaonZeroDefinition();
76 //G4KaonZeroLong::KaonZeroLongDefinition();
77 //G4KaonZeroShort::KaonZeroShortDefinition();
78 //
79 //// barions
80 //G4Proton::ProtonDefinition();
81 //G4AntiProton::AntiProtonDefinition();
82 //G4Neutron::NeutronDefinition();
83 //G4AntiNeutron::AntiNeutronDefinition();
84 //
85 //// hyperons
86 //G4Lambda::LambdaDefinition();
87 //G4SigmaPlus::SigmaPlusDefinition();
88 //G4SigmaZero::SigmaZeroDefinition();
89 //G4SigmaMinus::SigmaMinusDefinition();
90 //G4XiMinus::XiMinusDefinition();
91 //G4XiZero::XiZeroDefinition();
92 //G4OmegaMinus::OmegaMinusDefinition();
93 //G4AntiLambda::AntiLambdaDefinition();
94 //G4AntiSigmaPlus::AntiSigmaPlusDefinition();
95 //G4AntiSigmaZero::AntiSigmaZeroDefinition();
96 //G4AntiSigmaMinus::AntiSigmaMinusDefinition();
97 //G4AntiXiMinus::AntiXiMinusDefinition();
98 //G4AntiXiZero::AntiXiZeroDefinition();
99 //G4AntiOmegaMinus::AntiOmegaMinusDefinition();
100
101 // The same can be done shorter, using supported Geant4 particle initialization classes:
102 G4BosonConstructor consBos; consBos.ConstructParticle();
103 G4LeptonConstructor consLep; consLep.ConstructParticle();
104 G4MesonConstructor consMes; consMes.ConstructParticle();
105 G4BaryonConstructor consBar; consBar.ConstructParticle();
106 G4IonConstructor consIon; consIon.ConstructParticle();
107 G4ShortLivedConstructor consShL; consShL.ConstructParticle();
108}
109
110
111void PhotInPhysicsList::ConstructProcess()
112{
113#ifdef debug
114 G4cout<<"PhotInPhysicsList::ConstructProcess: is called "<<G4endl;
115#endif
116 AddTransportation(); // Transportation is a "process" and defined in the basic class
117
118 // Add Electromagnetic interaction Processes and Decays
119 G4Decay* theDecayProcess = new G4Decay(); // @@ When this class is decayed? (M.K.)
120
121 // For high energy hadronic (QGSC): this will be the model class for high energies
122 G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator;
123
124 // a cascade interface interface to chiral invriant phase space decay
125 // at particle level.
126 G4StringChipsParticleLevelInterface* theCascade =
127 new G4StringChipsParticleLevelInterface();
128
129 // here come the high energy parts
130 // the string model; still not quite according to design - Explicite use of the forseen
131 // interfaces will be tested and documented in this program by beta-02 at latest.
132 G4VPartonStringModel* theStringModel;
133 theStringModel = new G4QGSModel<G4QGSParticipants>;
134 theTheoModel->SetTransport(theCascade);
135 theTheoModel->SetHighEnergyGenerator(theStringModel);
136 theTheoModel->SetMinEnergy(19*GeV);
137 theTheoModel->SetMaxEnergy(100*TeV);
138
139 G4VLongitudinalStringDecay* theFragmentation = new G4LundStringFragmentation;
140 G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theFragmentation);
141 theStringModel->SetFragmentationModel(theStringDecay);
142
143 // done with the generator model (most of the above is also available as default)
144 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess; //@@
145 G4LElastic* theElasticModel = new G4LElastic; //@@
146 theElasticProcess->RegisterMe(theElasticModel);
147 G4HadronElasticProcess* theElasticProcess1 = new G4HadronElasticProcess;
148
149 theParticleIterator->reset();
150 while( (*theParticleIterator)() )
151 {
152 G4ParticleDefinition* particle = theParticleIterator->value();
153#ifdef debug
154 G4cout<<"PhotInPhysList::ConstructProcess: Part="<<particle->GetParticleName()<<G4endl;
155#endif
156 G4ProcessManager* pmanager = particle->GetProcessManager();
157
158 // ================================ Decays ===========================================
159
160 if (theDecayProcess->IsApplicable(*particle))
161 {
162 pmanager ->AddProcess(theDecayProcess);
163 // set ordering for PostStepDoIt and AtRestDoIt
164 pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
165 pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
166 }
167
168 G4String particleName = particle->GetParticleName();
169
170 // =========================== EM Interactions & Photo-nuclear =======================
171
172 if (particleName == "gamma") // gamma
173 {
174 // Construct processes for gamma
175 // ........The OLD (GHAR) usage of the CHIPS Photo-nuclear reactions .............
176 //G4PhotoNuclearProcess* thePhotoNuclearProcess = new G4PhotoNuclearProcess;
177 //G4QGSMFragmentation* theFragmentation = new G4QGSMFragmentation;
178 //G4QGSModel<G4GammaParticipants>*theStringModel=new G4QGSModel<G4GammaParticipants>;
179 //G4GammaNuclearReaction* theGammaReaction = new G4GammaNuclearReaction;
180 //G4TheoFSGenerator* theModel = new G4TheoFSGenerator;
181 //G4StringChipsParticleLevelInterface* theCascade =
182 // new G4StringChipsParticleLevelInterface;
183 //theModel->SetTransport(theCascade);
184 //theModel->SetHighEnergyGenerator(theStringModel);
185 //G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theFragmentation);
186 //theStringModel->SetFragmentationModel(theStringDecay);
187 //theGammaReaction->SetMaxEnergy(3.5*GeV);
188 //thePhotoNuclearProcess->RegisterMe(theGammaReaction);
189 //theModel->SetMinEnergy(3.*GeV);
190 //theModel->SetMaxEnergy(100*TeV);
191 //thePhotoNuclearProcess->RegisterMe(theModel);
192 //pmanager->AddDiscreteProcess(thePhotoNuclearProcess);
193 // ........The NEW (native) usage of the CHIPS Photo-nuclear reactions .............
194 pmanager->AddDiscreteProcess(new G4QCollision);
195
196 // Pure Electromagnetic Processes
197 pmanager->AddDiscreteProcess(new G4GammaConversion()); // Pair Production
198 pmanager->AddDiscreteProcess(new G4ComptonScattering()); // Compton Effect
199 pmanager->AddDiscreteProcess(new G4PhotoElectricEffect()); // Photo Effect
200 }
201 else if (particleName == "e-") //electron
202 {
203 // Construct processes for electron
204 // ........The OLD (GHAR) usage of the CHIPS Photo-nuclear reactions .............
205 //G4ElectronNuclearProcess* theElectronNuclearProcess = new G4ElectronNuclearProcess;
206 //G4ElectroNuclearReaction* theElectronReaction = new G4ElectroNuclearReaction;
207 //theElectronNuclearProcess->RegisterMe(theElectronReaction);
208 ////theElectronNuclearProcess->BiasCrossSectionByFactor(100); // Very dangerous
209 //pmanager->AddDiscreteProcess(theElectronNuclearProcess);
210 // ........The NEW (native) usage of the CHIPS electro-nuclear reactions ............
211 pmanager->AddDiscreteProcess(new G4QCollision);
212
213 // Pure Electromagnetic Processes
214 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);// ElectronMultipleScattering
215 pmanager->AddProcess(new G4eIonisation(),-1,2,2); // Electron Ionisation
216 pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3); // Electron BremsStrahlung
217 }
218 else if (particleName == "e+") //positron
219 {
220 // Construct processes for positron
221 // ........The OLD (GHAR) usage of the CHIPS Photo-nuclear reactions .............
222 //G4PositronNuclearProcess* thePositronNuclearProcess = new G4PositronNuclearProcess;
223 //G4ElectroNuclearReaction* thePositronReaction = new G4ElectroNuclearReaction;
224 //thePositronNuclearProcess->RegisterMe(thePositronReaction);
225 ////thePositronNuclearProcess->BiasCrossSectionByFactor(100); // Very dangerous
226 //pmanager->AddDiscreteProcess(thePositronNuclearProcess);
227 // ........The NEW (native) usage of the CHIPS positron-nuclear reactions ...........
228 pmanager->AddDiscreteProcess(new G4QCollision);
229
230 // Pure Electromagnetic Processes
231 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);// PositronMultipleScattering
232 pmanager->AddProcess(new G4eIonisation(),-1,2,2); // Positron Ionisation
233 pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3); // Positron BremsStrahlung
234 G4eplusAnnihilation* theAnnihilation = new G4eplusAnnihilation;
235 pmanager->AddDiscreteProcess(theAnnihilation); // Positron Annihilation on Flight
236 pmanager->AddRestProcess(theAnnihilation); // Positron Annihilation at Rest
237 }
238 else if( particleName == "mu+" || particleName == "mu-" ) //muon of both signs
239 {
240 // Construct processes for muon+/-
241 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1); // Mu Multiple Scattering
242 pmanager->AddProcess(new G4MuIonisation(),-1,2,2); // Mu Ionization
243 pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3); // Mu Bremsstrahlung
244 pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4); // Mu Pair Production
245 // ........The NEW (native) usage of the CHIPS Muo-nuclear reactions .............
246 pmanager->AddDiscreteProcess(new G4QCollision);
247 // ........Negative lepton stopping
248 if(particleName == "mu-") pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
249 }
250 else if( particleName == "tau+" || particleName == "tau-" )
251 {
252 // Construct processes for muon+/-
253 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1); // Tau Multiple Scattering
254 // ........The NEW (native) usage of the CHIPS Muo-nuclear reactions .............
255 pmanager->AddDiscreteProcess(new G4QCollision);
256 // ........Negative lepton stopping
257 if(particleName == "tau-")pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
258 }
259 else if( particleName == "GenericIon" )
260 {
261 // Construct processes for ions
262 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
263 pmanager->AddProcess(new G4hIonisation(),-1,2,2);
264 }
265 else
266 {
267 if (particle->GetPDGCharge() && (particle->GetParticleName() != "chargedgeantino") &&
268 !particle->IsShortLived() )
269 {
270 // short lived particles except geantino
271 pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
272 pmanager->AddProcess(new G4hIonisation(),-1,2,2);
273 }
274 }
275
276 // ============================== Hadronic Interactions ===============================
277
278 if (particleName == "pi+")
279 {
280 pmanager->AddDiscreteProcess(theElasticProcess);
281 G4PionPlusInelasticProcess* theInelasticProcess =
282 new G4PionPlusInelasticProcess("inelastic");
283 G4LEPionPlusInelastic* theInelasticModel = new G4LEPionPlusInelastic;
284 theInelasticProcess->RegisterMe(theInelasticModel);
285 theInelasticProcess->RegisterMe(theTheoModel);
286 pmanager->AddDiscreteProcess(theInelasticProcess);
287 }
288 else if (particleName == "pi-")
289 {
290 pmanager->AddDiscreteProcess(theElasticProcess);
291 G4PionMinusInelasticProcess* theInelasticProcess =
292 new G4PionMinusInelasticProcess("inelastic");
293 G4LEPionMinusInelastic* theInelasticModel = new G4LEPionMinusInelastic;
294 theInelasticProcess->RegisterMe(theInelasticModel);
295 theInelasticProcess->RegisterMe(theTheoModel);
296 pmanager->AddDiscreteProcess(theInelasticProcess);
297 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
298 }
299 else if (particleName == "kaon+")
300 {
301 pmanager->AddDiscreteProcess(theElasticProcess);
302 G4KaonPlusInelasticProcess* theInelasticProcess =
303 new G4KaonPlusInelasticProcess("inelastic");
304 G4LEKaonPlusInelastic* theInelasticModel = new G4LEKaonPlusInelastic;
305 theInelasticProcess->RegisterMe(theInelasticModel);
306 new G4HEKaonPlusInelastic; // @@
307 theInelasticProcess->RegisterMe(theTheoModel);
308 pmanager->AddDiscreteProcess(theInelasticProcess);
309 }
310 else if (particleName == "kaon0S")
311 {
312 pmanager->AddDiscreteProcess(theElasticProcess);
313 G4KaonZeroSInelasticProcess* theInelasticProcess =
314 new G4KaonZeroSInelasticProcess("inelastic");
315 G4LEKaonZeroSInelastic* theInelasticModel = new G4LEKaonZeroSInelastic;
316 theInelasticProcess->RegisterMe(theInelasticModel);
317 theInelasticProcess->RegisterMe(theTheoModel);
318 pmanager->AddDiscreteProcess(theInelasticProcess);
319 }
320 else if (particleName == "kaon0L")
321 {
322 pmanager->AddDiscreteProcess(theElasticProcess);
323 G4KaonZeroLInelasticProcess* theInelasticProcess =
324 new G4KaonZeroLInelasticProcess("inelastic");
325 G4LEKaonZeroLInelastic* theInelasticModel = new G4LEKaonZeroLInelastic;
326 theInelasticProcess->RegisterMe(theInelasticModel);
327 new G4HEKaonZeroInelastic;
328 theInelasticProcess->RegisterMe(theTheoModel);
329 pmanager->AddDiscreteProcess(theInelasticProcess);
330 }
331 else if (particleName == "kaon-")
332 {
333 pmanager->AddDiscreteProcess(theElasticProcess);
334 G4KaonMinusInelasticProcess* theInelasticProcess =
335 new G4KaonMinusInelasticProcess("inelastic");
336 G4LEKaonMinusInelastic* theInelasticModel = new G4LEKaonMinusInelastic;
337 theInelasticProcess->RegisterMe(theInelasticModel);
338 theInelasticProcess->RegisterMe(theTheoModel);
339 pmanager->AddDiscreteProcess(theInelasticProcess);
340 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
341 }
342 else if (particleName == "proton")
343 {
344 pmanager->AddDiscreteProcess(theElasticProcess);
345 G4ProtonInelasticProcess* theInelasticProcess =
346 new G4ProtonInelasticProcess("inelastic");
347 G4LEProtonInelastic* theInelasticModel = new G4LEProtonInelastic;
348 theInelasticProcess->RegisterMe(theInelasticModel);
349 theInelasticProcess->RegisterMe(theTheoModel);
350 pmanager->AddDiscreteProcess(theInelasticProcess);
351 }
352 else if (particleName == "anti_proton")
353 {
354 pmanager->AddDiscreteProcess(theElasticProcess);
355 G4AntiProtonInelasticProcess* theInelasticProcess =
356 new G4AntiProtonInelasticProcess("inelastic");
357 G4LEAntiProtonInelastic* theInelasticModel = new G4LEAntiProtonInelastic;
358 theInelasticProcess->RegisterMe(theInelasticModel);
359 G4HEAntiProtonInelastic* theHEInelasticModel = new G4HEAntiProtonInelastic;
360 theInelasticProcess->RegisterMe(theHEInelasticModel);
361 pmanager->AddDiscreteProcess(theInelasticProcess);
362 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
363 }
364 else if (particleName == "neutron")
365 {
366 // elastic scattering
367 G4LElastic* theElasticModel1 = new G4LElastic;
368 theElasticProcess1->RegisterMe(theElasticModel1);
369 pmanager->AddDiscreteProcess(theElasticProcess1);
370 // inelastic scattering
371 G4NeutronInelasticProcess* theInelasticProcess =
372 new G4NeutronInelasticProcess("inelastic");
373 G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
374 theInelasticProcess->RegisterMe(theInelasticModel);
375 theInelasticProcess->RegisterMe(theTheoModel);
376 pmanager->AddDiscreteProcess(theInelasticProcess);
377 // fission
378 G4HadronFissionProcess* theFissionProcess = new G4HadronFissionProcess;
379 G4LFission* theFissionModel = new G4LFission;
380 theFissionProcess->RegisterMe(theFissionModel);
381 pmanager->AddDiscreteProcess(theFissionProcess);
382 // capture
383 G4HadronCaptureProcess* theCaptureProcess = new G4HadronCaptureProcess;
384 G4LCapture* theCaptureModel = new G4LCapture;
385 theCaptureProcess->RegisterMe(theCaptureModel);
386 pmanager->AddDiscreteProcess(theCaptureProcess);
387 }
388 else if (particleName == "anti_neutron")
389 {
390 pmanager->AddDiscreteProcess(theElasticProcess);
391 G4AntiNeutronInelasticProcess* theInelasticProcess =
392 new G4AntiNeutronInelasticProcess("inelastic");
393 G4LEAntiNeutronInelastic* theInelasticModel = new G4LEAntiNeutronInelastic;
394 theInelasticProcess->RegisterMe(theInelasticModel);
395 G4HEAntiNeutronInelastic* theHEInelasticModel = new G4HEAntiNeutronInelastic;
396 theInelasticProcess->RegisterMe(theHEInelasticModel);
397 pmanager->AddDiscreteProcess(theInelasticProcess);
398 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
399 }
400 else if (particleName == "lambda")
401 {
402 pmanager->AddDiscreteProcess(theElasticProcess);
403 G4LambdaInelasticProcess* theInelasticProcess =
404 new G4LambdaInelasticProcess("inelastic");
405 G4LELambdaInelastic* theInelasticModel = new G4LELambdaInelastic;
406 theInelasticProcess->RegisterMe(theInelasticModel);
407 G4HELambdaInelastic* theHEInelasticModel = new G4HELambdaInelastic;
408 theInelasticProcess->RegisterMe(theHEInelasticModel);
409 pmanager->AddDiscreteProcess(theInelasticProcess);
410 }
411 else if (particleName == "anti_lambda")
412 {
413 pmanager->AddDiscreteProcess(theElasticProcess);
414 G4AntiLambdaInelasticProcess* theInelasticProcess =
415 new G4AntiLambdaInelasticProcess("inelastic");
416 G4LEAntiLambdaInelastic* theInelasticModel = new G4LEAntiLambdaInelastic;
417 theInelasticProcess->RegisterMe(theInelasticModel);
418 G4HEAntiLambdaInelastic* theHEInelasticModel = new G4HEAntiLambdaInelastic;
419 theInelasticProcess->RegisterMe(theHEInelasticModel);
420 pmanager->AddDiscreteProcess(theInelasticProcess);
421 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
422 }
423 else if (particleName == "sigma+")
424 {
425 pmanager->AddDiscreteProcess(theElasticProcess);
426 G4SigmaPlusInelasticProcess* theInelasticProcess =
427 new G4SigmaPlusInelasticProcess("inelastic");
428 G4LESigmaPlusInelastic* theInelasticModel = new G4LESigmaPlusInelastic;
429 theInelasticProcess->RegisterMe(theInelasticModel);
430 G4HESigmaPlusInelastic* theHEInelasticModel = new G4HESigmaPlusInelastic;
431 theInelasticProcess->RegisterMe(theHEInelasticModel);
432 pmanager->AddDiscreteProcess(theInelasticProcess);
433 }
434 else if (particleName == "sigma-")
435 {
436 pmanager->AddDiscreteProcess(theElasticProcess);
437 G4SigmaMinusInelasticProcess* theInelasticProcess =
438 new G4SigmaMinusInelasticProcess("inelastic");
439 G4LESigmaMinusInelastic* theInelasticModel = new G4LESigmaMinusInelastic;
440 theInelasticProcess->RegisterMe(theInelasticModel);
441 G4HESigmaMinusInelastic* theHEInelasticModel = new G4HESigmaMinusInelastic;
442 theInelasticProcess->RegisterMe(theHEInelasticModel);
443 pmanager->AddDiscreteProcess(theInelasticProcess);
444 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
445 }
446 else if (particleName == "anti_sigma+")
447 {
448 pmanager->AddDiscreteProcess(theElasticProcess);
449 G4AntiSigmaPlusInelasticProcess* theInelasticProcess =
450 new G4AntiSigmaPlusInelasticProcess("inelastic");
451 G4LEAntiSigmaPlusInelastic* theInelasticModel = new G4LEAntiSigmaPlusInelastic;
452 theInelasticProcess->RegisterMe(theInelasticModel);
453 G4HEAntiSigmaPlusInelastic* theHEInelasticModel = new G4HEAntiSigmaPlusInelastic;
454 theInelasticProcess->RegisterMe(theHEInelasticModel);
455 pmanager->AddDiscreteProcess(theInelasticProcess);
456 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
457 }
458 else if (particleName == "anti_sigma-")
459 {
460 pmanager->AddDiscreteProcess(theElasticProcess);
461 G4AntiSigmaMinusInelasticProcess* theInelasticProcess =
462 new G4AntiSigmaMinusInelasticProcess("inelastic");
463 G4LEAntiSigmaMinusInelastic* theInelasticModel = new G4LEAntiSigmaMinusInelastic;
464 theInelasticProcess->RegisterMe(theInelasticModel);
465 G4HEAntiSigmaMinusInelastic* theHEInelasticModel = new G4HEAntiSigmaMinusInelastic;
466 theInelasticProcess->RegisterMe(theHEInelasticModel);
467 pmanager->AddDiscreteProcess(theInelasticProcess);
468 }
469 else if (particleName == "xi0")
470 {
471 pmanager->AddDiscreteProcess(theElasticProcess);
472 G4XiZeroInelasticProcess* theInelasticProcess =
473 new G4XiZeroInelasticProcess("inelastic");
474 G4LEXiZeroInelastic* theInelasticModel = new G4LEXiZeroInelastic;
475 theInelasticProcess->RegisterMe(theInelasticModel);
476 G4HEXiZeroInelastic* theHEInelasticModel = new G4HEXiZeroInelastic;
477 theInelasticProcess->RegisterMe(theHEInelasticModel);
478 pmanager->AddDiscreteProcess(theInelasticProcess);
479 }
480 else if (particleName == "xi-")
481 {
482 pmanager->AddDiscreteProcess(theElasticProcess);
483 G4XiMinusInelasticProcess* theInelasticProcess =
484 new G4XiMinusInelasticProcess("inelastic");
485 G4LEXiMinusInelastic* theInelasticModel = new G4LEXiMinusInelastic;
486 theInelasticProcess->RegisterMe(theInelasticModel);
487 G4HEXiMinusInelastic* theHEInelasticModel = new G4HEXiMinusInelastic;
488 theInelasticProcess->RegisterMe(theHEInelasticModel);
489 pmanager->AddDiscreteProcess(theInelasticProcess);
490 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
491 }
492 else if (particleName == "anti_xi0")
493 {
494 pmanager->AddDiscreteProcess(theElasticProcess);
495 G4AntiXiZeroInelasticProcess* theInelasticProcess =
496 new G4AntiXiZeroInelasticProcess("inelastic");
497 G4LEAntiXiZeroInelastic* theInelasticModel = new G4LEAntiXiZeroInelastic;
498 theInelasticProcess->RegisterMe(theInelasticModel);
499 G4HEAntiXiZeroInelastic* theHEInelasticModel = new G4HEAntiXiZeroInelastic;
500 theInelasticProcess->RegisterMe(theHEInelasticModel);
501 pmanager->AddDiscreteProcess(theInelasticProcess);
502 }
503 else if (particleName == "anti_xi-")
504 {
505 pmanager->AddDiscreteProcess(theElasticProcess);
506 G4AntiXiMinusInelasticProcess* theInelasticProcess =
507 new G4AntiXiMinusInelasticProcess("inelastic");
508 G4LEAntiXiMinusInelastic* theInelasticModel = new G4LEAntiXiMinusInelastic;
509 theInelasticProcess->RegisterMe(theInelasticModel);
510 G4HEAntiXiMinusInelastic* theHEInelasticModel = new G4HEAntiXiMinusInelastic;
511 theInelasticProcess->RegisterMe(theHEInelasticModel);
512 pmanager->AddDiscreteProcess(theInelasticProcess);
513 }
514 else if (particleName == "deuteron")
515 {
516 pmanager->AddDiscreteProcess(theElasticProcess);
517 G4DeuteronInelasticProcess* theInelasticProcess =
518 new G4DeuteronInelasticProcess("inelastic");
519 G4LEDeuteronInelastic* theInelasticModel = new G4LEDeuteronInelastic;
520 theInelasticProcess->RegisterMe(theInelasticModel);
521 pmanager->AddDiscreteProcess(theInelasticProcess);
522 }
523 else if (particleName == "triton")
524 {
525 pmanager->AddDiscreteProcess(theElasticProcess);
526 G4TritonInelasticProcess* theInelasticProcess =
527 new G4TritonInelasticProcess("inelastic");
528 G4LETritonInelastic* theInelasticModel = new G4LETritonInelastic;
529 theInelasticProcess->RegisterMe(theInelasticModel);
530 pmanager->AddDiscreteProcess(theInelasticProcess);
531 }
532 else if (particleName == "alpha")
533 {
534 pmanager->AddDiscreteProcess(theElasticProcess);
535 G4AlphaInelasticProcess* theInelasticProcess =
536 new G4AlphaInelasticProcess("inelastic");
537 G4LEAlphaInelastic* theInelasticModel = new G4LEAlphaInelastic;
538 theInelasticProcess->RegisterMe(theInelasticModel);
539 pmanager->AddDiscreteProcess(theInelasticProcess);
540 }
541 else if (particleName == "omega-")
542 {
543 pmanager->AddDiscreteProcess(theElasticProcess);
544 G4OmegaMinusInelasticProcess* theInelasticProcess =
545 new G4OmegaMinusInelasticProcess("inelastic");
546 G4LEOmegaMinusInelastic* theInelasticModel = new G4LEOmegaMinusInelastic;
547 theInelasticProcess->RegisterMe(theInelasticModel);
548 G4HEOmegaMinusInelastic* theHEInelasticModel = new G4HEOmegaMinusInelastic;
549 theInelasticProcess->RegisterMe(theHEInelasticModel);
550 pmanager->AddDiscreteProcess(theInelasticProcess);
551 pmanager->AddRestProcess(new G4QCaptureAtRest, ordDefault);
552 }
553 else if (particleName == "anti_omega-")
554 {
555 pmanager->AddDiscreteProcess(theElasticProcess);
556 G4AntiOmegaMinusInelasticProcess* theInelasticProcess =
557 new G4AntiOmegaMinusInelasticProcess("inelastic");
558 G4LEAntiOmegaMinusInelastic* theInelasticModel = new G4LEAntiOmegaMinusInelastic;
559 theInelasticProcess->RegisterMe(theInelasticModel);
560 G4HEAntiOmegaMinusInelastic* theHEInelasticModel = new G4HEAntiOmegaMinusInelastic;
561 theInelasticProcess->RegisterMe(theHEInelasticModel);
562 pmanager->AddDiscreteProcess(theInelasticProcess);
563 }
564 }
565}
566
567void PhotInPhysicsList::SetCuts()
568{
569 if(verboseLevel>0) G4cout<<"PhotInPhysicsList::SetCuts: default cut length : "
570 <<G4BestUnit(defaultCutValue,"Length")<<G4endl;
571 // These values are used as the default production thresholds for the world volume.
572 SetCutsWithDefault();
573
574 // Production thresholds for detector regions
575
576 G4double fact = 1.; // Multiplicative factor for default cuts
577 for(G4int i=0; i< PhotInNumSections; i++)
578 {
579 G4Region* reg = G4RegionStore::GetInstance()-> GetRegion(PhotInRegName[i]);
580 G4ProductionCuts* cuts = new G4ProductionCuts;
581 cuts->SetProductionCut(defaultCutValue*fact);
582 reg->SetProductionCuts(cuts);
583 fact *= 10.; // @@ Increment the multiplicative factor by order of magnitude
584 }
585}
586
587
Note: See TracBrowser for help on using the repository browser.