source: trunk/examples/advanced/Tiara/source/tiara/src/TiaraPhysicsList.cc @ 1002

Last change on this file since 1002 was 807, checked in by garnier, 16 years ago

update

File size: 26.3 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: TiaraPhysicsList.cc,v 1.6 2006/06/29 15:45:20 gunter Exp $
28// GEANT4 tag $Name:  $
29//
30
31#include "globals.hh"
32#include <iomanip>                
33
34#include "TiaraPhysicsList.hh"
35#include "G4ParticleDefinition.hh"
36#include "G4ParticleWithCuts.hh"
37#include "G4ProcessManager.hh"
38#include "G4ProcessVector.hh"
39#include "G4ParticleTypes.hh"
40#include "G4ParticleTable.hh"
41#include "G4BosonConstructor.hh"
42#include "G4LeptonConstructor.hh"
43#include "G4MesonConstructor.hh"
44#include "G4BaryonConstructor.hh"
45#include "G4IonConstructor.hh"
46#include "G4ShortLivedConstructor.hh"
47#include "G4Material.hh"
48#include "G4MaterialTable.hh"
49
50TiaraPhysicsList::TiaraPhysicsList() : G4VUserPhysicsList()
51{
52  SetVerboseLevel(1);
53}
54
55TiaraPhysicsList::~TiaraPhysicsList()
56{}
57
58void TiaraPhysicsList::ConstructParticle()
59{
60  // In this method, static member functions should be called
61  // for all particles which you want to use.
62  // This ensures that objects of these particle types will be
63  // created in the program.
64
65  ConstructAllBosons();
66  ConstructAllLeptons();
67  ConstructAllMesons();
68  ConstructAllBaryons();
69  ConstructAllIons();
70  ConstructAllShortLiveds();
71}
72
73void TiaraPhysicsList::ConstructAllBosons()
74{
75  // Construct all bosons
76  G4BosonConstructor pConstructor;
77  pConstructor.ConstructParticle();
78}
79
80void TiaraPhysicsList::ConstructAllLeptons()
81{
82  // Construct all leptons
83  G4LeptonConstructor pConstructor;
84  pConstructor.ConstructParticle();
85}
86
87void TiaraPhysicsList::ConstructAllMesons()
88{
89  //  Construct all mesons
90  G4MesonConstructor pConstructor;
91  pConstructor.ConstructParticle();
92}
93
94void TiaraPhysicsList::ConstructAllBaryons()
95{
96  //  Construct all barions
97  G4BaryonConstructor pConstructor;
98  pConstructor.ConstructParticle();
99}
100
101void TiaraPhysicsList::ConstructAllIons()
102{
103  //  Construct light ions
104  G4IonConstructor pConstructor;
105  pConstructor.ConstructParticle(); 
106}
107
108void TiaraPhysicsList::ConstructAllShortLiveds()
109{
110  //  Construct  resonaces and quarks
111  G4ShortLivedConstructor pConstructor;
112  pConstructor.ConstructParticle(); 
113}
114
115void TiaraPhysicsList::ConstructProcess()
116{
117  AddTransportation();
118  ConstructEM();
119  ConstructLeptHad();
120  ConstructHad();
121  ConstructGeneral();
122}
123
124#include "G4ComptonScattering.hh"
125#include "G4GammaConversion.hh"
126#include "G4PhotoElectricEffect.hh"
127
128#include "G4MultipleScattering.hh"
129
130#include "G4eIonisation.hh"
131#include "G4eBremsstrahlung.hh"
132#include "G4eplusAnnihilation.hh"
133
134#include "G4MuIonisation.hh"
135#include "G4MuBremsstrahlung.hh"
136#include "G4MuPairProduction.hh"
137
138#include "G4hIonisation.hh"
139#include "G4ionIonisation.hh"
140
141void TiaraPhysicsList::ConstructEM()
142{
143  theParticleIterator->reset();
144  while( (*theParticleIterator)() ){
145    G4ParticleDefinition* particle = theParticleIterator->value();
146    G4ProcessManager* pmanager = particle->GetProcessManager();
147    G4String particleName = particle->GetParticleName();
148     
149    if (particleName == "gamma") {
150    // gamma
151      // Construct processes for gamma
152      pmanager->AddDiscreteProcess(new G4GammaConversion());
153      pmanager->AddDiscreteProcess(new G4ComptonScattering());     
154      pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
155
156    } else if (particleName == "e-") {
157    //electron
158      // Construct processes for electron
159      pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
160      pmanager->AddProcess(new G4eIonisation(),-1,2,2);
161      pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
162 
163    } else if (particleName == "e+") {
164    //positron
165      // Construct processes for positron
166     pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
167     
168     pmanager->AddProcess(new G4eIonisation(),-1,2,2);
169     pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);     
170     pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
171 
172    } else if( particleName == "mu+" || 
173               particleName == "mu-"    ) {
174    //muon 
175     // Construct processes for muon+
176     pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
177     pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
178     pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
179     pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);       
180     
181    } else if( particleName == "GenericIon" ) {
182      pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
183      pmanager->AddProcess(new G4ionIonisation(),-1,2,2); 
184    } else { 
185      if ((particle->GetPDGCharge() != 0.0) && 
186          (particle->GetParticleName() != "chargedgeantino")&&
187          (!particle->IsShortLived()) ) {
188     // all others charged particles except geantino
189       pmanager->AddProcess(new G4MultipleScattering(),-1,1,1);
190       pmanager->AddProcess(new G4hIonisation(),-1,2,2);       
191     }
192    }
193  }
194}
195
196// Hadron Processes
197
198#include "G4HadronElasticProcess.hh"
199#include "G4HadronFissionProcess.hh"
200#include "G4HadronCaptureProcess.hh"
201
202#include "G4PionPlusInelasticProcess.hh"
203#include "G4PionMinusInelasticProcess.hh"
204#include "G4KaonPlusInelasticProcess.hh"
205#include "G4KaonZeroSInelasticProcess.hh"
206#include "G4KaonZeroLInelasticProcess.hh"
207#include "G4KaonMinusInelasticProcess.hh"
208#include "G4ProtonInelasticProcess.hh"
209#include "G4AntiProtonInelasticProcess.hh"
210#include "G4NeutronInelasticProcess.hh"
211#include "G4AntiNeutronInelasticProcess.hh"
212#include "G4LambdaInelasticProcess.hh"
213#include "G4AntiLambdaInelasticProcess.hh"
214#include "G4SigmaPlusInelasticProcess.hh"
215#include "G4SigmaMinusInelasticProcess.hh"
216#include "G4AntiSigmaPlusInelasticProcess.hh"
217#include "G4AntiSigmaMinusInelasticProcess.hh"
218#include "G4XiZeroInelasticProcess.hh"
219#include "G4XiMinusInelasticProcess.hh"
220#include "G4AntiXiZeroInelasticProcess.hh"
221#include "G4AntiXiMinusInelasticProcess.hh"
222#include "G4DeuteronInelasticProcess.hh"
223#include "G4TritonInelasticProcess.hh"
224#include "G4AlphaInelasticProcess.hh"
225#include "G4OmegaMinusInelasticProcess.hh"
226#include "G4AntiOmegaMinusInelasticProcess.hh"
227
228// Low-energy Models
229
230#include "G4LElastic.hh"
231#include "G4LFission.hh"
232#include "G4LCapture.hh"
233
234#include "G4LEPionPlusInelastic.hh"
235#include "G4LEPionMinusInelastic.hh"
236#include "G4LEKaonPlusInelastic.hh"
237#include "G4LEKaonZeroSInelastic.hh"
238#include "G4LEKaonZeroLInelastic.hh"
239#include "G4LEKaonMinusInelastic.hh"
240#include "G4LEProtonInelastic.hh"
241#include "G4LEAntiProtonInelastic.hh"
242#include "G4LENeutronInelastic.hh"
243#include "G4LEAntiNeutronInelastic.hh"
244#include "G4LELambdaInelastic.hh"
245#include "G4LEAntiLambdaInelastic.hh"
246#include "G4LESigmaPlusInelastic.hh"
247#include "G4LESigmaMinusInelastic.hh"
248#include "G4LEAntiSigmaPlusInelastic.hh"
249#include "G4LEAntiSigmaMinusInelastic.hh"
250#include "G4LEXiZeroInelastic.hh"
251#include "G4LEXiMinusInelastic.hh"
252#include "G4LEAntiXiZeroInelastic.hh"
253#include "G4LEAntiXiMinusInelastic.hh"
254#include "G4LEDeuteronInelastic.hh"
255#include "G4LETritonInelastic.hh"
256#include "G4LEAlphaInelastic.hh"
257#include "G4LEOmegaMinusInelastic.hh"
258#include "G4LEAntiOmegaMinusInelastic.hh"
259
260// -- generator models
261#include "G4TheoFSGenerator.hh"
262#include "G4ExcitationHandler.hh"
263#include "G4Evaporation.hh"
264#include "G4CompetitiveFission.hh"
265#include "G4FermiBreakUp.hh"
266#include "G4StatMF.hh"
267#include "G4GeneratorPrecompoundInterface.hh"
268#include "G4Fancy3DNucleus.hh"
269#include "G4LEProtonInelastic.hh"
270#include "G4StringModel.hh"
271#include "G4PreCompoundModel.hh"
272#include "G4FTFModel.hh"
273#include "G4QGSMFragmentation.hh"
274#include "G4ExcitedStringDecay.hh"
275
276//
277// ConstructHad()
278//
279// Makes discrete physics processes for the hadrons, at present limited
280// to those particles with GHEISHA interactions (INTRC > 0).
281// The processes are: Elastic scattering, Inelastic scattering,
282// Fission (for neutron only), and Capture (neutron).
283//
284// F.W.Jones  06-JUL-1998
285//
286
287void TiaraPhysicsList::ConstructHad()
288{
289    // this will be the model class for high energies
290    G4TheoFSGenerator * theTheoModel = new G4TheoFSGenerator;
291       
292    // all models for treatment of thermal nucleus
293    G4Evaporation * theEvaporation = new G4Evaporation;
294    G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
295    G4StatMF * theMF = new G4StatMF;
296
297    // Evaporation logic
298    G4ExcitationHandler * theHandler = new G4ExcitationHandler;
299        theHandler->SetEvaporation(theEvaporation);
300        theHandler->SetFermiModel(theFermiBreakUp);
301        theHandler->SetMultiFragmentation(theMF);
302        theHandler->SetMaxAandZForFermiBreakUp(12, 6);
303        theHandler->SetMinEForMultiFrag(3*MeV);
304       
305    // Pre equilibrium stage
306    G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel(theHandler);
307
308   
309    // a no-cascade generator-precompound interaface
310    G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface;
311            theCascade->SetDeExcitation(thePreEquilib); 
312       
313    // here come the high energy parts
314    // the string model; still not quite according to design - Explicite use of the forseen interfaces
315    // will be tested and documented in this program by beta-02 at latest.
316    G4VPartonStringModel * theStringModel;
317    theStringModel = new G4FTFModel;
318    theTheoModel->SetTransport(theCascade);
319    theTheoModel->SetHighEnergyGenerator(theStringModel);
320    theTheoModel->SetMinEnergy(19*GeV);
321    theTheoModel->SetMaxEnergy(100*TeV);
322
323      G4VLongitudinalStringDecay * theFragmentation = new G4QGSMFragmentation;
324      G4ExcitedStringDecay * theStringDecay = new G4ExcitedStringDecay(theFragmentation);
325      theStringModel->SetFragmentationModel(theStringDecay);
326
327// done with the generator model (most of the above is also available as default)
328   G4HadronElasticProcess* theElasticProcess = 
329                                    new G4HadronElasticProcess;
330   G4LElastic* theElasticModel = new G4LElastic;
331   theElasticProcess->RegisterMe(theElasticModel);
332   G4HadronElasticProcess* theElasticProcess1 = 
333                                    new G4HadronElasticProcess;
334   theParticleIterator->reset();
335   while ((*theParticleIterator)()) {
336      G4ParticleDefinition* particle = theParticleIterator->value();
337      G4ProcessManager* pmanager = particle->GetProcessManager();
338      G4String particleName = particle->GetParticleName();
339     
340      if (particleName == "pi+") {
341         pmanager->AddDiscreteProcess(theElasticProcess);
342         G4PionPlusInelasticProcess* theInelasticProcess = 
343                                new G4PionPlusInelasticProcess("inelastic");
344         G4LEPionPlusInelastic* theInelasticModel = 
345                                new G4LEPionPlusInelastic;
346         theInelasticProcess->RegisterMe(theInelasticModel);
347         theInelasticProcess->RegisterMe(theTheoModel);
348         pmanager->AddDiscreteProcess(theInelasticProcess);
349      }
350      else if (particleName == "pi-") {
351         pmanager->AddDiscreteProcess(theElasticProcess);
352         G4PionMinusInelasticProcess* theInelasticProcess = 
353                                new G4PionMinusInelasticProcess("inelastic");
354         G4LEPionMinusInelastic* theInelasticModel = 
355                                new G4LEPionMinusInelastic;
356         theInelasticProcess->RegisterMe(theInelasticModel);
357         theInelasticProcess->RegisterMe(theTheoModel);
358         pmanager->AddDiscreteProcess(theInelasticProcess);
359      }
360      else if (particleName == "kaon+") {
361         pmanager->AddDiscreteProcess(theElasticProcess);
362         G4KaonPlusInelasticProcess* theInelasticProcess = 
363                                  new G4KaonPlusInelasticProcess("inelastic");
364         G4LEKaonPlusInelastic* theInelasticModel = new G4LEKaonPlusInelastic;
365         theInelasticProcess->RegisterMe(theInelasticModel);
366         theInelasticProcess->RegisterMe(theTheoModel);
367         pmanager->AddDiscreteProcess(theInelasticProcess);
368      }
369      else if (particleName == "kaon0S") {
370         pmanager->AddDiscreteProcess(theElasticProcess);
371         G4KaonZeroSInelasticProcess* theInelasticProcess = 
372                             new G4KaonZeroSInelasticProcess("inelastic");
373         G4LEKaonZeroSInelastic* theInelasticModel = 
374                             new G4LEKaonZeroSInelastic;
375         theInelasticProcess->RegisterMe(theInelasticModel);
376         theInelasticProcess->RegisterMe(theTheoModel);
377         pmanager->AddDiscreteProcess(theInelasticProcess);
378      }
379      else if (particleName == "kaon0L") {
380         pmanager->AddDiscreteProcess(theElasticProcess);
381         G4KaonZeroLInelasticProcess* theInelasticProcess = 
382                             new G4KaonZeroLInelasticProcess("inelastic");
383         G4LEKaonZeroLInelastic* theInelasticModel = 
384                             new G4LEKaonZeroLInelastic;
385         theInelasticProcess->RegisterMe(theInelasticModel);
386         theInelasticProcess->RegisterMe(theTheoModel);
387         pmanager->AddDiscreteProcess(theInelasticProcess);
388      }
389      else if (particleName == "kaon-") {
390         pmanager->AddDiscreteProcess(theElasticProcess);
391         G4KaonMinusInelasticProcess* theInelasticProcess = 
392                                 new G4KaonMinusInelasticProcess("inelastic");
393         G4LEKaonMinusInelastic* theInelasticModel = 
394                                 new G4LEKaonMinusInelastic;
395         theInelasticProcess->RegisterMe(theInelasticModel);
396         theInelasticProcess->RegisterMe(theTheoModel);
397         pmanager->AddDiscreteProcess(theInelasticProcess);
398      }
399      else if (particleName == "proton") {
400         pmanager->AddDiscreteProcess(theElasticProcess);
401         G4ProtonInelasticProcess* theInelasticProcess = 
402                                    new G4ProtonInelasticProcess("inelastic");
403         G4LEProtonInelastic* theInelasticModel = new G4LEProtonInelastic;
404         theInelasticProcess->RegisterMe(theInelasticModel);
405         theInelasticProcess->RegisterMe(theTheoModel);
406         pmanager->AddDiscreteProcess(theInelasticProcess);
407      }
408      else if (particleName == "anti_proton") {
409         pmanager->AddDiscreteProcess(theElasticProcess);
410         G4AntiProtonInelasticProcess* theInelasticProcess = 
411                                new G4AntiProtonInelasticProcess("inelastic");
412         G4LEAntiProtonInelastic* theInelasticModel = 
413                                new G4LEAntiProtonInelastic;
414         theInelasticProcess->RegisterMe(theInelasticModel);
415         theInelasticProcess->RegisterMe(theTheoModel);
416         pmanager->AddDiscreteProcess(theInelasticProcess);
417      }
418      else if (particleName == "neutron") {
419         
420          // elastic scattering
421         G4LElastic* theElasticModel1 = new G4LElastic;
422         theElasticProcess1->RegisterMe(theElasticModel1);
423         pmanager->AddDiscreteProcess(theElasticProcess1);
424          // inelastic scattering
425         G4NeutronInelasticProcess* theInelasticProcess = 
426                                    new G4NeutronInelasticProcess("inelastic");
427         G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
428         theInelasticProcess->RegisterMe(theInelasticModel);
429         theInelasticProcess->RegisterMe(theTheoModel);
430         pmanager->AddDiscreteProcess(theInelasticProcess);
431          // fission
432         G4HadronFissionProcess* theFissionProcess =
433                                    new G4HadronFissionProcess;
434         G4LFission* theFissionModel = new G4LFission;
435         theFissionProcess->RegisterMe(theFissionModel);
436         pmanager->AddDiscreteProcess(theFissionProcess);
437         // capture
438         G4HadronCaptureProcess* theCaptureProcess =
439                                    new G4HadronCaptureProcess;
440         G4LCapture* theCaptureModel = new G4LCapture;
441         theCaptureProcess->RegisterMe(theCaptureModel);
442         pmanager->AddDiscreteProcess(theCaptureProcess);
443      } 
444      else if (particleName == "anti_neutron") {
445         pmanager->AddDiscreteProcess(theElasticProcess);
446         G4AntiNeutronInelasticProcess* theInelasticProcess = 
447                               new G4AntiNeutronInelasticProcess("inelastic");
448         G4LEAntiNeutronInelastic* theInelasticModel = 
449                               new G4LEAntiNeutronInelastic;
450         theInelasticProcess->RegisterMe(theInelasticModel);
451         theInelasticProcess->RegisterMe(theTheoModel);
452         pmanager->AddDiscreteProcess(theInelasticProcess);
453      }
454      else if (particleName == "lambda") {
455         pmanager->AddDiscreteProcess(theElasticProcess);
456         G4LambdaInelasticProcess* theInelasticProcess = 
457                                    new G4LambdaInelasticProcess("inelastic");
458         G4LELambdaInelastic* theInelasticModel = new G4LELambdaInelastic;
459         theInelasticProcess->RegisterMe(theInelasticModel);
460         theInelasticProcess->RegisterMe(theTheoModel);
461         pmanager->AddDiscreteProcess(theInelasticProcess);
462      }
463      else if (particleName == "anti_lambda") {
464         pmanager->AddDiscreteProcess(theElasticProcess);
465         G4AntiLambdaInelasticProcess* theInelasticProcess = 
466                                new G4AntiLambdaInelasticProcess("inelastic");
467         G4LEAntiLambdaInelastic* theInelasticModel = 
468                                new G4LEAntiLambdaInelastic;
469         theInelasticProcess->RegisterMe(theInelasticModel);
470         theInelasticProcess->RegisterMe(theTheoModel);
471         pmanager->AddDiscreteProcess(theInelasticProcess);
472      }
473      else if (particleName == "sigma+") {
474         pmanager->AddDiscreteProcess(theElasticProcess);
475         G4SigmaPlusInelasticProcess* theInelasticProcess = 
476                                 new G4SigmaPlusInelasticProcess("inelastic");
477         G4LESigmaPlusInelastic* theInelasticModel = 
478                                 new G4LESigmaPlusInelastic;
479         theInelasticProcess->RegisterMe(theInelasticModel);
480         theInelasticProcess->RegisterMe(theTheoModel);
481         pmanager->AddDiscreteProcess(theInelasticProcess);
482      }
483      else if (particleName == "sigma-") {
484         pmanager->AddDiscreteProcess(theElasticProcess);
485         G4SigmaMinusInelasticProcess* theInelasticProcess = 
486                                 new G4SigmaMinusInelasticProcess("inelastic");
487         G4LESigmaMinusInelastic* theInelasticModel = 
488                                 new G4LESigmaMinusInelastic;
489         theInelasticProcess->RegisterMe(theInelasticModel);
490         theInelasticProcess->RegisterMe(theTheoModel);
491         pmanager->AddDiscreteProcess(theInelasticProcess);
492      }
493      else if (particleName == "anti_sigma+") {
494         pmanager->AddDiscreteProcess(theElasticProcess);
495         G4AntiSigmaPlusInelasticProcess* theInelasticProcess = 
496                             new G4AntiSigmaPlusInelasticProcess("inelastic");
497         G4LEAntiSigmaPlusInelastic* theInelasticModel = 
498                                 new G4LEAntiSigmaPlusInelastic;
499         theInelasticProcess->RegisterMe(theInelasticModel);
500         theInelasticProcess->RegisterMe(theTheoModel);
501         pmanager->AddDiscreteProcess(theInelasticProcess);
502      }
503      else if (particleName == "anti_sigma-") {
504         pmanager->AddDiscreteProcess(theElasticProcess);
505         G4AntiSigmaMinusInelasticProcess* theInelasticProcess = 
506                            new G4AntiSigmaMinusInelasticProcess("inelastic");
507         G4LEAntiSigmaMinusInelastic* theInelasticModel = 
508                                 new G4LEAntiSigmaMinusInelastic;
509         theInelasticProcess->RegisterMe(theInelasticModel);
510         theInelasticProcess->RegisterMe(theTheoModel);
511         pmanager->AddDiscreteProcess(theInelasticProcess);
512      }
513      else if (particleName == "xi0") {
514         pmanager->AddDiscreteProcess(theElasticProcess);
515         G4XiZeroInelasticProcess* theInelasticProcess = 
516                            new G4XiZeroInelasticProcess("inelastic");
517         G4LEXiZeroInelastic* theInelasticModel = 
518                                 new G4LEXiZeroInelastic;
519         theInelasticProcess->RegisterMe(theInelasticModel);
520         theInelasticProcess->RegisterMe(theTheoModel);
521         pmanager->AddDiscreteProcess(theInelasticProcess);
522      }
523      else if (particleName == "xi-") {
524         pmanager->AddDiscreteProcess(theElasticProcess);
525         G4XiMinusInelasticProcess* theInelasticProcess = 
526                            new G4XiMinusInelasticProcess("inelastic");
527         G4LEXiMinusInelastic* theInelasticModel = 
528                                 new G4LEXiMinusInelastic;
529         theInelasticProcess->RegisterMe(theInelasticModel);
530         theInelasticProcess->RegisterMe(theTheoModel);
531         pmanager->AddDiscreteProcess(theInelasticProcess);
532      }
533      else if (particleName == "anti_xi0") {
534         pmanager->AddDiscreteProcess(theElasticProcess);
535         G4AntiXiZeroInelasticProcess* theInelasticProcess = 
536                            new G4AntiXiZeroInelasticProcess("inelastic");
537         G4LEAntiXiZeroInelastic* theInelasticModel = 
538                                 new G4LEAntiXiZeroInelastic;
539         theInelasticProcess->RegisterMe(theInelasticModel);
540         theInelasticProcess->RegisterMe(theTheoModel);
541         pmanager->AddDiscreteProcess(theInelasticProcess);
542      }
543      else if (particleName == "anti_xi-") {
544         pmanager->AddDiscreteProcess(theElasticProcess);
545         G4AntiXiMinusInelasticProcess* theInelasticProcess = 
546                            new G4AntiXiMinusInelasticProcess("inelastic");
547         G4LEAntiXiMinusInelastic* theInelasticModel = 
548                                 new G4LEAntiXiMinusInelastic;
549         theInelasticProcess->RegisterMe(theInelasticModel);
550         theInelasticProcess->RegisterMe(theTheoModel);
551         pmanager->AddDiscreteProcess(theInelasticProcess);
552      }
553      else if (particleName == "deuteron") {
554         pmanager->AddDiscreteProcess(theElasticProcess);
555         G4DeuteronInelasticProcess* theInelasticProcess = 
556                            new G4DeuteronInelasticProcess("inelastic");
557         G4LEDeuteronInelastic* theInelasticModel = 
558                                 new G4LEDeuteronInelastic;
559         theInelasticProcess->RegisterMe(theInelasticModel);
560         theInelasticProcess->RegisterMe(theTheoModel);
561         pmanager->AddDiscreteProcess(theInelasticProcess);
562      }
563      else if (particleName == "triton") {
564         pmanager->AddDiscreteProcess(theElasticProcess);
565         G4TritonInelasticProcess* theInelasticProcess = 
566                            new G4TritonInelasticProcess("inelastic");
567         G4LETritonInelastic* theInelasticModel = 
568                                 new G4LETritonInelastic;
569         theInelasticProcess->RegisterMe(theInelasticModel);
570         theInelasticProcess->RegisterMe(theTheoModel);
571         pmanager->AddDiscreteProcess(theInelasticProcess);
572      }
573      else if (particleName == "alpha") {
574         pmanager->AddDiscreteProcess(theElasticProcess);
575         G4AlphaInelasticProcess* theInelasticProcess = 
576                            new G4AlphaInelasticProcess("inelastic");
577         G4LEAlphaInelastic* theInelasticModel = 
578                                 new G4LEAlphaInelastic;
579         theInelasticProcess->RegisterMe(theInelasticModel);
580         theInelasticProcess->RegisterMe(theTheoModel);
581         pmanager->AddDiscreteProcess(theInelasticProcess);
582      }
583      else if (particleName == "omega-") {
584         pmanager->AddDiscreteProcess(theElasticProcess);
585         G4OmegaMinusInelasticProcess* theInelasticProcess = 
586                            new G4OmegaMinusInelasticProcess("inelastic");
587         G4LEOmegaMinusInelastic* theInelasticModel = 
588                                 new G4LEOmegaMinusInelastic;
589         theInelasticProcess->RegisterMe(theInelasticModel);
590         theInelasticProcess->RegisterMe(theTheoModel);
591         pmanager->AddDiscreteProcess(theInelasticProcess);
592      }
593      else if (particleName == "anti_omega-") {
594         pmanager->AddDiscreteProcess(theElasticProcess);
595         G4AntiOmegaMinusInelasticProcess* theInelasticProcess = 
596                            new G4AntiOmegaMinusInelasticProcess("inelastic");
597         G4LEAntiOmegaMinusInelastic* theInelasticModel = 
598                                 new G4LEAntiOmegaMinusInelastic;
599         theInelasticProcess->RegisterMe(theInelasticModel);
600         theInelasticProcess->RegisterMe(theTheoModel);
601         pmanager->AddDiscreteProcess(theInelasticProcess);
602      }
603   }
604}
605
606void TiaraPhysicsList::ConstructLeptHad()
607{;}
608
609#include "G4Decay.hh"
610void TiaraPhysicsList::ConstructGeneral()
611{
612  G4Decay* theDecayProcess = new G4Decay();
613  theParticleIterator->reset();
614  while( (*theParticleIterator)() ){
615    G4ParticleDefinition* particle = theParticleIterator->value();
616    G4ProcessManager* pmanager = particle->GetProcessManager();
617    if (theDecayProcess->IsApplicable(*particle)) { 
618      pmanager ->AddProcess(theDecayProcess);
619      pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
620      pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
621    }
622  }
623}
624
625void TiaraPhysicsList::SetCuts()
626{
627  if (verboseLevel >0)
628  {
629    G4cout << "TiaraPhysicsList::SetCuts:";
630    G4cout << "CutLength : " << defaultCutValue/mm << " (mm)" << G4endl;
631  } 
632  //   "G4VUserPhysicsList::SetCutsWithDefault" method sets
633  //   the default cut value for all particle types
634  SetCutsWithDefault();   
635}
Note: See TracBrowser for help on using the repository browser.