source: trunk/source/run/src/G4VUserPhysicsList.cc @ 1168

Last change on this file since 1168 was 1058, checked in by garnier, 15 years ago

file release beta

File size: 23.9 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: G4VUserPhysicsList.cc,v 1.64 2008/05/09 13:00:42 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//
31// ------------------------------------------------------------
32//      GEANT 4 class header file
33//
34// ------------------------------------------------------------
35//      History
36//       first version                    09 Jan 1998 by H.Kurashige
37//       Added SetEnergyRange             18 Jun 1998 by H.Kurashige
38//       Change for short lived particles 27 Jun 1998 by H.Kurashige
39//       G4BestUnit on output             12 nov 1998 by M.Maire
40//       Added RemoveProcessManager        9 Feb 1999 by H.Kurashige
41//       Fixed RemoveProcessManager       15 Apr 1999 by H.Kurashige
42//       Removed ConstructAllParticles()  15 Apr 1999 by H.Kurashige
43//       Modified for CUTS per REGION     10 Oct 2002 by H.Kurashige
44//       Check if particle IsShortLived   18 Jun 2003 by V.Ivanchenko
45//       Modify PreparePhysicsList        18 Jan 2006 by H.Kurashige
46// ------------------------------------------------------------
47
48#include "globals.hh"
49#include "G4VUserPhysicsList.hh"
50#include "G4ParticleDefinition.hh"
51#include "G4ProcessManager.hh"
52#include "G4ParticleTable.hh"
53#include "G4ProductionCutsTable.hh"
54#include "G4Material.hh"
55#include "G4UserPhysicsListMessenger.hh"
56#include "G4UImanager.hh"
57#include "G4UnitsTable.hh"
58#include "G4RegionStore.hh"
59#include "G4Region.hh"
60#include "G4ProductionCutsTable.hh"
61#include "G4ProductionCuts.hh"
62#include "G4MaterialCutsCouple.hh"
63
64#include "G4ios.hh"
65#include <iomanip>
66#include <fstream>
67
68////////////////////////////////////////////////////////
69G4VUserPhysicsList::G4VUserPhysicsList()
70                   :verboseLevel(1),
71                    fRetrievePhysicsTable(false),
72                    fStoredInAscii(true),
73                    fIsCheckedForRetrievePhysicsTable(false),
74                    fIsRestoredCutValues(false),
75                    directoryPhysicsTable("."),
76                    fDisplayThreshold(0),
77                    useCoupledTransportation(false)
78{
79  // default cut value  (1.0mm)
80  defaultCutValue = 1.0*mm;
81
82  // pointer to the particle table
83  theParticleTable = G4ParticleTable::GetParticleTable();
84  theParticleIterator = theParticleTable->GetIterator();
85
86  // pointer to the cuts table
87  fCutsTable =  G4ProductionCutsTable::GetProductionCutsTable();
88
89  // set energy range for SetCut calcuration
90  fCutsTable->SetEnergyRange(0.99*keV, 100*TeV);
91
92  // UI Messenger
93  theMessenger = new G4UserPhysicsListMessenger(this);
94}
95
96////////////////////////////////////////////////////////
97G4VUserPhysicsList::~G4VUserPhysicsList()
98{
99  if (theMessenger != 0) {
100    delete theMessenger;
101    theMessenger = 0;
102  }
103  RemoveProcessManager();
104
105  // invoke DeleteAllParticle
106  theParticleTable->DeleteAllParticles();
107
108}
109
110////////////////////////////////////////////////////////
111void G4VUserPhysicsList::SetVerboseLevel(G4int value)
112{
113  verboseLevel = value;
114  // set verboseLevel for G4ProductionCutsTable same as one for G4VUserPhysicsList:
115  fCutsTable->SetVerboseLevel(value);
116
117#ifdef G4VERBOSE
118  if (verboseLevel >1){
119    G4cout << "G4VUserPhysicsList::SetVerboseLevel  :";
120    G4cout << " Verbose level is set to " << verboseLevel << G4endl;
121  }
122#endif
123}
124
125////////////////////////////////////////////////////////
126void G4VUserPhysicsList::AddProcessManager(G4ParticleDefinition* newParticle,
127                                           G4ProcessManager*     newManager)
128{
129  if (newParticle == 0) return;
130  if (newParticle->GetProcessManager() != 0) {
131#ifdef G4VERBOSE
132    if (verboseLevel >1){
133      G4cout << "G4VUserPhysicsList::AddProcessManager: ";
134      G4cout  << newParticle->GetParticleName();
135      G4cout << " already has ProcessManager " << G4endl;
136    }
137#endif
138    return;
139  }
140
141  // create new process manager if newManager  == 0
142  if (newManager  == 0){
143    // Add ProcessManager
144    if (newParticle->GetParticleType() == "nucleus") {
145      // Create a copy of the process manager of "GenericIon" in case of "nucleus"
146      G4ParticleDefinition* genericIon =
147           (G4ParticleTable::GetParticleTable())->FindParticle("GenericIon");
148
149      if (genericIon != 0) {
150        G4ProcessManager* ionMan = genericIon->GetProcessManager();
151        if (ionMan != 0) {
152          newManager = new G4ProcessManager(*ionMan);
153        } else {
154          // no process manager has been registered yet
155          newManager = new G4ProcessManager(newParticle);
156        }
157      } else {
158        // "GenericIon" does not exist
159        newManager = new G4ProcessManager(newParticle);
160      }
161
162    } else {
163      // create process manager for particles other than "nucleus"
164      newManager = new G4ProcessManager(newParticle);
165    }
166  }
167
168  // set particle type
169  newManager->SetParticleType(newParticle);
170
171  // add the process manager
172  newParticle->SetProcessManager(newManager);
173
174#ifdef G4VERBOSE
175 if (verboseLevel >2){
176    G4cout << "G4VUserPhysicsList::AddProcessManager: ";
177    G4cout  << "adds ProcessManager to ";
178    G4cout  << newParticle->GetParticleName() << G4endl;
179    newManager->DumpInfo();
180  }
181#endif
182  if (newParticle->GetParticleType() == "nucleus") {
183    PreparePhysicsTable(newParticle);
184    BuildPhysicsTable(newParticle);
185  }
186}
187
188
189////////////////////////////////////////////////////////
190void G4VUserPhysicsList::InitializeProcessManager()
191{
192  // loop over all particles in G4ParticleTable
193  theParticleIterator->reset();
194  while( (*theParticleIterator)() ){
195    G4ParticleDefinition* particle = theParticleIterator->value();
196    G4ProcessManager* pmanager = particle->GetProcessManager();
197    if  (pmanager==0) {
198      // create process manager if the particle has no its one
199      pmanager = new G4ProcessManager(particle);
200      particle->SetProcessManager(pmanager);
201    }
202  }
203}
204
205////////////////////////////////////////////////////////
206void G4VUserPhysicsList::RemoveProcessManager()
207{
208  // loop over all particles in G4ParticleTable
209  theParticleIterator->reset();
210  while( (*theParticleIterator)() ){
211    G4ParticleDefinition* particle = theParticleIterator->value();
212    G4ProcessManager* pmanager = particle->GetProcessManager();
213    if  (pmanager!=0) delete pmanager;
214    particle->SetProcessManager(0);
215#ifdef G4VERBOSE
216    if (verboseLevel >2){
217      G4cout << "G4VUserPhysicsList::RemoveProcessManager: ";
218      G4cout  << "remove ProcessManager from ";
219      G4cout  << particle->GetParticleName() << G4endl;
220    }
221#endif
222  }
223}
224
225
226////////////////////////////////////////////////////////
227#include "G4Transportation.hh"
228#include "G4CoupledTransportation.hh"
229#include "G4RunManagerKernel.hh"
230#include "G4ScoringManager.hh"
231
232void G4VUserPhysicsList::AddTransportation()
233{
234  G4int verboseLevelTransport = 0;
235  G4VProcess* theTransportationProcess = 0;
236  G4int nParaWorld = G4RunManagerKernel::GetRunManagerKernel()->GetNumberOfParallelWorld();
237
238  if(nParaWorld || useCoupledTransportation || G4ScoringManager::GetScoringManagerIfExist())
239    {
240      theTransportationProcess = new G4CoupledTransportation(verboseLevelTransport);
241      G4cout << "#############################################################################" << G4endl
242             << " G4VUserPhysicsList::AddTransportation() --- G4CoupledTransportation is used " << G4endl
243             << "#############################################################################" << G4endl;
244    }
245  else
246    { theTransportationProcess = new G4Transportation(verboseLevelTransport); }
247 
248#ifdef G4VERBOSE
249    if (verboseLevel >2){
250      G4cout << "G4VUserPhysicsList::AddTransportation()  "<< G4endl;
251    }
252#endif
253
254  // loop over all particles in G4ParticleTable
255  theParticleIterator->reset();
256  while( (*theParticleIterator)() ){
257    G4ParticleDefinition* particle = theParticleIterator->value();
258    G4ProcessManager* pmanager = particle->GetProcessManager();
259    if (!particle->IsShortLived()) {
260      // Add transportation process for all particles other than  "shortlived"
261      if ( pmanager == 0) {
262        // Error !! no process manager
263        G4String particleName = particle->GetParticleName();
264        G4Exception("G4VUserPhysicsList::AddTransportation","No process manager",
265                    RunMustBeAborted, particleName );
266      } else {
267        // add transportation with ordering = ( -1, "first", "first" )
268        pmanager ->AddProcess(theTransportationProcess);
269        pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxAlongStep);
270        pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxPostStep);
271      }
272    } else {
273      // shortlived particle case
274     // Add transportation process for all particles other than  "shortlived"
275      if ( pmanager == 0) {
276        // Error !! no process manager
277        G4String particleName = particle->GetParticleName();
278        G4Exception("G4VUserPhysicsList::AddTransportation","No process manager",
279                    RunMustBeAborted, particleName );
280      } else {
281        // add transportation with ordering = ( -1, "first", "first" )
282        pmanager ->AddProcess(theTransportationProcess);
283        pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxAlongStep);
284        pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxPostStep);
285      }
286
287    }
288  }
289}
290
291
292////////////////////////////////////////////////////////
293void G4VUserPhysicsList::SetDefaultCutValue(G4double value)
294{
295   if (value<=0.0) {
296#ifdef G4VERBOSE
297     if (verboseLevel >0){
298       G4cout << "G4VUserPhysicsList::SetDefaultCutValue: negative cut values";
299       G4cout << "  :" << value/mm << "[mm]" << G4endl;
300     }
301#endif
302   } else {
303#ifdef G4VERBOSE
304     if (verboseLevel >1){
305       G4cout << "G4VUserPhysicsList::SetDefaultCutValue:";
306       G4cout << "default cut value is changed to   :" ;
307       G4cout << value/mm << "[mm]" << G4endl;
308     }
309#endif
310     defaultCutValue = value;
311
312   }
313}
314
315
316////////////////////////////////////////////////////////
317void G4VUserPhysicsList::SetCutValue(G4double aCut, const G4String& name)
318{
319  G4ParticleDefinition* particle = theParticleTable->FindParticle(name);
320  if (particle != 0){
321    if (!particle->IsShortLived()) {
322      //set cut value
323      SetParticleCuts( aCut ,particle );
324    }
325  }
326}
327
328////////////////////////////////////////////////////////
329void G4VUserPhysicsList::SetCutValue
330(G4double aCut, const G4String& pname, const G4String& rname)
331{
332  G4ParticleDefinition* particle = theParticleTable->FindParticle(pname);
333  G4Region* region = G4RegionStore::GetInstance()->GetRegion(rname);
334  if (particle != 0 && region != 0){
335    if (!particle->IsShortLived()) {
336      //set cut value
337      SetParticleCuts( aCut ,particle, region );
338    }
339  }
340}
341
342
343
344////////////////////////////////////////////////////////
345void G4VUserPhysicsList::SetCutsWithDefault()
346{
347  // default cut value
348  G4double cut = defaultCutValue;
349
350#ifdef G4VERBOSE
351  if (verboseLevel >1){
352    G4cout << "G4VUserPhysicsList::SetCutsWithDefault:";
353    G4cout << "CutLength : " << cut/mm << " (mm)" << G4endl;
354  }
355#endif
356
357  // set cut values for gamma at first and for e- and e+
358  SetCutValue(cut, "gamma");
359  SetCutValue(cut, "e-");
360  SetCutValue(cut, "e+");
361
362  // dump Cut values if verboseLevel==3
363  if (verboseLevel>2) {
364    DumpCutValuesTable();
365  }
366}
367
368
369////////////////////////////////////////////////////////
370void G4VUserPhysicsList::SetCutsForRegion(G4double aCut, const G4String& rname)
371{
372  // set cut values for gamma at first and for e- and e+
373  SetCutValue(aCut, "gamma", rname);
374  SetCutValue(aCut, "e-", rname);
375  SetCutValue(aCut, "e+", rname);
376}
377
378
379
380////////////////////////////////////////////////////////
381////////////////////////////////////////////////////////
382void G4VUserPhysicsList::SetParticleCuts( G4double cut, G4ParticleDefinition* particle, G4Region* region)
383{
384  if(!region) region = (*(G4RegionStore::GetInstance()))[0];
385  G4ProductionCuts* pcuts = region->GetProductionCuts();
386  pcuts->SetProductionCut(cut,particle);
387}
388
389///////////////////////////////////////////////////////////////
390void G4VUserPhysicsList::BuildPhysicsTable()
391{
392  //Prepare Physics table for all particles
393  theParticleIterator->reset();
394  while( (*theParticleIterator)() ){
395    G4ParticleDefinition* particle = theParticleIterator->value();
396    PreparePhysicsTable(particle); 
397  }
398
399  // ask processes to prepare physics table
400  if (fRetrievePhysicsTable) {
401    fIsRestoredCutValues = fCutsTable->RetrieveCutsTable(directoryPhysicsTable, fStoredInAscii);
402    // check if retrieve Cut Table successfully
403    if (!fIsRestoredCutValues) {
404#ifdef G4VERBOSE
405      if (verboseLevel>0){
406        G4cout << "G4VUserPhysicsList::BuildPhysicsTable";
407        G4cout << " Retrieve Cut Table failed !!" << G4endl;
408      }
409#endif 
410      G4Exception("G4VUserPhysicsList::BuildPhysicsTable","Fail to Retreive",
411                  RunMustBeAborted,"Production Cut Table can not be retreived");
412    } else {
413#ifdef G4VERBOSE
414      if (verboseLevel>2){
415        G4cout << "G4VUserPhysicsList::BuildPhysicsTable";
416        G4cout << "  Retrieve Cut Table successfully " << G4endl;
417      }
418#endif
419    }     
420  } else {
421#ifdef G4VERBOSE
422    if (verboseLevel>2){
423      G4cout << "G4VUserPhysicsList::BuildPhysicsTable";
424      G4cout << " does not retrieve Cut Table but calculate " << G4endl;
425    } 
426#endif     
427  }
428
429  // Sets a value to particle
430  // set cut values for gamma at first and for e- and e+
431  G4String particleName;
432  G4ParticleDefinition* GammaP = theParticleTable->FindParticle("gamma");
433  if(GammaP) BuildPhysicsTable(GammaP);
434  G4ParticleDefinition* EMinusP = theParticleTable->FindParticle("e-");
435  if(EMinusP) BuildPhysicsTable(EMinusP);
436  G4ParticleDefinition* EPlusP = theParticleTable->FindParticle("e+");
437  if(EPlusP) BuildPhysicsTable(EPlusP);
438  G4ParticleDefinition* ProtonP = theParticleTable->FindParticle("proton");
439  if(ProtonP) BuildPhysicsTable(ProtonP);
440
441
442  theParticleIterator->reset();
443  while( (*theParticleIterator)() ){
444    G4ParticleDefinition* particle = theParticleIterator->value();
445    if( particle!=GammaP && 
446        particle!=EMinusP && 
447        particle!=EPlusP && 
448        particle!=ProtonP   ){
449      BuildPhysicsTable(particle); 
450    }
451  }
452}
453///////////////////////////////////////////////////////////////
454void G4VUserPhysicsList::BuildPhysicsTable(G4ParticleDefinition* particle)
455{
456  if (fRetrievePhysicsTable) {
457    if ( !fIsRestoredCutValues){
458      // fail to retreive cut tables
459#ifdef G4VERBOSE
460      if (verboseLevel>0){
461        G4cout << "G4VUserPhysicsList::BuildPhysicsTable  ";
462        G4cout << "Physics table can not be retreived and will be calculated "<< G4endl;
463      }
464#endif
465      fRetrievePhysicsTable = false; 
466
467    } else {
468#ifdef G4VERBOSE
469      if (verboseLevel>2){
470        G4cout << "G4VUserPhysicsList::BuildPhysicsTable  ";
471        G4cout << " Retrieve Physics Table for ";
472        G4cout << particle->GetParticleName() << G4endl;
473      }
474#endif
475      //  Retrieve PhysicsTable from files for proccesses
476      RetrievePhysicsTable(particle, directoryPhysicsTable, fStoredInAscii);
477    }
478  }
479
480#ifdef G4VERBOSE
481  if (verboseLevel>2){
482    G4cout << "G4VUserPhysicsList::BuildPhysicsTable  ";
483    G4cout << "Calculate Physics Table for " << particle->GetParticleName() << G4endl;
484  }
485#endif
486  // Rebuild the physics tables for every process for this particle type
487  // if particle is not ShortLived
488  if(!particle->IsShortLived()) {
489    G4ProcessVector* pVector = particle->GetProcessManager()->GetProcessList();
490    for (G4int j=0; j < pVector->size(); ++j) {
491      (*pVector)[j]->BuildPhysicsTable(*particle);
492    }
493  }
494}
495
496///////////////////////////////////////////////////////////////
497void G4VUserPhysicsList::PreparePhysicsTable(G4ParticleDefinition* particle)
498{
499  // Prepare the physics tables for every process for this particle type
500  // if particle is not ShortLived
501  if(!particle->IsShortLived()) {
502    G4ProcessManager* pManager =  particle->GetProcessManager();
503    if (!pManager) {
504      G4cerr << "G4VUserPhysicsList::PreparePhysicsTable  : No Process Manager for " 
505             << particle->GetParticleName() <<G4endl;
506      G4cerr << particle->GetParticleName() << " should be created in your PhysicsList" <<G4endl;
507      G4Exception("G4VUserPhysicsList::PreparePhysicsTable","No process manager",
508                    RunMustBeAborted,  particle->GetParticleName() );
509    }
510
511    G4ProcessVector* pVector = pManager->GetProcessList();
512    for (G4int j=0; j < pVector->size(); ++j) {
513      (*pVector)[j]->PreparePhysicsTable(*particle);
514    }
515  }
516}
517
518///////////////////////////////////////////////////////////////
519void  G4VUserPhysicsList::BuildIntegralPhysicsTable(G4VProcess* process,
520                                                    G4ParticleDefinition* particle)
521{
522  //*********************************************************************
523  // temporary addition to make the integral schema of electromagnetic
524  // processes work.
525  //
526
527  if ( (process->GetProcessName() == "Imsc") ||
528       (process->GetProcessName() == "IeIoni") ||
529       (process->GetProcessName() == "IeBrems") ||
530       (process->GetProcessName() == "Iannihil") ||
531       (process->GetProcessName() == "IhIoni") ||
532       (process->GetProcessName() == "IMuIoni") ||
533       (process->GetProcessName() == "IMuBrems") ||
534       (process->GetProcessName() == "IMuPairProd")  ) {
535#ifdef G4VERBOSE
536    if (verboseLevel>2){
537      G4cout << "G4VUserPhysicsList::BuildIntegralPhysicsTable  ";
538      G4cout << " BuildPhysicsTable is invoked for ";
539      G4cout << process->GetProcessName();
540      G4cout << "(" << particle->GetParticleName() << ")" << G4endl;
541    }
542#endif
543    process->BuildPhysicsTable(*particle);
544  }
545}
546
547///////////////////////////////////////////////////////////////
548void G4VUserPhysicsList::DumpList() const
549{
550  theParticleIterator->reset();
551  G4int idx = 0;
552  while( (*theParticleIterator)() ){
553    G4ParticleDefinition* particle = theParticleIterator->value();
554    G4cout << particle->GetParticleName();
555    if ((idx++ % 4) == 3) {
556      G4cout << G4endl;
557    } else {
558      G4cout << ", ";
559    }     
560  }
561  G4cout << G4endl;
562}
563
564
565
566///////////////////////////////////////////////////////////////
567void G4VUserPhysicsList::DumpCutValuesTable(G4int nParticles) 
568{ 
569  fDisplayThreshold = nParticles; 
570}
571
572///////////////////////////////////////////////////////////////
573void G4VUserPhysicsList::DumpCutValuesTableIfRequested()
574{
575  if(fDisplayThreshold==0) return;
576  G4ProductionCutsTable::GetProductionCutsTable()->DumpCouples();
577  fDisplayThreshold = 0;
578}
579
580
581///////////////////////////////////////////////////////////////
582///////////////////////////////////////////////////////////////
583G4bool G4VUserPhysicsList::StorePhysicsTable(const G4String& directory)
584{
585  G4bool   ascii = fStoredInAscii;
586  G4String dir   = directory;
587  if (dir.isNull()) dir = directoryPhysicsTable; 
588  else directoryPhysicsTable = dir; 
589
590  // store CutsTable info
591  if (!fCutsTable->StoreCutsTable(dir, ascii)) {
592    G4Exception("G4VUserPhysicsList::StorePhysicsTable","Faile to store ",
593                JustWarning,"Cut Table can not be stored");     
594    return false;
595  }
596#ifdef G4VERBOSE 
597  if (verboseLevel>2){
598    G4cout << "G4VUserPhysicsList::StorePhysicsTable   ";
599    G4cout << " Store material and cut values successfully" << G4endl;
600  }
601#endif
602
603  G4bool success= true;
604
605  // loop over all particles in G4ParticleTable
606  theParticleIterator->reset();
607  while( (*theParticleIterator)() ){
608    G4ParticleDefinition* particle = theParticleIterator->value();
609    // Store physics tables for every process for this particle type
610    G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
611    G4int  j;
612    for ( j=0; j < pVector->size(); ++j) {
613      if (!(*pVector)[j]->StorePhysicsTable(particle,dir,ascii)){   
614        G4String comment =  "Fail to store for " + (*pVector)[j]->GetProcessName();
615        comment += "(" + particle->GetParticleName()  + ")";
616        G4Exception("G4VUserPhysicsList::StorePhysicsTable","Faile to store ",
617                    JustWarning,comment);       
618        success = false;
619      }
620    }
621    // end loop over processes
622  }
623  // end loop over particles
624  return success;
625}
626
627
628
629///////////////////////////////////////////////////////////////
630void  G4VUserPhysicsList::SetPhysicsTableRetrieved(const G4String& directory)
631{
632  fRetrievePhysicsTable = true;
633  if(!directory.isNull()) {
634    directoryPhysicsTable = directory;
635  }
636  fIsCheckedForRetrievePhysicsTable=false;
637  fIsRestoredCutValues = false;
638}
639
640///////////////////////////////////////////////////////////////
641void G4VUserPhysicsList::RetrievePhysicsTable(G4ParticleDefinition* particle, 
642                                              const G4String& directory,
643                                              G4bool          ascii)
644{
645  G4int  j;
646  G4bool success[100]; 
647  // Retrieve physics tables for every process for this particle type
648  G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
649  for ( j=0; j < pVector->size(); ++j) {
650    success[j] = 
651       (*pVector)[j]->RetrievePhysicsTable(particle,directory,ascii);
652
653    if (!success[j]) {
654#ifdef G4VERBOSE 
655      if (verboseLevel>2){
656        G4cout << "G4VUserPhysicsList::RetrievePhysicsTable   ";
657        G4cout << " Fail to retrieve Physics Table for ";
658        G4cout << (*pVector)[j]->GetProcessName() << G4endl;
659        G4cout << "Calculate Physics Table for " << particle->GetParticleName() << G4endl;
660      }
661#endif
662      (*pVector)[j]->BuildPhysicsTable(*particle);
663    }
664  }
665  for ( j=0; j < pVector->size(); ++j) {
666    // temporary addition to make the integral schema
667    if (!success[j]) BuildIntegralPhysicsTable((*pVector)[j], particle); 
668  }
669}
670
671void G4VUserPhysicsList::ResetCuts()
672{
673#ifdef G4VERBOSE 
674  if (verboseLevel>0){
675    G4cout << "G4VUserPhysicsList::ResetCuts() is obsolete.";
676    G4cout << " This method gives no effect and you can remove it. "<< G4endl;
677  }
678#endif
679}
680
681void G4VUserPhysicsList::SetApplyCuts(G4bool value, const G4String& name)
682{
683#ifdef G4VERBOSE 
684  if (verboseLevel>2){
685    G4cout << "G4VUserPhysicsList::SetApplyCuts for " << name << G4endl;
686  }
687#endif
688  if(name=="all") {
689    theParticleTable->FindParticle("gamma")->SetApplyCutsFlag(value);
690    theParticleTable->FindParticle("e-")->SetApplyCutsFlag(value);
691    theParticleTable->FindParticle("e+")->SetApplyCutsFlag(value);
692  } else {
693    theParticleTable->FindParticle(name)->SetApplyCutsFlag(value);
694  }
695}
696
697G4bool G4VUserPhysicsList::GetApplyCuts(const G4String& name) const
698{
699  return theParticleTable->FindParticle(name)->GetApplyCutsFlag();
700}
701
702
703
704
705/// obsolete methods
706
707
708void G4VUserPhysicsList::DumpCutValues(const G4String &particle_name)
709{
710  G4cerr << "WARNING !" << G4endl;
711  G4cerr << " Obsolete DumpCutValues() method is invoked for " << particle_name << G4endl;
712  G4cerr << " Please use DumpCutValuesTable() instead." << G4endl;
713  G4cerr << " This dummy method implementation will be removed soon." << G4endl;
714  DumpCutValuesTable();
715}
716
717void G4VUserPhysicsList::DumpCutValues(G4ParticleDefinition* )
718{
719  G4cerr << "WARNING !" << G4endl;
720  G4cerr << " DumpCutValues() became obsolete." << G4endl;
721  G4cerr << " Please use DumpCutValuesTable() instead." << G4endl;
722  G4cerr << " This dummy method implementation will be removed soon." << G4endl;
723  DumpCutValuesTable();
724}
Note: See TracBrowser for help on using the repository browser.