]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4PhysicsList.cxx
a9e0dd712c9197f2240d86632ace19c20fe49a82
[u/mrichter/AliRoot.git] / TGeant4 / TG4PhysicsList.cxx
1 // $Id$
2 // Category: physics
3 //
4 // According to:
5 // ExN04PhysicsList.cc,v 1.7 1999/12/15 14:49:26 gunter
6 // GEANT4 tag Name: geant4-01-01
7
8 #include "TG4PhysicsList.h"
9 #include "TG4PhysicsListMessenger.h"
10 #include "TG4PhysicsManager.h"
11 #include "TG4FlagVector.h"
12 #include "TG4SpecialCuts.h"
13 #include "TG4SpecialFlags.h"
14
15 #include <G4ParticleDefinition.hh>
16 #include <G4ProcessManager.hh>
17 #include <G4ProcessVector.hh>
18 #include <G4BosonConstructor.hh>
19 #include <G4LeptonConstructor.hh>
20 #include <G4MesonConstructor.hh>
21 #include <G4BaryonConstructor.hh>
22 #include <G4IonConstructor.hh>
23 #include <G4ShortLivedConstructor.hh>
24 #include <G4ProcessTable.hh>
25
26
27 TG4PhysicsList::TG4PhysicsList()
28   : G4VUserPhysicsList(),
29     fSetOptical(false),
30     fSetHadron(false),
31     fSetSpecialCuts(false),
32     fSetSpecialFlags(false)
33 {
34   // default cut value  (1.0mm) 
35   defaultCutValue = 1.0*mm;
36
37   // messenger
38   fMessenger = new TG4PhysicsListMessenger(this);
39
40   SetVerboseLevel(1);
41 }
42
43 TG4PhysicsList::TG4PhysicsList(const TG4PhysicsList& right)
44   : G4VUserPhysicsList(right)
45 {
46   // messenger
47   fMessenger = new TG4PhysicsListMessenger(this);
48   
49   fSetOptical = right.fSetOptical;
50   fSetHadron = right.fSetHadron;
51   fSetSpecialCuts = right.fSetSpecialCuts;
52   fSetSpecialFlags = right.fSetSpecialFlags;
53 }
54
55 TG4PhysicsList::~TG4PhysicsList() {
56 //
57 }
58
59 // operators
60
61 TG4PhysicsList& TG4PhysicsList::operator=(const TG4PhysicsList& right)
62 {
63   // check assignement to self
64   if (this == &right) return *this;
65
66   // base class assignement
67   G4VUserPhysicsList::operator=(right);
68
69   fSetOptical = right.fSetOptical;
70   fSetHadron = right.fSetHadron;
71   fSetSpecialCuts = right.fSetSpecialCuts;
72   fSetSpecialFlags = right.fSetSpecialFlags;
73   
74   return *this;
75 }  
76
77 // public methods
78
79 void TG4PhysicsList::ConstructParticle()
80 {
81 // In this method, static member functions should be called
82 // for all particles which you want to use.
83 // This ensures that objects of these particle types will be
84 // created in the program. 
85 // ---
86
87   // lock physics manager
88   TG4PhysicsManager* physicsManager = TG4PhysicsManager::Instance();
89   physicsManager->Lock();  
90   physicsManager->SetPhysicsList(this);
91  
92   // create all particles
93   ConstructAllBosons();
94   ConstructAllLeptons();
95   ConstructAllMesons();
96   ConstructAllBaryons();
97   ConstructAllIons();
98   ConstructAllShortLiveds();
99 }
100
101 void TG4PhysicsList::ConstructProcess()
102 {
103 // Constructs all processes.
104 // ---
105
106   AddTransportation();
107
108   ConstructEM();
109   if (fSetHadron)  ConstructHad();
110   if (fSetOptical) ConstructOp();
111   if (fSetSpecialCuts)  ConstructSpecialCuts();
112   if (fSetSpecialFlags) ConstructSpecialFlags();  
113   ConstructGeneral();
114   if (verboseLevel>1) PrintAllProcesses();
115   // InActivateEM();
116 }
117
118 void TG4PhysicsList::SetProcessActivation()
119 {
120 // (In)Activates built processes according
121 // to the setup in TG4PhysicsManager::fFlagVector.
122 // ---
123
124   G4cout << "TG4PhysicsList::SetProcessActivation() start" << G4endl;
125   
126   TG4PhysicsManager* physicsManager = TG4PhysicsManager::Instance();
127   TG4FlagVector* flagVector = physicsManager->GetFlagVector();
128
129   // uncomment following lines to print
130   // the flagVector values
131   //for (G4int i=0; i<kNoG3Flags; i++)
132   //{ cout << i << " flag: " << (*flagVector)[i] << G4endl; }
133
134   if (flagVector) {
135     theParticleIterator->reset();
136     while ((*theParticleIterator)())
137     {
138       G4ParticleDefinition* particle = theParticleIterator->value();
139       G4ProcessManager* processManager = particle->GetProcessManager(); 
140       G4ProcessVector* processVector = processManager->GetProcessList();
141     
142       // set processes flags
143       for (G4int i=0; i<processManager->GetProcessListLength(); i++) {
144         G4int flag = flagVector->GetFlag((*processVector)[i]);
145         if ((flag == kInActivate) && 
146             (processManager->GetProcessActivation(i))) {
147           if (verboseLevel>1) {
148              G4cout << "Set process inactivation for " 
149                     << (*processVector)[i]->GetProcessName() << G4endl;
150           }
151           processManager->SetProcessActivation(i,false);
152         }  
153         else if (((flag == kActivate) || (flag == kActivate2)) &&
154                  (!processManager->GetProcessActivation(i))) {
155           if (verboseLevel>1) {
156              G4cout << "Set process activation for " 
157                     << (*processVector)[i]->GetProcessName() << G4endl;
158           }
159           processManager->SetProcessActivation(i,true);
160         } 
161       }
162     }
163   }
164   else {
165     G4String text = "TG4PhysicsList::SetProcessActivation: \n";
166     text = text + "    Vector of processes flags is not set.";
167     TG4Globals::Warning(text);
168   }    
169   G4cout << "TG4PhysicsList::SetProcessActivation() end" << G4endl;
170 }
171
172 void TG4PhysicsList::PrintAllProcesses() const
173 {
174 // Prints all processes.
175 // ---
176
177   G4cout << "TG4PhysicsList processes: " << G4endl;
178   G4cout << "========================= " << G4endl;
179  
180   G4ProcessTable* processTable = G4ProcessTable::GetProcessTable();
181   G4ProcessTable::G4ProcNameVector* processNameList 
182     = processTable->GetNameList();
183
184   for (G4int i=0; i <processNameList->size(); i++){
185     G4cout << "   " << (*processNameList)[i] << G4endl;
186   }  
187 }
188
189 // protected methods
190
191 #include <G4ComptonScattering.hh>
192 #include <G4GammaConversion.hh>
193 #include <G4PhotoElectricEffect.hh>
194
195 #include <G4MultipleScattering.hh>
196
197 #include <G4eIonisation.hh>
198 #include <G4eBremsstrahlung.hh>
199 #include <G4eplusAnnihilation.hh>
200
201 #include <G4MuIonisation.hh>
202 #include <G4MuBremsstrahlung.hh>
203 #include <G4MuPairProduction.hh>
204
205 #include <G4hIonisation.hh>
206 #include <G4ionIonisation.hh>
207
208 void TG4PhysicsList::ConstructEM()
209 {
210 // Constructs electromagnetic processes.
211 // ---
212
213   theParticleIterator->reset();
214   while( (*theParticleIterator)() ){
215     G4ParticleDefinition* particle = theParticleIterator->value();
216     G4ProcessManager* pmanager = particle->GetProcessManager();
217     G4String particleName = particle->GetParticleName();
218      
219     if (particleName == "gamma") {
220     // gamma
221       // Construct processes for gamma
222       pmanager->AddDiscreteProcess(new G4GammaConversion());
223       pmanager->AddDiscreteProcess(new G4ComptonScattering());      
224       pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
225
226     } else if (particleName == "e-") {
227     //electron
228       // Construct processes for electron
229       G4VProcess* theeminusMultipleScattering = new G4MultipleScattering();
230       G4VProcess* theeminusIonisation = new G4eIonisation();
231       G4VProcess* theeminusBremsstrahlung = new G4eBremsstrahlung();
232       // add processes
233       pmanager->AddProcess(theeminusMultipleScattering);
234       pmanager->AddProcess(theeminusIonisation);
235       pmanager->AddProcess(theeminusBremsstrahlung);      
236       // set ordering for AlongStepDoIt
237       pmanager->SetProcessOrdering(theeminusMultipleScattering, idxAlongStep,  1);
238       pmanager->SetProcessOrdering(theeminusIonisation, idxAlongStep,  2);
239       // set ordering for PostStepDoIt
240       pmanager->SetProcessOrdering(theeminusMultipleScattering, idxPostStep, 1);
241       pmanager->SetProcessOrdering(theeminusIonisation, idxPostStep, 2);
242       pmanager->SetProcessOrdering(theeminusBremsstrahlung, idxPostStep, 3);
243
244     } else if (particleName == "e+") {
245     //positron
246       // Construct processes for positron
247       G4VProcess* theeplusMultipleScattering = new G4MultipleScattering();
248       G4VProcess* theeplusIonisation = new G4eIonisation();
249       G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
250       G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
251       // add processes
252       pmanager->AddProcess(theeplusMultipleScattering);
253       pmanager->AddProcess(theeplusIonisation);
254       pmanager->AddProcess(theeplusBremsstrahlung);
255       pmanager->AddProcess(theeplusAnnihilation);
256       // set ordering for AtRestDoIt
257       pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
258       // set ordering for AlongStepDoIt
259       pmanager->SetProcessOrdering(theeplusMultipleScattering, idxAlongStep,  1);
260       pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep,  2);
261       // set ordering for PostStepDoIt
262       pmanager->SetProcessOrdering(theeplusMultipleScattering, idxPostStep, 1);
263       pmanager->SetProcessOrdering(theeplusIonisation, idxPostStep, 2);
264       pmanager->SetProcessOrdering(theeplusBremsstrahlung, idxPostStep, 3);
265       pmanager->SetProcessOrdering(theeplusAnnihilation, idxPostStep, 4);
266   
267     } else if( particleName == "mu+" || 
268                particleName == "mu-"    ) {
269     //muon  
270      // Construct processes for muon+
271      G4VProcess* aMultipleScattering = new G4MultipleScattering();
272      G4VProcess* aBremsstrahlung = new G4MuBremsstrahlung();
273      G4VProcess* aPairProduction = new G4MuPairProduction();
274      G4VProcess* anIonisation = new G4MuIonisation();
275       // add processes
276      pmanager->AddProcess(anIonisation);
277      pmanager->AddProcess(aMultipleScattering);
278      pmanager->AddProcess(aBremsstrahlung);
279      pmanager->AddProcess(aPairProduction);
280      // set ordering for AlongStepDoIt
281      pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep,  1);
282      pmanager->SetProcessOrdering(anIonisation, idxAlongStep,  2);
283      // set ordering for PostStepDoIt
284      pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep, 1);
285      pmanager->SetProcessOrdering(anIonisation, idxPostStep, 2);
286      pmanager->SetProcessOrdering(aBremsstrahlung, idxPostStep, 3);
287      pmanager->SetProcessOrdering(aPairProduction, idxPostStep, 4);
288      
289     } else if( particleName == "GenericIon" ) {
290      G4VProcess* aionIonization = new G4ionIonisation;
291      G4VProcess* aMultipleScattering = new G4MultipleScattering();
292      pmanager->AddProcess(aionIonization);
293      pmanager->AddProcess(aMultipleScattering);
294      // set ordering for AlongStepDoIt
295      pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep,  1);
296      pmanager->SetProcessOrdering(aionIonization, idxAlongStep,  2);
297      // set ordering for PostStepDoIt
298      pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep, 1);
299      pmanager->SetProcessOrdering(aionIonization, idxPostStep, 2);
300
301    } else if ((!particle->IsShortLived()) &&
302               (particle->GetPDGCharge() != 0.0) && 
303               (particle->GetParticleName() != "chargedgeantino")) {
304      // all others charged particles except geantino
305      G4VProcess* aMultipleScattering = new G4MultipleScattering();
306      G4VProcess* anIonisation = new G4hIonisation();
307      // add processes
308      pmanager->AddProcess(anIonisation);
309      pmanager->AddProcess(aMultipleScattering);
310      // set ordering for AlongStepDoIt
311      pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep,  1);
312      pmanager->SetProcessOrdering(anIonisation, idxAlongStep,  2);
313      // set ordering for PostStepDoIt
314      pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep, 1);
315      pmanager->SetProcessOrdering(anIonisation, idxPostStep, 2);
316     }
317   }
318 }
319
320
321 // Hadron Processes
322
323 #include <G4HadronElasticProcess.hh>
324
325 #include <G4PionPlusInelasticProcess.hh>
326 #include <G4PionMinusInelasticProcess.hh>
327 #include <G4KaonPlusInelasticProcess.hh>
328 #include <G4KaonZeroSInelasticProcess.hh>
329 #include <G4KaonZeroLInelasticProcess.hh>
330 #include <G4KaonMinusInelasticProcess.hh>
331 #include <G4ProtonInelasticProcess.hh>
332 #include <G4AntiProtonInelasticProcess.hh>
333 #include <G4NeutronInelasticProcess.hh>
334 #include <G4AntiNeutronInelasticProcess.hh>
335 #include <G4LambdaInelasticProcess.hh>
336 #include <G4AntiLambdaInelasticProcess.hh>
337 #include <G4SigmaPlusInelasticProcess.hh>
338 #include <G4SigmaMinusInelasticProcess.hh>
339 #include <G4AntiSigmaPlusInelasticProcess.hh>
340 #include <G4AntiSigmaMinusInelasticProcess.hh>
341 #include <G4XiZeroInelasticProcess.hh>
342 #include <G4XiMinusInelasticProcess.hh>
343 #include <G4AntiXiZeroInelasticProcess.hh>
344 #include <G4AntiXiMinusInelasticProcess.hh>
345 #include <G4DeuteronInelasticProcess.hh>
346 #include <G4TritonInelasticProcess.hh>
347 #include <G4AlphaInelasticProcess.hh>
348 #include <G4OmegaMinusInelasticProcess.hh>
349 #include <G4AntiOmegaMinusInelasticProcess.hh>
350
351 // Low-energy Models
352
353 #include <G4LElastic.hh>
354
355 #include <G4LEPionPlusInelastic.hh>
356 #include <G4LEPionMinusInelastic.hh>
357 #include <G4LEKaonPlusInelastic.hh>
358 #include <G4LEKaonZeroSInelastic.hh>
359 #include <G4LEKaonZeroLInelastic.hh>
360 #include <G4LEKaonMinusInelastic.hh>
361 #include <G4LEProtonInelastic.hh>
362 #include <G4LEAntiProtonInelastic.hh>
363 #include <G4LENeutronInelastic.hh>
364 #include <G4LEAntiNeutronInelastic.hh>
365 #include <G4LELambdaInelastic.hh>
366 #include <G4LEAntiLambdaInelastic.hh>
367 #include <G4LESigmaPlusInelastic.hh>
368 #include <G4LESigmaMinusInelastic.hh>
369 #include <G4LEAntiSigmaPlusInelastic.hh>
370 #include <G4LEAntiSigmaMinusInelastic.hh>
371 #include <G4LEXiZeroInelastic.hh>
372 #include <G4LEXiMinusInelastic.hh>
373 #include <G4LEAntiXiZeroInelastic.hh>
374 #include <G4LEAntiXiMinusInelastic.hh>
375 #include <G4LEDeuteronInelastic.hh>
376 #include <G4LETritonInelastic.hh>
377 #include <G4LEAlphaInelastic.hh>
378 #include <G4LEOmegaMinusInelastic.hh>
379 #include <G4LEAntiOmegaMinusInelastic.hh>
380
381 // High-energy Models
382
383 #include <G4HEPionPlusInelastic.hh>
384 #include <G4HEPionMinusInelastic.hh>
385 #include <G4HEKaonPlusInelastic.hh>
386 #include <G4HEKaonZeroInelastic.hh>
387 #include <G4HEKaonZeroInelastic.hh>
388 #include <G4HEKaonMinusInelastic.hh>
389 #include <G4HEProtonInelastic.hh>
390 #include <G4HEAntiProtonInelastic.hh>
391 #include <G4HENeutronInelastic.hh>
392 #include <G4HEAntiNeutronInelastic.hh>
393 #include <G4HELambdaInelastic.hh>
394 #include <G4HEAntiLambdaInelastic.hh>
395 #include <G4HESigmaPlusInelastic.hh>
396 #include <G4HESigmaMinusInelastic.hh>
397 #include <G4HEAntiSigmaPlusInelastic.hh>
398 #include <G4HEAntiSigmaMinusInelastic.hh>
399 #include <G4HEXiZeroInelastic.hh>
400 #include <G4HEXiMinusInelastic.hh>
401 #include <G4HEAntiXiZeroInelastic.hh>
402 #include <G4HEAntiXiMinusInelastic.hh>
403 #include <G4HEOmegaMinusInelastic.hh>
404 #include <G4HEAntiOmegaMinusInelastic.hh>
405
406 // Stopping processes
407
408 #ifdef TRIUMF_STOP_PIMINUS
409 #include <G4PionMinusAbsorptionAtRest.hh>
410 #else
411 #include <G4PiMinusAbsorptionAtRest.hh>
412 #endif
413 #ifdef TRIUMF_STOP_KMINUS
414 #include <G4KaonMinusAbsorption.hh>
415 #else
416 #include <G4KaonMinusAbsorptionAtRest.hh>
417 #endif
418
419 void TG4PhysicsList::ConstructHad()
420 {
421 //
422 // ConstructHad()
423 //
424 // Makes discrete physics processes for the hadrons, at present limited
425 // to those particles with GHEISHA interactions (INTRC > 0).
426 // The processes are: Elastic scattering and Inelastic scattering.
427 //
428 // F.W.Jones  09-JUL-1998
429 // ---
430
431    G4cout << "### TG4PhysicsList::ConstructHad()" << G4endl;
432
433    G4HadronElasticProcess* theElasticProcess = 
434                                     new G4HadronElasticProcess;
435    G4LElastic* theElasticModel = new G4LElastic;
436    theElasticProcess->RegisterMe(theElasticModel);
437
438    theParticleIterator->reset();
439    while ((*theParticleIterator)()) {
440       G4ParticleDefinition* particle = theParticleIterator->value();
441       G4ProcessManager* pmanager = particle->GetProcessManager();
442       G4String particleName = particle->GetParticleName();
443      
444       if (particleName == "pi+") {
445          pmanager->AddDiscreteProcess(theElasticProcess);
446          G4PionPlusInelasticProcess* theInelasticProcess = 
447                                 new G4PionPlusInelasticProcess("inelastic");
448          G4LEPionPlusInelastic* theLEInelasticModel = 
449                                 new G4LEPionPlusInelastic;
450          theInelasticProcess->RegisterMe(theLEInelasticModel);
451          G4HEPionPlusInelastic* theHEInelasticModel = 
452                                 new G4HEPionPlusInelastic;
453          theInelasticProcess->RegisterMe(theHEInelasticModel);
454          pmanager->AddDiscreteProcess(theInelasticProcess);
455       }
456       else if (particleName == "pi-") {
457          pmanager->AddDiscreteProcess(theElasticProcess);
458          G4PionMinusInelasticProcess* theInelasticProcess = 
459                                 new G4PionMinusInelasticProcess("inelastic");
460          G4LEPionMinusInelastic* theLEInelasticModel = 
461                                 new G4LEPionMinusInelastic;
462          theInelasticProcess->RegisterMe(theLEInelasticModel);
463          G4HEPionMinusInelastic* theHEInelasticModel = 
464                                 new G4HEPionMinusInelastic;
465          theInelasticProcess->RegisterMe(theHEInelasticModel);
466          pmanager->AddDiscreteProcess(theInelasticProcess);
467 #ifdef TRIUMF_STOP_PIMINUS
468          pmanager->AddRestProcess(new G4PionMinusAbsorptionAtRest, ordDefault);
469 #else
470          G4String prcNam;
471          pmanager->AddRestProcess(
472            new G4PiMinusAbsorptionAtRest(
473                 prcNam="PiMinusAbsorptionAtRest"), ordDefault);
474 #endif
475       }
476       else if (particleName == "kaon+") {
477          pmanager->AddDiscreteProcess(theElasticProcess);
478          G4KaonPlusInelasticProcess* theInelasticProcess = 
479                                   new G4KaonPlusInelasticProcess("inelastic");
480          G4LEKaonPlusInelastic* theLEInelasticModel = 
481                                   new G4LEKaonPlusInelastic;
482          theInelasticProcess->RegisterMe(theLEInelasticModel);
483          G4HEKaonPlusInelastic* theHEInelasticModel = 
484                                   new G4HEKaonPlusInelastic;
485          theInelasticProcess->RegisterMe(theHEInelasticModel);
486          pmanager->AddDiscreteProcess(theInelasticProcess);
487       }
488       else if (particleName == "kaon0S") {
489          pmanager->AddDiscreteProcess(theElasticProcess);
490          G4KaonZeroSInelasticProcess* theInelasticProcess = 
491                              new G4KaonZeroSInelasticProcess("inelastic");
492          G4LEKaonZeroSInelastic* theLEInelasticModel = 
493                              new G4LEKaonZeroSInelastic;
494          theInelasticProcess->RegisterMe(theLEInelasticModel);
495          G4HEKaonZeroInelastic* theHEInelasticModel = 
496                              new G4HEKaonZeroInelastic;
497          theInelasticProcess->RegisterMe(theHEInelasticModel);
498          pmanager->AddDiscreteProcess(theInelasticProcess);
499       }
500       else if (particleName == "kaon0L") {
501          pmanager->AddDiscreteProcess(theElasticProcess);
502          G4KaonZeroLInelasticProcess* theInelasticProcess = 
503                              new G4KaonZeroLInelasticProcess("inelastic");
504          G4LEKaonZeroLInelastic* theLEInelasticModel = 
505                              new G4LEKaonZeroLInelastic;
506          theInelasticProcess->RegisterMe(theLEInelasticModel);
507          G4HEKaonZeroInelastic* theHEInelasticModel = 
508                              new G4HEKaonZeroInelastic;
509          theInelasticProcess->RegisterMe(theHEInelasticModel);
510          pmanager->AddDiscreteProcess(theInelasticProcess);
511       }
512       else if (particleName == "kaon-") {
513          pmanager->AddDiscreteProcess(theElasticProcess);
514          G4KaonMinusInelasticProcess* theInelasticProcess = 
515                                  new G4KaonMinusInelasticProcess("inelastic");
516          G4LEKaonMinusInelastic* theLEInelasticModel = 
517                                  new G4LEKaonMinusInelastic;
518          theInelasticProcess->RegisterMe(theLEInelasticModel);
519          G4HEKaonMinusInelastic* theHEInelasticModel = 
520                                  new G4HEKaonMinusInelastic;
521          theInelasticProcess->RegisterMe(theHEInelasticModel);
522          pmanager->AddDiscreteProcess(theInelasticProcess);
523 #ifdef TRIUMF_STOP_KMINUS
524          pmanager->AddRestProcess(new G4KaonMinusAbsorption, ordDefault);
525 #else
526          pmanager->AddRestProcess(new G4KaonMinusAbsorptionAtRest, ordDefault);
527 #endif
528       }
529       else if (particleName == "proton") {
530          pmanager->AddDiscreteProcess(theElasticProcess);
531          G4ProtonInelasticProcess* theInelasticProcess = 
532                                     new G4ProtonInelasticProcess("inelastic");
533          G4LEProtonInelastic* theLEInelasticModel = new G4LEProtonInelastic;
534          theInelasticProcess->RegisterMe(theLEInelasticModel);
535          G4HEProtonInelastic* theHEInelasticModel = new G4HEProtonInelastic;
536          theInelasticProcess->RegisterMe(theHEInelasticModel);
537          pmanager->AddDiscreteProcess(theInelasticProcess);
538       }
539       else if (particleName == "anti_proton") {
540          pmanager->AddDiscreteProcess(theElasticProcess);
541          G4AntiProtonInelasticProcess* theInelasticProcess = 
542                                 new G4AntiProtonInelasticProcess("inelastic");
543          G4LEAntiProtonInelastic* theLEInelasticModel = 
544                                 new G4LEAntiProtonInelastic;
545          theInelasticProcess->RegisterMe(theLEInelasticModel);
546          G4HEAntiProtonInelastic* theHEInelasticModel = 
547                                 new G4HEAntiProtonInelastic;
548          theInelasticProcess->RegisterMe(theHEInelasticModel);
549          pmanager->AddDiscreteProcess(theInelasticProcess);
550       }
551       else if (particleName == "neutron") {
552          pmanager->AddDiscreteProcess(theElasticProcess);
553          G4NeutronInelasticProcess* theInelasticProcess = 
554                                     new G4NeutronInelasticProcess("inelastic");
555          G4LENeutronInelastic* theLEInelasticModel = 
556                                     new G4LENeutronInelastic;
557          theInelasticProcess->RegisterMe(theLEInelasticModel);
558          G4HENeutronInelastic* theHEInelasticModel = 
559                                     new G4HENeutronInelastic;
560          theInelasticProcess->RegisterMe(theHEInelasticModel);
561          pmanager->AddDiscreteProcess(theInelasticProcess);
562       }  
563       else if (particleName == "anti_neutron") {
564          pmanager->AddDiscreteProcess(theElasticProcess);
565          G4AntiNeutronInelasticProcess* theInelasticProcess = 
566                                new G4AntiNeutronInelasticProcess("inelastic");
567          G4LEAntiNeutronInelastic* theLEInelasticModel = 
568                                new G4LEAntiNeutronInelastic;
569          theInelasticProcess->RegisterMe(theLEInelasticModel);
570          G4HEAntiNeutronInelastic* theHEInelasticModel = 
571                                new G4HEAntiNeutronInelastic;
572          theInelasticProcess->RegisterMe(theHEInelasticModel);
573          pmanager->AddDiscreteProcess(theInelasticProcess);
574       }
575       else if (particleName == "lambda") {
576          pmanager->AddDiscreteProcess(theElasticProcess);
577          G4LambdaInelasticProcess* theInelasticProcess = 
578                                     new G4LambdaInelasticProcess("inelastic");
579          G4LELambdaInelastic* theLEInelasticModel = new G4LELambdaInelastic;
580          theInelasticProcess->RegisterMe(theLEInelasticModel);
581          G4HELambdaInelastic* theHEInelasticModel = new G4HELambdaInelastic;
582          theInelasticProcess->RegisterMe(theHEInelasticModel);
583          pmanager->AddDiscreteProcess(theInelasticProcess);
584       }
585       else if (particleName == "anti_lambda") {
586          pmanager->AddDiscreteProcess(theElasticProcess);
587          G4AntiLambdaInelasticProcess* theInelasticProcess = 
588                                 new G4AntiLambdaInelasticProcess("inelastic");
589          G4LEAntiLambdaInelastic* theLEInelasticModel = 
590                                 new G4LEAntiLambdaInelastic;
591          theInelasticProcess->RegisterMe(theLEInelasticModel);
592          G4HEAntiLambdaInelastic* theHEInelasticModel = 
593                                 new G4HEAntiLambdaInelastic;
594          theInelasticProcess->RegisterMe(theHEInelasticModel);
595          pmanager->AddDiscreteProcess(theInelasticProcess);
596       }
597       else if (particleName == "sigma+") {
598          pmanager->AddDiscreteProcess(theElasticProcess);
599          G4SigmaPlusInelasticProcess* theInelasticProcess = 
600                                  new G4SigmaPlusInelasticProcess("inelastic");
601          G4LESigmaPlusInelastic* theLEInelasticModel = 
602                                  new G4LESigmaPlusInelastic;
603          theInelasticProcess->RegisterMe(theLEInelasticModel);
604          G4HESigmaPlusInelastic* theHEInelasticModel = 
605                                  new G4HESigmaPlusInelastic;
606          theInelasticProcess->RegisterMe(theHEInelasticModel);
607          pmanager->AddDiscreteProcess(theInelasticProcess);
608       }
609       else if (particleName == "sigma-") {
610          pmanager->AddDiscreteProcess(theElasticProcess);
611          G4SigmaMinusInelasticProcess* theInelasticProcess = 
612                                  new G4SigmaMinusInelasticProcess("inelastic");
613          G4LESigmaMinusInelastic* theLEInelasticModel = 
614                                  new G4LESigmaMinusInelastic;
615          theInelasticProcess->RegisterMe(theLEInelasticModel);
616          G4HESigmaMinusInelastic* theHEInelasticModel = 
617                                  new G4HESigmaMinusInelastic;
618          theInelasticProcess->RegisterMe(theHEInelasticModel);
619          pmanager->AddDiscreteProcess(theInelasticProcess);
620       }
621       else if (particleName == "anti_sigma+") {
622          pmanager->AddDiscreteProcess(theElasticProcess);
623          G4AntiSigmaPlusInelasticProcess* theInelasticProcess = 
624                              new G4AntiSigmaPlusInelasticProcess("inelastic");
625          G4LEAntiSigmaPlusInelastic* theLEInelasticModel = 
626                                  new G4LEAntiSigmaPlusInelastic;
627          theInelasticProcess->RegisterMe(theLEInelasticModel);
628          G4HEAntiSigmaPlusInelastic* theHEInelasticModel = 
629                                  new G4HEAntiSigmaPlusInelastic;
630          theInelasticProcess->RegisterMe(theHEInelasticModel);
631          pmanager->AddDiscreteProcess(theInelasticProcess);
632       }
633       else if (particleName == "anti_sigma-") {
634          pmanager->AddDiscreteProcess(theElasticProcess);
635          G4AntiSigmaMinusInelasticProcess* theInelasticProcess = 
636                             new G4AntiSigmaMinusInelasticProcess("inelastic");
637          G4LEAntiSigmaMinusInelastic* theLEInelasticModel = 
638                                  new G4LEAntiSigmaMinusInelastic;
639          theInelasticProcess->RegisterMe(theLEInelasticModel);
640          G4HEAntiSigmaMinusInelastic* theHEInelasticModel = 
641                                  new G4HEAntiSigmaMinusInelastic;
642          theInelasticProcess->RegisterMe(theHEInelasticModel);
643          pmanager->AddDiscreteProcess(theInelasticProcess);
644       }
645       else if (particleName == "xi0") {
646          pmanager->AddDiscreteProcess(theElasticProcess);
647          G4XiZeroInelasticProcess* theInelasticProcess = 
648                             new G4XiZeroInelasticProcess("inelastic");
649          G4LEXiZeroInelastic* theLEInelasticModel = 
650                                  new G4LEXiZeroInelastic;
651          theInelasticProcess->RegisterMe(theLEInelasticModel);
652          G4HEXiZeroInelastic* theHEInelasticModel = 
653                                  new G4HEXiZeroInelastic;
654          theInelasticProcess->RegisterMe(theHEInelasticModel);
655          pmanager->AddDiscreteProcess(theInelasticProcess);
656       }
657       else if (particleName == "xi-") {
658          pmanager->AddDiscreteProcess(theElasticProcess);
659          G4XiMinusInelasticProcess* theInelasticProcess = 
660                             new G4XiMinusInelasticProcess("inelastic");
661          G4LEXiMinusInelastic* theLEInelasticModel = 
662                                  new G4LEXiMinusInelastic;
663          theInelasticProcess->RegisterMe(theLEInelasticModel);
664          G4HEXiMinusInelastic* theHEInelasticModel = 
665                                  new G4HEXiMinusInelastic;
666          theInelasticProcess->RegisterMe(theHEInelasticModel);
667          pmanager->AddDiscreteProcess(theInelasticProcess);
668       }
669       else if (particleName == "anti_xi0") {
670          pmanager->AddDiscreteProcess(theElasticProcess);
671          G4AntiXiZeroInelasticProcess* theInelasticProcess = 
672                             new G4AntiXiZeroInelasticProcess("inelastic");
673          G4LEAntiXiZeroInelastic* theLEInelasticModel = 
674                                  new G4LEAntiXiZeroInelastic;
675          theInelasticProcess->RegisterMe(theLEInelasticModel);
676          G4HEAntiXiZeroInelastic* theHEInelasticModel = 
677                                  new G4HEAntiXiZeroInelastic;
678          theInelasticProcess->RegisterMe(theHEInelasticModel);
679          pmanager->AddDiscreteProcess(theInelasticProcess);
680       }
681       else if (particleName == "anti_xi-") {
682          pmanager->AddDiscreteProcess(theElasticProcess);
683          G4AntiXiMinusInelasticProcess* theInelasticProcess = 
684                             new G4AntiXiMinusInelasticProcess("inelastic");
685          G4LEAntiXiMinusInelastic* theLEInelasticModel = 
686                                  new G4LEAntiXiMinusInelastic;
687          theInelasticProcess->RegisterMe(theLEInelasticModel);
688          G4HEAntiXiMinusInelastic* theHEInelasticModel = 
689                                  new G4HEAntiXiMinusInelastic;
690          theInelasticProcess->RegisterMe(theHEInelasticModel);
691          pmanager->AddDiscreteProcess(theInelasticProcess);
692       }
693       else if (particleName == "deuteron") {
694          pmanager->AddDiscreteProcess(theElasticProcess);
695          G4DeuteronInelasticProcess* theInelasticProcess = 
696                             new G4DeuteronInelasticProcess("inelastic");
697          G4LEDeuteronInelastic* theLEInelasticModel = 
698                                  new G4LEDeuteronInelastic;
699          theInelasticProcess->RegisterMe(theLEInelasticModel);
700          pmanager->AddDiscreteProcess(theInelasticProcess);
701       }
702       else if (particleName == "triton") {
703          pmanager->AddDiscreteProcess(theElasticProcess);
704          G4TritonInelasticProcess* theInelasticProcess = 
705                             new G4TritonInelasticProcess("inelastic");
706          G4LETritonInelastic* theLEInelasticModel = 
707                                  new G4LETritonInelastic;
708          theInelasticProcess->RegisterMe(theLEInelasticModel);
709          pmanager->AddDiscreteProcess(theInelasticProcess);
710       }
711       else if (particleName == "alpha") {
712          pmanager->AddDiscreteProcess(theElasticProcess);
713          G4AlphaInelasticProcess* theInelasticProcess = 
714                             new G4AlphaInelasticProcess("inelastic");
715          G4LEAlphaInelastic* theLEInelasticModel = 
716                                  new G4LEAlphaInelastic;
717          theInelasticProcess->RegisterMe(theLEInelasticModel);
718          pmanager->AddDiscreteProcess(theInelasticProcess);
719       }
720       else if (particleName == "omega-") {
721          pmanager->AddDiscreteProcess(theElasticProcess);
722          G4OmegaMinusInelasticProcess* theInelasticProcess = 
723                             new G4OmegaMinusInelasticProcess("inelastic");
724          G4LEOmegaMinusInelastic* theLEInelasticModel = 
725                                  new G4LEOmegaMinusInelastic;
726          theInelasticProcess->RegisterMe(theLEInelasticModel);
727          G4HEOmegaMinusInelastic* theHEInelasticModel = 
728                                  new G4HEOmegaMinusInelastic;
729          theInelasticProcess->RegisterMe(theHEInelasticModel);
730          pmanager->AddDiscreteProcess(theInelasticProcess);
731       }
732       else if (particleName == "anti_omega-") {
733          pmanager->AddDiscreteProcess(theElasticProcess);
734          G4AntiOmegaMinusInelasticProcess* theInelasticProcess = 
735                             new G4AntiOmegaMinusInelasticProcess("inelastic");
736          G4LEAntiOmegaMinusInelastic* theLEInelasticModel = 
737                                  new G4LEAntiOmegaMinusInelastic;
738          theInelasticProcess->RegisterMe(theLEInelasticModel);
739          G4HEAntiOmegaMinusInelastic* theHEInelasticModel = 
740                                  new G4HEAntiOmegaMinusInelastic;
741          theInelasticProcess->RegisterMe(theHEInelasticModel);
742          pmanager->AddDiscreteProcess(theInelasticProcess);
743       }
744    }
745
746    G4cout << "### TG4PhysicsList::ConstructHad() finished." << G4endl;
747
748 }
749
750
751 #include <G4Cerenkov.hh>
752 #include <G4OpAbsorption.hh>
753 #include <G4OpRayleigh.hh>
754 #include <G4OpBoundaryProcess.hh>
755
756 void TG4PhysicsList::ConstructOp()
757 {
758 // Constructs optical processes.
759 // According to ExN06PhysicsList.cc.
760 // (geant4 1.1)
761 // ---
762
763   G4Cerenkov*     theCerenkovProcess = new G4Cerenkov("Cerenkov");
764   G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
765   G4OpRayleigh*   theRayleighScatteringProcess = new G4OpRayleigh();
766   G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
767
768   theCerenkovProcess->DumpPhysicsTable();
769   //theAbsorptionProcess->DumpPhysicsTable();
770   //theRayleighScatteringProcess->DumpPhysicsTable();
771
772   // add verbose 
773   //theCerenkovProcess->SetVerboseLevel(1);
774   //theAbsorptionProcess->SetVerboseLevel(1);
775   //theRayleighScatteringProcess->SetVerboseLevel(1);
776   //theBoundaryProcess->SetVerboseLevel(1);
777
778   G4int maxNumPhotons = 300;
779
780   theCerenkovProcess->SetTrackSecondariesFirst(true);
781   theCerenkovProcess->SetMaxNumPhotonsPerStep(maxNumPhotons);
782
783   //G4OpticalSurfaceModel themodel = unified;   
784   // model from GEANT3
785   G4OpticalSurfaceModel themodel = glisur;
786   theBoundaryProcess->SetModel(themodel);
787
788   theParticleIterator->reset();
789   while( (*theParticleIterator)() ){
790     G4ParticleDefinition* particle = theParticleIterator->value();
791     G4ProcessManager* processManager = particle->GetProcessManager();
792     G4String particleName = particle->GetParticleName();
793     if (theCerenkovProcess->IsApplicable(*particle)) {
794       processManager->AddContinuousProcess(theCerenkovProcess);
795     }
796     if (particleName == "opticalphoton") {
797       G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl;
798       processManager->AddDiscreteProcess(theAbsorptionProcess);
799       processManager->AddDiscreteProcess(theRayleighScatteringProcess);
800       processManager->AddDiscreteProcess(theBoundaryProcess);
801     }
802   }
803 }
804
805 #include <G4Decay.hh>
806
807 void TG4PhysicsList::ConstructGeneral()
808 {
809 // Constructs general processes.
810 // ---
811
812   // Add Decay Process
813   G4Decay* theDecayProcess = new G4Decay();
814   theParticleIterator->reset();
815   while( (*theParticleIterator)() ){
816     G4ParticleDefinition* particle = theParticleIterator->value();
817     G4ProcessManager* pmanager = particle->GetProcessManager();
818     if (theDecayProcess->IsApplicable(*particle)) { 
819       pmanager ->AddProcess(theDecayProcess);
820       // set ordering for PostStepDoIt and AtRestDoIt
821       pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
822       pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
823     }
824   }
825 }
826
827 void TG4PhysicsList::ConstructSpecialCuts()
828 {
829 // Adds TG4SpecialCuts "process" that activates
830 // the kinetic energy cuts defined in 
831 // the vector of cuts (PhysicsManager::fCutVector) or in TG4Limits.
832 // ---
833
834   TG4PhysicsManager* physicsManager 
835     = TG4PhysicsManager::Instance();
836
837   if (physicsManager->IsSpecialCuts())
838   {
839     TG4CutVector* cutVector
840       = physicsManager->GetCutVector(); 
841     TG4boolVector* isCutVector 
842       = physicsManager->GetIsCutVector(); 
843
844     theParticleIterator->reset();
845     while ((*theParticleIterator)())
846     {
847       G4ParticleDefinition* particle = theParticleIterator->value();
848       TG3ParticleWSP particleWSP 
849         = physicsManager->GetG3ParticleWSP(particle);
850       G4String name;
851       physicsManager->GetG3ParticleWSPName(particleWSP, name);
852       
853       // uncomment this to see all particles "WSP"
854       //G4cout << "Iterating particle: " 
855       //       << particle->GetParticleName() << " " << particleWSP << " "
856       //       << name << G4endl;
857
858       // special process is created in case
859       // cutVector (vector of kinetic energy cuts) is set
860       // or the special cut is set by TG4Limits
861       if ((particleWSP !=kNofParticlesWSP) && 
862           ((*isCutVector)[particleWSP])) {
863         // check if process already exists
864         G4String processName = "specialCutFor" + name;
865         G4VProcess* process = FindProcess(processName);
866         if (!process) {
867           process = new TG4SpecialCuts(particleWSP, cutVector, processName);
868         }  
869         //particle->GetProcessManager()->AddProcess(process, 0, -1, 1);
870         particle->GetProcessManager()->AddDiscreteProcess(process);
871       }
872     }
873
874     if (verboseLevel>0) {
875       G4cout << "TG4PhysicsList::ConstructSpecialCuts: " << G4endl;
876       if (cutVector)
877         G4cout << "   Global kinetic energy cuts are set." << G4endl;
878       G4cout << "   Special cuts process is defined for: " << G4endl 
879              << "   ";
880       for (G4int i=0; i<kAny; i++) {
881         G4String name;
882         physicsManager->GetG3ParticleWSPName(i, name);
883         if ((*isCutVector)[i]) G4cout << name << " ";
884       }  
885       G4cout << G4endl;
886     }  
887   }
888 }
889
890 void TG4PhysicsList::ConstructSpecialFlags()
891 {
892 // Adds TG4SpecialFlags "process" that activates
893 // the control process flags defined in TG4Limits.
894 // ---
895
896   TG4PhysicsManager* physicsManager 
897     = TG4PhysicsManager::Instance();
898
899   if (physicsManager->IsSpecialFlags())
900   {
901     G4cout << "IsSpecialFlags started" << G4endl;
902     TG4boolVector* isFlagVector 
903       = physicsManager->GetIsFlagVector(); 
904
905     theParticleIterator->reset();
906     while ((*theParticleIterator)())
907     {
908       G4ParticleDefinition* particle = theParticleIterator->value();
909       TG3ParticleWSP particleWSP 
910         = physicsManager->GetG3ParticleWSP(particle);
911       //G4String name;
912       //GetG3ParticleWSPName(particleWSP, name);
913
914       // special process is set in case
915       // the special flag is set by TG4Limits    
916       if ((particleWSP !=kNofParticlesWSP) && 
917           ((*isFlagVector)[particleWSP])) {
918         // check if process already exists
919         G4String processName = "specialFlag";
920         G4VProcess* process = FindProcess(processName);
921         if (!process) {
922           process = new TG4SpecialFlags(processName);
923         }  
924         //particle->GetProcessManager()->AddProcess(process, 0, -1, 1);
925         particle->GetProcessManager()->AddDiscreteProcess(process);
926       }
927     }
928
929     if (verboseLevel>0) {
930       G4cout << "TG4PhysicsList::ConstructSpecialFlagss: " << G4endl;
931       G4cout << "   Special flags process is defined for: " << G4endl
932              << "   ";
933       for (G4int i=0; i<kNofParticlesWSP; i++) {
934         G4String name;
935         physicsManager->GetG3ParticleWSPName(i, name);
936         if ((*isFlagVector)[i]) G4cout << name << " ";
937       }  
938       G4cout << G4endl;
939     }  
940   }
941 }
942
943 void TG4PhysicsList::SetCuts()
944 {
945 // Sets the default cut value for all particle types
946 // other then e-/e+. 
947 // The cut value for e-/e+ is high in oredr to supress
948 // tracking of delta electrons.
949 // ---
950
951   // SetCutsWithDefault();   
952          // "G4VUserPhysicsList::SetCutsWithDefault" method sets 
953          // the default cut value for all particle types.
954
955   // default cut value
956   G4double cut  = defaultCutValue;
957   G4double ecut = 10.*m; 
958   //G4double ecut = cut; 
959
960 #ifdef G4VERBOSE    
961   if (verboseLevel >1){
962     G4cout << "G4VUserPhysicsList::SetCutsWithDefault:";
963     G4cout << "CutLength : " << cut/mm << " (mm)" << G4endl;
964   }  
965 #endif
966
967   // set cut values for gamma at first and for e- second and next for e+,
968   // because some processes for e+/e- need cut values for gamma 
969   SetCutValue(cut, "gamma");
970   SetCutValue(ecut, "e-");
971   SetCutValue(ecut, "e+");
972  
973   // set cut values for proton and anti_proton before all other hadrons
974   // because some processes for hadrons need cut values for proton/anti_proton 
975   SetCutValue(cut, "proton");
976   SetCutValue(cut, "anti_proton");
977   
978   SetCutValueForOthers(cut);
979
980   if (verboseLevel>1) {
981     DumpCutValuesTable();
982   }
983 }
984
985 void TG4PhysicsList::ConstructAllBosons()
986 {
987 // Construct all bosons
988 // ---
989
990   G4BosonConstructor pConstructor;
991   pConstructor.ConstructParticle();
992 }
993
994 void TG4PhysicsList::ConstructAllLeptons()
995 {
996 // Construct all leptons
997 // ---
998
999   G4LeptonConstructor pConstructor;
1000   pConstructor.ConstructParticle();
1001 }
1002
1003 void TG4PhysicsList::ConstructAllMesons()
1004 {
1005 // Construct all mesons
1006 // ---
1007
1008   G4MesonConstructor pConstructor;
1009   pConstructor.ConstructParticle();
1010 }
1011
1012 void TG4PhysicsList::ConstructAllBaryons()
1013 {
1014 // Construct all barions
1015 // ---
1016
1017   G4BaryonConstructor pConstructor;
1018   pConstructor.ConstructParticle();
1019 }
1020
1021 void TG4PhysicsList::ConstructAllIons()
1022 {
1023 // Construct light ions
1024 // ---
1025
1026   G4IonConstructor pConstructor;
1027   pConstructor.ConstructParticle();  
1028 }
1029
1030 void TG4PhysicsList::ConstructAllShortLiveds()
1031 {
1032 // Construct  resonaces and quarks
1033 // ---
1034
1035   G4ShortLivedConstructor pConstructor;
1036   pConstructor.ConstructParticle();  
1037 }
1038
1039 // private methods
1040
1041 G4VProcess* TG4PhysicsList::FindProcess(G4String processName) const
1042 {
1043 // Finds G4VProcess with specified name.
1044 // ---
1045
1046   G4ProcessTable* processTable = G4ProcessTable::GetProcessTable();
1047
1048   G4ProcessVector* processVector 
1049     = processTable->FindProcesses(processName);
1050   G4VProcess* firstFoundProcess = 0;
1051   if (processVector->entries()>0) firstFoundProcess= (*processVector)[0];
1052
1053   processVector->clear();
1054   delete processVector;
1055   
1056   return firstFoundProcess;
1057 }
1058
1059 void TG4PhysicsList::InActivateProcess(G4String processName, 
1060                         G4ParticleDefinition* particle)              
1061 {
1062 // Activates the process specified by name for the specified 
1063 // particle.
1064 // Only for tests - to be removed.
1065 // ---
1066
1067   G4ProcessManager* processManager = particle->GetProcessManager();
1068   G4ProcessVector* processVector = processManager->GetProcessList();
1069   for (G4int i=0; i<processVector->entries(); i++) {
1070     if ((*processVector)[i]->GetProcessName() == processName) {
1071       processManager->SetProcessActivation((*processVector)[i], false);
1072       return;
1073     }
1074   }
1075  
1076   G4String text = "TG4PhysicsList::InActivateProcess: ";
1077   text = text + processName + " is not set for ";
1078   text = text + particle->GetParticleName();
1079   TG4Globals::Exception(text);
1080 }
1081
1082 void TG4PhysicsList::InActivateEM()
1083 {
1084 // Inactivates specified electromagnetic processes.
1085 // Only for tests - to be removed.
1086 // !! This method must be called after all Construct methods.
1087 // Uncomment the selected line(s) to inactivate desired processes.
1088 // ---
1089
1090  theParticleIterator->reset();
1091   while ((*theParticleIterator)())
1092   {
1093     G4ParticleDefinition* particle = theParticleIterator->value();
1094     G4String name = particle->GetParticleName();     
1095     
1096     if (name == "gamma") {
1097       // gamma
1098       //InActivateProcess("phot",  particle);      
1099       //InActivateProcess("compt", particle);      
1100       //InActivateProcess("conv", particle);      
1101     } 
1102     else if (name == "e-") {
1103       //electron
1104       InActivateProcess("msc",  particle);      
1105       G4cout << "msc inactivated." << G4endl;
1106       //InActivateProcess("eIoni",  particle);      
1107       //G4cout << "eIoni inactivated." << G4endl;
1108       InActivateProcess("eBrem",  particle);      
1109       G4cout << "eBrem inactivated." << G4endl;
1110     } 
1111     else if (name == "e+") {
1112       //positron
1113       //InActivateProcess("msc",  particle);      
1114       //InActivateProcess("eIoni",  particle);      
1115       //InActivateProcess("eBrem",  particle);      
1116       //InActivateProcess("annihil",  particle);      
1117     } 
1118     else if (name == "mu+" || name == "mu-") {
1119       //muon  
1120       //InActivateProcess("msc",  particle);      
1121       //InActivateProcess("MuIoni",  particle);      
1122       //InActivateProcess("MuBrem",  particle);      
1123       //InActivateProcess("MuPairProd",  particle);      
1124     } 
1125     else if ((!particle->IsShortLived()) &&
1126              (particle->GetPDGCharge() != 0.0) && 
1127              (particle->GetParticleName() != "chargedgeantino")) 
1128     {
1129       // all others charged particles except geantino
1130       //InActivateProcess("msc",  particle);      
1131       //InActivateProcess("hIoni",  particle);      
1132     }
1133   }
1134 }
1135