Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VUserPhysicsList.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4VUserPhysicsList.cc 101155 2016-11-08 08:21:41Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // ------------------------------------------------------------
34 // History
35 // first version 09 Jan 1998 by H.Kurashige
36 // Added SetEnergyRange 18 Jun 1998 by H.Kurashige
37 // Change for short lived particles 27 Jun 1998 by H.Kurashige
38 // G4BestUnit on output 12 nov 1998 by M.Maire
39 // Added RemoveProcessManager 9 Feb 1999 by H.Kurashige
40 // Fixed RemoveProcessManager 15 Apr 1999 by H.Kurashige
41 // Removed ConstructAllParticles() 15 Apr 1999 by H.Kurashige
42 // Modified for CUTS per REGION 10 Oct 2002 by H.Kurashige
43 // Check if particle IsShortLived 18 Jun 2003 by V.Ivanchenko
44 // Modify PreparePhysicsList 18 Jan 2006 by H.Kurashige
45 // Added PhysicsListHelper 29 APr. 2011 H.Kurashige
46 // Added default impelmentation of SetCuts 10 June 2011 H.Kurashige
47 // SetCuts is not 'pure virtual' any more
48 // Transformation for G4MT 26 Mar 2013 A. Dotti
49 // PL is shared by threads. Adding a method for workers
50 // To initialize thread specific data
51 // ------------------------------------------------------------
52 
53 #include <iomanip>
54 #include <fstream>
55 
56 #include "G4PhysicsListHelper.hh"
57 #include "G4VUserPhysicsList.hh"
58 
59 //Andrea Dotti (Jan 13, 2013), transformation for G4MT
60 #include "G4VMultipleScattering.hh"
61 #include "G4VEnergyLossProcess.hh"
62 
63 
64 #include "globals.hh"
65 #include "G4SystemOfUnits.hh"
66 #include "G4ios.hh"
67 #include "G4ParticleDefinition.hh"
68 #include "G4ProcessManager.hh"
69 #include "G4ParticleTable.hh"
70 #include "G4ProductionCutsTable.hh"
71 #include "G4Material.hh"
73 #include "G4UImanager.hh"
74 #include "G4UnitsTable.hh"
75 #include "G4RegionStore.hh"
76 #include "G4Region.hh"
77 #include "G4ProductionCutsTable.hh"
78 #include "G4ProductionCuts.hh"
79 #include "G4MaterialCutsCouple.hh"
80 
81 // This static member is thread local. For each thread, it holds the array
82 // size of G4VUPLData instances.
83 //
84 #define G4MT_theMessenger ((this->subInstanceManager.offset[this->g4vuplInstanceID])._theMessenger)
85 #define G4MT_thePLHelper ((this->subInstanceManager.offset[this->g4vuplInstanceID])._thePLHelper)
86 #define fIsPhysicsTableBuilt ((this->subInstanceManager.offset[this->g4vuplInstanceID])._fIsPhysicsTableBuilt)
87 #define fDisplayThreshold ((this->subInstanceManager.offset[this->g4vuplInstanceID])._fDisplayThreshold)
88 #define theParticleIterator ((this->subInstanceManager.offset[this->g4vuplInstanceID])._theParticleIterator)
89 
90 // This field helps to use the class G4VUPLManager
91 //
93 
95 {
97  _theMessenger = 0;
99  _fIsPhysicsTableBuilt = false;
100  _fDisplayThreshold = 0;
101 }
102 
105  :verboseLevel(1),
106  defaultCutValue(1.0 * mm),
107  isSetDefaultCutValue(false),
108  fRetrievePhysicsTable(false),
109  fStoredInAscii(true),
110  fIsCheckedForRetrievePhysicsTable(false),
111  fIsRestoredCutValues(false),
112  directoryPhysicsTable("."),
113  //fDisplayThreshold(0),
114  //fIsPhysicsTableBuilt(false),
115  fDisableCheckParticleList(false)
116 {
118  // default cut value (1.0mm)
119  defaultCutValue = 1.0*mm;
120 
121  // pointer to the particle table
123  //theParticleIterator = theParticleTable->GetIterator();
124 
125  // pointer to the cuts table
127 
128  // set energy range for SetCut calcuration
129  fCutsTable->SetEnergyRange(0.99*keV, 100*TeV);
130 
131  // UI Messenger
132  //theMessenger = new G4UserPhysicsListMessenger(this);
134 
135  // PhysicsListHelper
136  //thePLHelper = G4PhysicsListHelper::GetPhysicsListHelper();
137  //thePLHelper->SetVerboseLevel(verboseLevel);
138  //G4MT_thePLHelper = G4PhysicsListHelper::GetPhysicsListHelper(); //AND
139  G4MT_thePLHelper->SetVerboseLevel(verboseLevel); //AND
140 
141  fIsPhysicsTableBuilt = false;
142  fDisplayThreshold = 0;
143 
144 }
145 
147 {
148  //Remember messengers are per-thread, so this needs to be done by each worker
149  //and due to the presence of "this" cannot be done in G4VUPLData::initialize()
151 }
152 
155 {
156  if (G4MT_theMessenger != 0) {
157  delete G4MT_theMessenger;
158  G4MT_theMessenger = 0;
159  }
161 
162  // invoke DeleteAllParticle
164 
165 }
166 
169  :verboseLevel(right.verboseLevel),
170  defaultCutValue(right.defaultCutValue),
171  isSetDefaultCutValue(right.isSetDefaultCutValue),
172  fRetrievePhysicsTable(right.fRetrievePhysicsTable),
173  fStoredInAscii(right.fStoredInAscii),
174  fIsCheckedForRetrievePhysicsTable(right.fIsCheckedForRetrievePhysicsTable),
175  fIsRestoredCutValues(right.fIsRestoredCutValues),
176  directoryPhysicsTable(right.directoryPhysicsTable),
177  //fDisplayThreshold(right.fDisplayThreshold),
178  //fIsPhysicsTableBuilt(right.fIsPhysicsTableBuilt),
179  fDisableCheckParticleList(right.fDisableCheckParticleList)
180 {
182  // pointer to the particle table
185  // pointer to the cuts table
187 
188  // UI Messenger
189  //theMessenger = new G4UserPhysicsListMessenger(this);
191 
192  // PhysicsListHelper
193  //thePLHelper = G4PhysicsListHelper::GetPhysicsListHelper();
194  //thePLHelper->SetVerboseLevel(verboseLevel);
196  G4MT_thePLHelper->SetVerboseLevel(verboseLevel); //AND
197 
198  fIsPhysicsTableBuilt = right.GetSubInstanceManager().offset[right.GetInstanceID()]._fIsPhysicsTableBuilt;
199  fDisplayThreshold = right.GetSubInstanceManager().offset[right.GetInstanceID()]._fDisplayThreshold;
200 }
201 
202 
205 {
206  if (this != &right) {
207  verboseLevel = right.verboseLevel;
215  //fDisplayThreshold = right.fDisplayThreshold;
216  fIsPhysicsTableBuilt = right.GetSubInstanceManager().offset[right.GetInstanceID()]._fIsPhysicsTableBuilt;
217  fDisplayThreshold = right.GetSubInstanceManager().offset[right.GetInstanceID()]._fDisplayThreshold;
218  //fIsPhysicsTableBuilt = right.fIsPhysicsTableBuilt;
220  }
221  return *this;
222 }
223 
227 {
228  if (newParticle == 0) return;
229  G4Exception("G4VUserPhysicsList::AddProcessManager",
230  "Run0252", JustWarning,
231  "This method is obsolete");
232 }
233 
234 
237 {
238  //Request lock for particle table accesses. Some changes are inside
239  //this critical region.
240 #ifdef G4MULTITHREADED
241  G4MUTEXLOCK(&G4ParticleTable::particleTableMutex);
242  G4ParticleTable::lockCount++;
243 #endif
245 
246  // loop over all particles in G4ParticleTable
247  theParticleIterator->reset();
248  while( (*theParticleIterator)() ){
249  G4ParticleDefinition* particle = theParticleIterator->value();
250  G4ProcessManager* pmanager = particle->GetProcessManager();
251 
252  if (pmanager==0) {
253  // create process manager if the particle does not have its own.
254  pmanager = new G4ProcessManager(particle);
255  particle->SetProcessManager(pmanager);
256  if( particle->GetMasterProcessManager() == 0 ) particle->SetMasterProcessManager(pmanager);
257 #ifdef G4VERBOSE
258  if (verboseLevel >2){
259  G4cout << "G4VUserPhysicsList::InitializeProcessManager: creating ProcessManager to "
260  << particle->GetParticleName() << G4endl;
261  }
262 #endif
263  }
264  }
265 
266  if(gion)
267  {
268  G4ProcessManager* gionPM = gion->GetProcessManager();
269  // loop over all particles once again (this time, with all general ions)
270  theParticleIterator->reset(false);
271  while( (*theParticleIterator)() ){
272  G4ParticleDefinition* particle = theParticleIterator->value();
273  if(particle->IsGeneralIon())
274  {
275  particle->SetProcessManager(gionPM);
276 #ifdef G4VERBOSE
277  if (verboseLevel >2){
278  G4cout << "G4VUserPhysicsList::InitializeProcessManager: copying ProcessManager to "
279  << particle->GetParticleName() << G4endl;
280  }
281 #endif
282  }
283  }
284  }
285 
286  //release lock for particle table accesses.
287 #ifdef G4MULTITHREADED
288  G4MUTEXUNLOCK(&G4ParticleTable::particleTableMutex);
289 #endif
290 // G4cout << "Particle table is released by G4VUserPhysicsList::InitializeProcessManager" << G4endl;
291 
292 }
293 
296 {
297  //Request lock for particle table accesses. Some changes are inside
298  //this critical region.
299 #ifdef G4MULTITHREADED
300  G4MUTEXLOCK(&G4ParticleTable::particleTableMutex);
301  G4ParticleTable::lockCount++;
302 #endif
303 // G4cout << "Particle table is held by G4VUserPhysicsList::InitializeProcessManager" << G4endl;
304 
305  // loop over all particles in G4ParticleTable
306  theParticleIterator->reset();
307  while( (*theParticleIterator)() ){
308  G4ParticleDefinition* particle = theParticleIterator->value();
310  {
311  if(particle->GetParticleSubType()!="generic" || particle->GetParticleName()=="GenericIon")
312  {
313  G4ProcessManager* pmanager = particle->GetProcessManager();
314  if (pmanager!=0) delete pmanager;
315 #ifdef G4VERBOSE
316  if (verboseLevel >2){
317  G4cout << "G4VUserPhysicsList::RemoveProcessManager: ";
318  G4cout << "remove ProcessManager from ";
319  G4cout << particle->GetParticleName() << G4endl;
320  }
321 #endif
322  }
323  particle->SetProcessManager(0);
324  }
325  }
326 
327  //release lock for particle table accesses.
328 #ifdef G4MULTITHREADED
329  G4MUTEXUNLOCK(&G4ParticleTable::particleTableMutex);
330 #endif
331 // G4cout << "Particle table is released by G4VUserPhysicsList::InitializeProcessManager" << G4endl;
332 
333 }
334 
337 {
338  if ( !isSetDefaultCutValue ){
340  }
341 
342 #ifdef G4VERBOSE
343  if (verboseLevel >1){
344  G4cout << "G4VUserPhysicsList::SetCuts: " << G4endl;
345  G4cout << "Cut for gamma: " << GetCutValue("gamma")/mm
346  << "[mm]" << G4endl;
347  G4cout << "Cut for e-: " << GetCutValue("e-")/mm
348  << "[mm]" << G4endl;
349  G4cout << "Cut for e+: " << GetCutValue("e+")/mm
350  << "[mm]" << G4endl;
351  G4cout << "Cut for proton: " << GetCutValue("proton")/mm
352  << "[mm]" << G4endl;
353  }
354 #endif
355 
356  // dump Cut values if verboseLevel==3
357  if (verboseLevel>2) {
359  }
360 }
361 
362 
365 {
366  if (value<0.0) {
367 #ifdef G4VERBOSE
368  if (verboseLevel >0){
369  G4cout << "G4VUserPhysicsList::SetDefaultCutValue: negative cut values"
370  << " :" << value/mm << "[mm]" << G4endl;
371  }
372 #endif
373  return;
374  }
375 
377  isSetDefaultCutValue = true;
378 
379  // set cut values for gamma at first and for e- and e+
380  SetCutValue(defaultCutValue, "gamma");
383  SetCutValue(defaultCutValue, "proton");
384 
385 #ifdef G4VERBOSE
386  if (verboseLevel >1){
387  G4cout << "G4VUserPhysicsList::SetDefaultCutValue:"
388  << "default cut value is changed to :"
389  << defaultCutValue/mm << "[mm]" << G4endl;
390  }
391 #endif
392  }
393 
394 
397 {
398  size_t nReg = (G4RegionStore::GetInstance())->size();
399  if (nReg==0) {
400 #ifdef G4VERBOSE
401  if (verboseLevel>0){
402  G4cout << "G4VUserPhysicsList::GetCutValue "
403  <<" : No Default Region " <<G4endl;
404  }
405 #endif
406  G4Exception("G4VUserPhysicsList::GetCutValue",
407  "Run0253", FatalException,
408  "No Default Region");
409  return -1.*mm;
410  }
411  G4Region* region = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld", false);
412  return region->GetProductionCuts()->GetProductionCut(name);
413 }
414 
417 {
418  SetParticleCuts( aCut ,name );
419 }
420 
423 (G4double aCut, const G4String& pname, const G4String& rname)
424 {
425  G4Region* region = G4RegionStore::GetInstance()->GetRegion(rname);
426  if (region != 0){
427  //set cut value
428  SetParticleCuts( aCut ,pname, region );
429  } else {
430 #ifdef G4VERBOSE
431  if (verboseLevel>0){
432  G4cout << "G4VUserPhysicsList::SetCutValue "
433  <<" : No Region of " << rname << G4endl;
434  }
435 #endif
436  }
437 }
438 
439 
442 {
445 }
446 
449 {
450  // set cut values for gamma at first and for e- and e+
451  SetCutValue(aCut, "gamma", rname);
452  SetCutValue(aCut, "e-", rname);
453  SetCutValue(aCut, "e+", rname);
454  SetCutValue(aCut, "proton", rname);
455 }
456 
457 
458 
461 {
462  SetParticleCuts(cut, particle->GetParticleName(), region);
463 }
464 
466 void G4VUserPhysicsList::SetParticleCuts( G4double cut, const G4String& particleName, G4Region* region)
467 {
468  if (cut<0.0) {
469 #ifdef G4VERBOSE
470  if (verboseLevel >0){
471  G4cout << "G4VUserPhysicsList::SetParticleCuts: negative cut values"
472  << " :" << cut/mm << "[mm]"
473  << " for "<< particleName << G4endl;
474  }
475 #endif
476  return;
477  }
478 
479  G4Region* world_region = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld", false);
480  if(!region){
481  size_t nReg = (G4RegionStore::GetInstance())->size();
482  if (nReg==0) {
483 #ifdef G4VERBOSE
484  if (verboseLevel>0){
485  G4cout << "G4VUserPhysicsList::SetParticleCuts "
486  <<" : No Default Region " <<G4endl;
487  }
488 #endif
489  G4Exception("G4VUserPhysicsList::SetParticleCuts ",
490  "Run0254", FatalException,
491  "No Default Region");
492  return;
493  }
494  region = world_region;
495  }
496 
497  if ( !isSetDefaultCutValue ){
499  }
500 
501  G4ProductionCuts* pcuts = region->GetProductionCuts();
502  if(region != world_region &&
503  pcuts==G4ProductionCutsTable::GetProductionCutsTable()->GetDefaultProductionCuts())
504  { // This region had no unique cuts yet but shares the default cuts.
505  // Need to create a new object before setting the value.
506  pcuts = new G4ProductionCuts(
507  *(G4ProductionCutsTable::GetProductionCutsTable()->GetDefaultProductionCuts()));
508  region->SetProductionCuts(pcuts);
509  }
510  pcuts->SetProductionCut(cut,particleName);
511 #ifdef G4VERBOSE
512  if (verboseLevel>2){
513  G4cout << "G4VUserPhysicsList::SetParticleCuts: "
514  << " :" << cut/mm << "[mm]"
515  << " for "<< particleName << G4endl;
516  }
517 #endif
518 }
519 
522 {
523  //Prepare Physics table for all particles
524  theParticleIterator->reset();
525  while( (*theParticleIterator)() ){
526  G4ParticleDefinition* particle = theParticleIterator->value();
527  PreparePhysicsTable(particle);
528  }
529 
530  // ask processes to prepare physics table
531  if (fRetrievePhysicsTable) {
533  // check if retrieve Cut Table successfully
534  if (!fIsRestoredCutValues) {
535 #ifdef G4VERBOSE
536  if (verboseLevel>0){
537  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
538  << " Retrieve Cut Table failed !!" << G4endl;
539  }
540 #endif
541  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
542  "Run0255", RunMustBeAborted,
543  "Fail to retrieve Production Cut Table");
544  } else {
545 #ifdef G4VERBOSE
546  if (verboseLevel>2){
547  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
548  << " Retrieve Cut Table successfully " << G4endl;
549  }
550 #endif
551  }
552  } else {
553 #ifdef G4VERBOSE
554  if (verboseLevel>2){
555  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
556  << " does not retrieve Cut Table but calculate " << G4endl;
557  }
558 #endif
559  }
560 
561  // Sets a value to particle
562  // set cut values for gamma at first and for e- and e+
563  G4String particleName;
565  if(GammaP) BuildPhysicsTable(GammaP);
567  if(EMinusP) BuildPhysicsTable(EMinusP);
569  if(EPlusP) BuildPhysicsTable(EPlusP);
570  G4ParticleDefinition* ProtonP = theParticleTable->FindParticle("proton");
571  if(ProtonP) BuildPhysicsTable(ProtonP);
572 
573  theParticleIterator->reset();
574  while( (*theParticleIterator)() ){
575  G4ParticleDefinition* particle = theParticleIterator->value();
576  if( particle!=GammaP &&
577  particle!=EMinusP &&
578  particle!=EPlusP &&
579  particle!=ProtonP ){
580  BuildPhysicsTable(particle);
581  }
582  }
583 
584  // Set flag
585  fIsPhysicsTableBuilt = true;
586 
587 }
589 //Change in order to share physics tables for two kind of process.
591 {
592  if(!(particle->GetMasterProcessManager())) {
593  G4cout << "#### G4VUserPhysicsList::BuildPhysicsTable() - BuildPhysicsTable("
594  << particle->GetParticleName() << ") skipped..." << G4endl;
595  return;
596  }
597  if (fRetrievePhysicsTable) {
598  if ( !fIsRestoredCutValues){
599  // fail to retreive cut tables
600 #ifdef G4VERBOSE
601  if (verboseLevel>0){
602  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
603  << "Physics table can not be retreived and will be calculated "
604  << G4endl;
605  }
606 #endif
607  fRetrievePhysicsTable = false;
608 
609  } else {
610 #ifdef G4VERBOSE
611  if (verboseLevel>2){
612  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
613  << " Retrieve Physics Table for "
614  << particle->GetParticleName() << G4endl;
615  }
616 #endif
617  // Retrieve PhysicsTable from files for proccesses
619  }
620  }
621 
622 #ifdef G4VERBOSE
623  if (verboseLevel>2){
624  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
625  << "Calculate Physics Table for "
626  << particle->GetParticleName() << G4endl;
627  }
628 #endif
629  // Rebuild the physics tables for every process for this particle type
630  // if particle is not ShortLived
631  if(!particle->IsShortLived()) {
632  G4ProcessManager* pManager = particle->GetProcessManager();
633  if (!pManager) {
634 #ifdef G4VERBOSE
635  if (verboseLevel>0){
636  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
637  <<" : No Process Manager for "
638  << particle->GetParticleName() << G4endl;
639  G4cout << particle->GetParticleName()
640  << " should be created in your PhysicsList" <<G4endl;
641  }
642 #endif
643  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
644  "Run0271", FatalException,
645  "No process manager");
646  return;
647  }
648 
649  //Get processes from master thread;
650  G4ProcessManager* pManagerShadow = particle->GetMasterProcessManager();
651 
652  G4ProcessVector* pVector = pManager->GetProcessList();
653  if (!pVector) {
654 #ifdef G4VERBOSE
655  if (verboseLevel>0){
656  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
657  <<" : No Process Vector for "
658  << particle->GetParticleName() <<G4endl;
659  }
660 #endif
661  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
662  "Run0272", FatalException,
663  "No process Vector");
664  return;
665  }
666 #ifdef G4VERBOSE
667  if (verboseLevel>2){
668  G4cout << "G4VUserPhysicsList::BuildPhysicsTable %%%%%% " << particle->GetParticleName() << G4endl;
669  G4cout << " ProcessManager : " << pManager << " ProcessManagerShadow : " << pManagerShadow << G4endl;
670  for(G4int iv1=0;iv1<pVector->size();iv1++)
671  { G4cout << " " << iv1 << " - " << (*pVector)[iv1]->GetProcessName() << G4endl; }
672  G4cout << "--------------------------------------------------------------" << G4endl;
673  G4ProcessVector* pVectorShadow = pManagerShadow->GetProcessList();
674 
675  for(G4int iv2=0;iv2<pVectorShadow->size();iv2++)
676  { G4cout << " " << iv2 << " - " << (*pVectorShadow)[iv2]->GetProcessName() << G4endl; }
677  }
678 #endif
679  for (G4int j=0; j < pVector->size(); ++j) {
680  //Andrea July 16th 2013 : migration to new interface...
681  //Infer if we are in a worker thread or master thread
682  //Master thread is the one in which the process manager
683  // and process manager shadow pointers are the same
684  if ( pManagerShadow == pManager )
685  {
686  (*pVector)[j]->BuildPhysicsTable(*particle);
687  }
688  else
689  {
690  (*pVector)[j]->BuildWorkerPhysicsTable(*particle);
691  }
692 
693  } //End loop on processes vector
694  } //End if short-lived
695 }
696 
699 {
700  if(!(particle->GetMasterProcessManager())) {
703  return;
704  }
705  // Prepare the physics tables for every process for this particle type
706  // if particle is not ShortLived
707  if(!particle->IsShortLived()) {
708  G4ProcessManager* pManager = particle->GetProcessManager();
709  if (!pManager) {
710 #ifdef G4VERBOSE
711  if (verboseLevel>0) {
712  G4cout<< "G4VUserPhysicsList::PreparePhysicsTable "
713  << ": No Process Manager for "
714  << particle->GetParticleName() <<G4endl;
715  G4cout << particle->GetParticleName()
716  << " should be created in your PhysicsList" <<G4endl;
717  }
718 #endif
719  G4Exception("G4VUserPhysicsList::PreparePhysicsTable",
720  "Run0273", FatalException,
721  "No process manager");
722  return;
723  }
724 
725  //Get processes from master thread
726  G4ProcessManager* pManagerShadow = particle->GetMasterProcessManager();
727  //Andrea Dotti 15 Jan 2013: Change of interface of MSC
728  //G4ProcessVector* pVectorShadow = pManagerShadow->GetProcessList();
729 
730  G4ProcessVector* pVector = pManager->GetProcessList();
731  if (!pVector) {
732 #ifdef G4VERBOSE
733  if (verboseLevel>0) {
734  G4cout << "G4VUserPhysicsList::PreparePhysicsTable "
735  << ": No Process Vector for "
736  << particle->GetParticleName() <<G4endl;
737  }
738 #endif
739  G4Exception("G4VUserPhysicsList::PreparePhysicsTable",
740  "Run0274", FatalException,
741  "No process Vector");
742  return;
743  }
744  for (G4int j=0; j < pVector->size(); ++j) {
745 
746  //Andrea July 16th 2013 : migration to new interface...
747  //Infer if we are in a worker thread or master thread
748  //Master thread is the one in which the process manager
749  // and process manager shadow pointers are the same
750  if ( pManagerShadow == pManager )
751  {
752  (*pVector)[j]->PreparePhysicsTable(*particle);
753  }
754  else
755  {
756  (*pVector)[j]->PrepareWorkerPhysicsTable(*particle);
757  }
758  } //End loop on processes vector
759  } //End if pn ShortLived
760 }
761 
762 //TODO Should we change this function?
765  G4ParticleDefinition* particle)
766 {
767  //*********************************************************************
768  // temporary addition to make the integral schema of electromagnetic
769  // processes work.
770  //
771 
772  if ( (process->GetProcessName() == "Imsc") ||
773  (process->GetProcessName() == "IeIoni") ||
774  (process->GetProcessName() == "IeBrems") ||
775  (process->GetProcessName() == "Iannihil") ||
776  (process->GetProcessName() == "IhIoni") ||
777  (process->GetProcessName() == "IMuIoni") ||
778  (process->GetProcessName() == "IMuBrems") ||
779  (process->GetProcessName() == "IMuPairProd") ) {
780 #ifdef G4VERBOSE
781  if (verboseLevel>2){
782  G4cout << "G4VUserPhysicsList::BuildIntegralPhysicsTable "
783  << " BuildPhysicsTable is invoked for "
784  << process->GetProcessName()
785  << "(" << particle->GetParticleName() << ")" << G4endl;
786  }
787 #endif
788  process->BuildPhysicsTable(*particle);
789  }
790 }
791 
794 {
795  theParticleIterator->reset();
796  G4int idx = 0;
797  while( (*theParticleIterator)() ){
798  G4ParticleDefinition* particle = theParticleIterator->value();
799  G4cout << particle->GetParticleName();
800  if ((idx++ % 4) == 3) {
801  G4cout << G4endl;
802  } else {
803  G4cout << ", ";
804  }
805  }
806  G4cout << G4endl;
807 }
808 
809 
812 {
813  fDisplayThreshold = flag;
814 }
815 
818 {
819  if(fDisplayThreshold==0) return;
821  fDisplayThreshold = 0;
822 }
823 
824 
827 {
828  G4bool ascii = fStoredInAscii;
829  G4String dir = directory;
830  if (dir.isNull()) dir = directoryPhysicsTable;
831  else directoryPhysicsTable = dir;
832 
833  // store CutsTable info
834  if (!fCutsTable->StoreCutsTable(dir, ascii)) {
835  G4Exception("G4VUserPhysicsList::StorePhysicsTable",
836  "Run0281", JustWarning,
837  "Fail to store Cut Table");
838  return false;
839  }
840 #ifdef G4VERBOSE
841  if (verboseLevel>2){
842  G4cout << "G4VUserPhysicsList::StorePhysicsTable "
843  << " Store material and cut values successfully" << G4endl;
844  }
845 #endif
846 
847  G4bool success= true;
848 
849  // loop over all particles in G4ParticleTable
850  theParticleIterator->reset();
851  while( (*theParticleIterator)() ){
852  G4ParticleDefinition* particle = theParticleIterator->value();
853  // Store physics tables for every process for this particle type
854  G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
855  G4int j;
856  for ( j=0; j < pVector->size(); ++j) {
857  if (!(*pVector)[j]->StorePhysicsTable(particle,dir,ascii)){
858  G4String comment = "Fail to store physics table for ";
859  comment += (*pVector)[j]->GetProcessName();
860  comment += "(" + particle->GetParticleName() + ")";
861  G4Exception("G4VUserPhysicsList::StorePhysicsTable",
862  "Run0282", JustWarning,
863  comment);
864  success = false;
865  }
866  }
867  // end loop over processes
868  }
869  // end loop over particles
870  return success;
871 }
872 
873 
874 
877 {
878  fRetrievePhysicsTable = true;
879  if(!directory.isNull()) {
880  directoryPhysicsTable = directory;
881  }
883  fIsRestoredCutValues = false;
884 }
885 
888  const G4String& directory,
889  G4bool ascii)
890 {
891  G4int j;
892  G4bool success[100];
893  // Retrieve physics tables for every process for this particle type
894  G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
895  for ( j=0; j < pVector->size(); ++j) {
896  success[j] =
897  (*pVector)[j]->RetrievePhysicsTable(particle,directory,ascii);
898 
899  if (!success[j]) {
900 #ifdef G4VERBOSE
901  if (verboseLevel>2){
902  G4cout << "G4VUserPhysicsList::RetrievePhysicsTable "
903  << " Fail to retrieve Physics Table for "
904  << (*pVector)[j]->GetProcessName() << G4endl;
905  G4cout << "Calculate Physics Table for "
906  << particle->GetParticleName() << G4endl;
907  }
908 #endif
909  (*pVector)[j]->BuildPhysicsTable(*particle);
910  }
911  }
912  for ( j=0; j < pVector->size(); ++j) {
913  // temporary addition to make the integral schema
914  if (!success[j]) BuildIntegralPhysicsTable((*pVector)[j], particle);
915  }
916 }
917 
918 
921 {
922 #ifdef G4VERBOSE
923  if (verboseLevel>2){
924  G4cout << "G4VUserPhysicsList::SetApplyCuts for " << name << G4endl;
925  }
926 #endif
927  if(name=="all") {
932  } else {
934  }
935 }
936 
939 {
941 }
942 
943 
946 {
948  G4MT_thePLHelper->CheckParticleList();
949  }
950 }
951 
954 {
955  G4MT_thePLHelper->AddTransportation();
956 }
957 
960 {
961  G4MT_thePLHelper->UseCoupledTransportation(vl);
962 }
963 
966  G4ParticleDefinition* particle)
967 {
968  return G4MT_thePLHelper->RegisterProcess(process, particle);
969 }
970 
973 {
974  return (subInstanceManager.offset[g4vuplInstanceID])._theParticleIterator;
975 }
976 
979 {
981  // set verboseLevel for G4ProductionCutsTable same as one for G4VUserPhysicsList:
983 
984  G4MT_thePLHelper->SetVerboseLevel(verboseLevel);
985 
986 #ifdef G4VERBOSE
987  if (verboseLevel >1){
988  G4cout << "G4VUserPhysicsList::SetVerboseLevel :"
989  << " Verbose level is set to " << verboseLevel << G4endl;
990  }
991 #endif
992 }
993 
994 
997 
1000 {
1001 #ifdef G4VERBOSE
1002  if (verboseLevel>0){
1003  G4cout << "G4VUserPhysicsList::ResetCuts() is obsolete."
1004  << " This method gives no effect and you can remove it. "<< G4endl;
1005  }
1006 #endif
1007 }
#define G4MUTEXUNLOCK
Definition: G4Threading.hh:180
G4double GetCutValue(const G4String &pname) const
void SetDefaultCutValue(G4double newCutValue)
void BuildIntegralPhysicsTable(G4VProcess *, G4ParticleDefinition *)
const XML_Char * name
Definition: expat.h:151
G4ProductionCuts * GetProductionCuts() const
void SetApplyCuts(G4bool value, const G4String &name)
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4PhysicsListHelper * _thePLHelper
static constexpr double mm
Definition: G4SIunits.hh:115
void SetProcessManager(G4ProcessManager *aProcessManager)
void SetCutValue(G4double aCut, const G4String &pname)
void SetEnergyRange(G4double lowedge, G4double highedge)
void PreparePhysicsTable(G4ParticleDefinition *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
G4ProductionCutsTable * fCutsTable
G4double GetProductionCut(G4int index) const
G4bool fIsCheckedForRetrievePhysicsTable
void SetProductionCut(G4double cut, G4int index=-1)
G4bool GetApplyCuts(const G4String &name) const
void SetCutsForRegion(G4double aCut, const G4String &rname)
G4ParticleTable * theParticleTable
const G4String & GetParticleSubType() const
G4ParticleDefinition * GetGenericIon() const
G4UserPhysicsListMessenger * _theMessenger
int G4int
Definition: G4Types.hh:78
G4int GetInstanceID() const
const G4String & GetParticleName() const
G4bool IsGeneralIon() const
G4RUN_DLL G4ThreadLocalStatic T * offset
static constexpr double TeV
Definition: G4SIunits.hh:218
static G4RegionStore * GetInstance()
G4bool _fIsPhysicsTableBuilt
#define fDisplayThreshold
void SetPhysicsTableRetrieved(const G4String &directory="")
void SetParticleCuts(G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
void DumpCutValuesTable(G4int flag=1)
void AddProcessManager(G4ParticleDefinition *newParticle, G4ProcessManager *newManager=0)
G4GLOB_DLL std::ostream G4cout
G4ProcessManager * GetMasterProcessManager() const
const XML_Char int const XML_Char * value
Definition: expat.h:331
void SetVerboseLevel(G4int value)
G4bool RetrieveCutsTable(const G4String &directory, G4bool ascii=false)
bool G4bool
Definition: G4Types.hh:79
G4int GetInstanceID() const
void SetVerboseLevel(G4int value)
#define theParticleIterator
G4ParticleTable::G4PTblDicIterator * _theParticleIterator
#define G4MUTEXLOCK
Definition: G4Threading.hh:179
G4VUserPhysicsList & operator=(const G4VUserPhysicsList &)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
string pname
Definition: eplot.py:33
G4bool StorePhysicsTable(const G4String &directory=".")
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int size() const
static G4ProductionCutsTable * GetProductionCutsTable()
void UseCoupledTransportation(G4bool vl=true)
G4PART_DLL G4ThreadLocalStatic G4int slavetotalspace
static G4ParticleTable * GetParticleTable()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
G4ProcessManager * GetProcessManager() const
static G4PhysicsListHelper * GetPhysicsListHelper()
G4bool StoreCutsTable(const G4String &directory, G4bool ascii=false)
static G4RUN_DLL G4VUPLManager subInstanceManager
static const G4VUPLManager & GetSubInstanceManager()
void ResetCuts()
obsolete methods
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
#define G4endl
Definition: G4ios.hh:61
G4bool GetApplyCutsFlag() const
#define G4MT_theMessenger
void SetProductionCuts(G4ProductionCuts *cut)
G4int CreateSubInstance()
#define G4MT_thePLHelper
double G4double
Definition: G4Types.hh:76
#define fIsPhysicsTableBuilt
virtual void RetrievePhysicsTable(G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
static constexpr double keV
Definition: G4SIunits.hh:216
G4PTblDicIterator * GetIterator() const
G4bool isNull() const
void SetMasterProcessManager(G4ProcessManager *aNewPM)
G4ProcessVector * GetProcessList() const