Geant4  9.6.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$
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 // ------------------------------------------------------------
49 
50 #include "G4VUserPhysicsList.hh"
51 
52 #include "globals.hh"
53 #include "G4SystemOfUnits.hh"
54 #include "G4PhysicsListHelper.hh"
55 #include "G4ParticleDefinition.hh"
56 #include "G4ProcessManager.hh"
57 #include "G4ParticleTable.hh"
58 #include "G4ProductionCutsTable.hh"
59 #include "G4Material.hh"
61 #include "G4UImanager.hh"
62 #include "G4UnitsTable.hh"
63 #include "G4RegionStore.hh"
64 #include "G4Region.hh"
65 #include "G4ProductionCutsTable.hh"
66 #include "G4ProductionCuts.hh"
67 #include "G4MaterialCutsCouple.hh"
68 
69 #include "G4ios.hh"
70 #include <iomanip>
71 #include <fstream>
72 
75  :verboseLevel(1),
76  defaultCutValue(1.0 * mm),
77  isSetDefaultCutValue(false),
78  fRetrievePhysicsTable(false),
79  fStoredInAscii(true),
80  fIsCheckedForRetrievePhysicsTable(false),
81  fIsRestoredCutValues(false),
82  directoryPhysicsTable("."),
83  fDisplayThreshold(0),
84  fIsPhysicsTableBuilt(false),
85  fDisableCheckParticleList(false)
86 {
87  // default cut value (1.0mm)
88  defaultCutValue = 1.0*mm;
89 
90  // pointer to the particle table
93 
94  // pointer to the cuts table
96 
97  // set energy range for SetCut calcuration
98  fCutsTable->SetEnergyRange(0.99*keV, 100*TeV);
99 
100  // UI Messenger
102 
103  // PhysicsListHelper
106 
107 }
108 
111 {
112  if (theMessenger != 0) {
113  delete theMessenger;
114  theMessenger = 0;
115  }
117 
118  // invoke DeleteAllParticle
120 
121 }
122 
125  :verboseLevel(right.verboseLevel),
126  defaultCutValue(right.defaultCutValue),
127  isSetDefaultCutValue(right.isSetDefaultCutValue),
128  fRetrievePhysicsTable(right.fRetrievePhysicsTable),
129  fStoredInAscii(right.fStoredInAscii),
130  fIsCheckedForRetrievePhysicsTable(right.fIsCheckedForRetrievePhysicsTable),
131  fIsRestoredCutValues(right.fIsRestoredCutValues),
132  directoryPhysicsTable(right.directoryPhysicsTable),
133  fDisplayThreshold(right.fDisplayThreshold),
134  fIsPhysicsTableBuilt(right.fIsPhysicsTableBuilt),
135  fDisableCheckParticleList(right.fDisableCheckParticleList)
136 {
137  // pointer to the particle table
140 
141  // pointer to the cuts table
143 
144  // UI Messenger
146 
147  // PhysicsListHelper
150 
151 }
152 
153 
156 {
157  if (this != &right) {
158  verboseLevel = right.verboseLevel;
169  }
170  return *this;
171 }
172 
175  G4ProcessManager* newManager)
176 {
177  if (newParticle == 0) return;
178  if (newParticle->GetProcessManager() != 0) {
179 #ifdef G4VERBOSE
180  if (verboseLevel >1){
181  G4cout << "G4VUserPhysicsList::AddProcessManager: "
182  << newParticle->GetParticleName()
183  << " already has ProcessManager " << G4endl;
184  }
185 #endif
186  return;
187  }
188 
189  // create new process manager if newManager == 0
190  if (newManager == 0){
191  // Add ProcessManager
192  if (newParticle->GetParticleType() == "nucleus") {
193  // Create a copy of the process manager of "GenericIon" in case of "nucleus"
194  G4ParticleDefinition* genericIon =
195  (G4ParticleTable::GetParticleTable())->FindParticle("GenericIon");
196 
197  if (genericIon != 0) {
198  G4ProcessManager* ionMan = genericIon->GetProcessManager();
199  if (ionMan != 0) {
200  newManager = new G4ProcessManager(*ionMan);
201  } else {
202  // no process manager has been registered yet
203  newManager = new G4ProcessManager(newParticle);
204  G4Exception("G4VUserPhysicsList::AddProcessManager",
205  "Run0251", RunMustBeAborted,
206  "GenericIon has no ProcessMamanger");
207  }
208  } else {
209  // "GenericIon" does not exist
210  newManager = new G4ProcessManager(newParticle);
211  G4Exception("G4VUserPhysicsList::AddProcessManager",
212  "Run0252", RunMustBeAborted,
213  "GenericIon does not exist");
214  }
215 
216  } else {
217  // create process manager for particles other than "nucleus"
218  newManager = new G4ProcessManager(newParticle);
219  }
220  }
221 
222  // set particle type
223  newManager->SetParticleType(newParticle);
224 
225  // add the process manager
226  newParticle->SetProcessManager(newManager);
227 
228 #ifdef G4VERBOSE
229  if (verboseLevel >2){
230  G4cout << "G4VUserPhysicsList::AddProcessManager: "
231  << "adds ProcessManager to "
232  << newParticle->GetParticleName() << G4endl;
233  newManager->DumpInfo();
234  }
235 #endif
237  && (newParticle->GetParticleType() == "nucleus")) {
238  PreparePhysicsTable(newParticle);
239  BuildPhysicsTable(newParticle);
240  }
241 }
242 
243 
246 {
247  // loop over all particles in G4ParticleTable
249  while( (*theParticleIterator)() ){
251  G4ProcessManager* pmanager = particle->GetProcessManager();
252  if (pmanager==0) {
253  // create process manager if the particle has no its one
254  pmanager = new G4ProcessManager(particle);
255  particle->SetProcessManager(pmanager);
256  }
257  }
258 }
259 
262 {
263  // loop over all particles in G4ParticleTable
265  while( (*theParticleIterator)() ){
267  G4ProcessManager* pmanager = particle->GetProcessManager();
268  if (pmanager!=0) delete pmanager;
269  particle->SetProcessManager(0);
270 #ifdef G4VERBOSE
271  if (verboseLevel >2){
272  G4cout << "G4VUserPhysicsList::RemoveProcessManager: "
273  << "remove ProcessManager from "
274  << particle->GetParticleName() << G4endl;
275  }
276 #endif
277  }
278 }
279 
280 
283 {
284  if ( !isSetDefaultCutValue ){
286  }
287 
288 #ifdef G4VERBOSE
289  if (verboseLevel >1){
290  G4cout << "G4VUserPhysicsList::SetCuts: " << G4endl;
291  G4cout << "Cut for gamma: " << GetCutValue("gamma")/mm
292  << "[mm]" << G4endl;
293  G4cout << "Cut for e-: " << GetCutValue("e-")/mm
294  << "[mm]" << G4endl;
295  G4cout << "Cut for e+: " << GetCutValue("e+")/mm
296  << "[mm]" << G4endl;
297  G4cout << "Cut for proton: " << GetCutValue("proton")/mm
298  << "[mm]" << G4endl;
299  }
300 #endif
301 
302  // dump Cut values if verboseLevel==3
303  if (verboseLevel>2) {
305  }
306 }
307 
308 
311 {
312  if (value<0.0) {
313 #ifdef G4VERBOSE
314  if (verboseLevel >0){
315  G4cout << "G4VUserPhysicsList::SetDefaultCutValue: negative cut values"
316  << " :" << value/mm << "[mm]" << G4endl;
317  }
318 #endif
319  return;
320  }
321 
323  isSetDefaultCutValue = true;
324 
325  // set cut values for gamma at first and for e- and e+
326  SetCutValue(defaultCutValue, "gamma");
329  SetCutValue(defaultCutValue, "proton");
330 
331 #ifdef G4VERBOSE
332  if (verboseLevel >1){
333  G4cout << "G4VUserPhysicsList::SetDefaultCutValue:"
334  << "default cut value is changed to :"
335  << defaultCutValue/mm << "[mm]" << G4endl;
336  }
337 #endif
338  }
339 
340 
343 {
344  size_t nReg = (G4RegionStore::GetInstance())->size();
345  if (nReg==0) {
346 #ifdef G4VERBOSE
347  if (verboseLevel>0){
348  G4cout << "G4VUserPhysicsList::GetCutValue "
349  <<" : No Default Region " <<G4endl;
350  }
351 #endif
352  G4Exception("G4VUserPhysicsList::GetCutValue",
353  "Run0253", FatalException,
354  "No Default Region");
355  return -1.*mm;
356  }
357  G4Region* region = (*(G4RegionStore::GetInstance()))[0];
358  return region->GetProductionCuts()->GetProductionCut(name);
359 }
360 
363 {
364  SetParticleCuts( aCut ,name );
365 }
366 
369 (G4double aCut, const G4String& pname, const G4String& rname)
370 {
371  G4Region* region = G4RegionStore::GetInstance()->GetRegion(rname);
372  if (region != 0){
373  //set cut value
374  SetParticleCuts( aCut ,pname, region );
375  } else {
376 #ifdef G4VERBOSE
377  if (verboseLevel>0){
378  G4cout << "G4VUserPhysicsList::SetCutValue "
379  <<" : No Region of " << rname << G4endl;
380  }
381 #endif
382  }
383 }
384 
385 
388 {
391 }
392 
395 {
396  // set cut values for gamma at first and for e- and e+
397  SetCutValue(aCut, "gamma", rname);
398  SetCutValue(aCut, "e-", rname);
399  SetCutValue(aCut, "e+", rname);
400  SetCutValue(aCut, "proton", rname);
401 }
402 
403 
404 
407 {
408  SetParticleCuts(cut, particle->GetParticleName(), region);
409 }
410 
412 void G4VUserPhysicsList::SetParticleCuts( G4double cut, const G4String& particleName, G4Region* region)
413 {
414  if (cut<0.0) {
415 #ifdef G4VERBOSE
416  if (verboseLevel >0){
417  G4cout << "G4VUserPhysicsList::SetParticleCuts: negative cut values"
418  << " :" << cut/mm << "[mm]"
419  << " for "<< particleName << G4endl;
420  }
421 #endif
422  return;
423  }
424 
425  if(!region){
426  size_t nReg = (G4RegionStore::GetInstance())->size();
427  if (nReg==0) {
428 #ifdef G4VERBOSE
429  if (verboseLevel>0){
430  G4cout << "G4VUserPhysicsList::SetParticleCuts "
431  <<" : No Default Region " <<G4endl;
432  }
433 #endif
434  G4Exception("G4VUserPhysicsList::SetParticleCuts ",
435  "Run0254", FatalException,
436  "No Default Region");
437  return;
438  }
439  region = (*(G4RegionStore::GetInstance()))[0];
440  }
441 
442  if ( !isSetDefaultCutValue ){
444  }
445 
446  G4ProductionCuts* pcuts = region->GetProductionCuts();
447  pcuts->SetProductionCut(cut,particleName);
448 #ifdef G4VERBOSE
449  if (verboseLevel>2){
450  G4cout << "G4VUserPhysicsList::SetParticleCuts: "
451  << " :" << cut/mm << "[mm]"
452  << " for "<< particleName << G4endl;
453  }
454 #endif
455 }
456 
459 {
460  //Prepare Physics table for all particles
462  while( (*theParticleIterator)() ){
464  PreparePhysicsTable(particle);
465  }
466 
467  // ask processes to prepare physics table
468  if (fRetrievePhysicsTable) {
470  // check if retrieve Cut Table successfully
471  if (!fIsRestoredCutValues) {
472 #ifdef G4VERBOSE
473  if (verboseLevel>0){
474  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
475  << " Retrieve Cut Table failed !!" << G4endl;
476  }
477 #endif
478  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
479  "Run0255", RunMustBeAborted,
480  "Fail to retrieve Production Cut Table");
481  } else {
482 #ifdef G4VERBOSE
483  if (verboseLevel>2){
484  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
485  << " Retrieve Cut Table successfully " << G4endl;
486  }
487 #endif
488  }
489  } else {
490 #ifdef G4VERBOSE
491  if (verboseLevel>2){
492  G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
493  << " does not retrieve Cut Table but calculate " << G4endl;
494  }
495 #endif
496  }
497 
498  // Sets a value to particle
499  // set cut values for gamma at first and for e- and e+
500  G4String particleName;
502  if(GammaP) BuildPhysicsTable(GammaP);
504  if(EMinusP) BuildPhysicsTable(EMinusP);
506  if(EPlusP) BuildPhysicsTable(EPlusP);
507  G4ParticleDefinition* ProtonP = theParticleTable->FindParticle("proton");
508  if(ProtonP) BuildPhysicsTable(ProtonP);
509 
510 
512  while( (*theParticleIterator)() ){
514  if( particle!=GammaP &&
515  particle!=EMinusP &&
516  particle!=EPlusP &&
517  particle!=ProtonP ){
518  BuildPhysicsTable(particle);
519  }
520  }
521 
522  // Set flag
523  fIsPhysicsTableBuilt = true;
524 
525 }
528 {
529  if (fRetrievePhysicsTable) {
530  if ( !fIsRestoredCutValues){
531  // fail to retreive cut tables
532 #ifdef G4VERBOSE
533  if (verboseLevel>0){
534  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
535  << "Physics table can not be retreived and will be calculated "
536  << G4endl;
537  }
538 #endif
539  fRetrievePhysicsTable = false;
540 
541  } else {
542 #ifdef G4VERBOSE
543  if (verboseLevel>2){
544  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
545  << " Retrieve Physics Table for "
546  << particle->GetParticleName() << G4endl;
547  }
548 #endif
549  // Retrieve PhysicsTable from files for proccesses
551  }
552  }
553 
554 #ifdef G4VERBOSE
555  if (verboseLevel>2){
556  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
557  << "Calculate Physics Table for "
558  << particle->GetParticleName() << G4endl;
559  }
560 #endif
561  // Rebuild the physics tables for every process for this particle type
562  // if particle is not ShortLived
563  if(!particle->IsShortLived()) {
564  G4ProcessManager* pManager = particle->GetProcessManager();
565  if (!pManager) {
566 #ifdef G4VERBOSE
567  if (verboseLevel>0){
568  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
569  <<" : No Process Manager for "
570  << particle->GetParticleName() << G4endl;
571  G4cout << particle->GetParticleName()
572  << " should be created in your PhysicsList" <<G4endl;
573  }
574 #endif
575  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
576  "Run0271", FatalException,
577  "No process manager");
578  return;
579  }
580  G4ProcessVector* pVector = pManager->GetProcessList();
581  if (!pVector) {
582 #ifdef G4VERBOSE
583  if (verboseLevel>0){
584  G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
585  <<" : No Process Vector for "
586  << particle->GetParticleName() <<G4endl;
587  }
588 #endif
589  G4Exception("G4VUserPhysicsList::BuildPhysicsTable",
590  "Run0272", FatalException,
591  "No process Vector");
592  return;
593  }
594  for (G4int j=0; j < pVector->size(); ++j) {
595  (*pVector)[j]->BuildPhysicsTable(*particle);
596  }
597  }
598 }
599 
602 {
603  // Prepare the physics tables for every process for this particle type
604  // if particle is not ShortLived
605  if(!particle->IsShortLived()) {
606  G4ProcessManager* pManager = particle->GetProcessManager();
607  if (!pManager) {
608 #ifdef G4VERBOSE
609  if (verboseLevel>0) {
610  G4cout<< "G4VUserPhysicsList::PreparePhysicsTable "
611  << ": No Process Manager for "
612  << particle->GetParticleName() <<G4endl;
613  G4cout << particle->GetParticleName()
614  << " should be created in your PhysicsList" <<G4endl;
615  }
616 #endif
617  G4Exception("G4VUserPhysicsList::PreparePhysicsTable",
618  "Run0273", FatalException,
619  "No process manager");
620  return;
621  }
622 
623  G4ProcessVector* pVector = pManager->GetProcessList();
624  if (!pVector) {
625 #ifdef G4VERBOSE
626  if (verboseLevel>0) {
627  G4cout << "G4VUserPhysicsList::PreparePhysicsTable "
628  << ": No Process Vector for "
629  << particle->GetParticleName() <<G4endl;
630  }
631 #endif
632  G4Exception("G4VUserPhysicsList::PreparePhysicsTable",
633  "Run0274", FatalException,
634  "No process Vector");
635  return;
636  }
637  for (G4int j=0; j < pVector->size(); ++j) {
638  (*pVector)[j]->PreparePhysicsTable(*particle);
639  }
640  }
641 }
642 
645  G4ParticleDefinition* particle)
646 {
647  //*********************************************************************
648  // temporary addition to make the integral schema of electromagnetic
649  // processes work.
650  //
651 
652  if ( (process->GetProcessName() == "Imsc") ||
653  (process->GetProcessName() == "IeIoni") ||
654  (process->GetProcessName() == "IeBrems") ||
655  (process->GetProcessName() == "Iannihil") ||
656  (process->GetProcessName() == "IhIoni") ||
657  (process->GetProcessName() == "IMuIoni") ||
658  (process->GetProcessName() == "IMuBrems") ||
659  (process->GetProcessName() == "IMuPairProd") ) {
660 #ifdef G4VERBOSE
661  if (verboseLevel>2){
662  G4cout << "G4VUserPhysicsList::BuildIntegralPhysicsTable "
663  << " BuildPhysicsTable is invoked for "
664  << process->GetProcessName()
665  << "(" << particle->GetParticleName() << ")" << G4endl;
666  }
667 #endif
668  process->BuildPhysicsTable(*particle);
669  }
670 }
671 
674 {
676  G4int idx = 0;
677  while( (*theParticleIterator)() ){
679  G4cout << particle->GetParticleName();
680  if ((idx++ % 4) == 3) {
681  G4cout << G4endl;
682  } else {
683  G4cout << ", ";
684  }
685  }
686  G4cout << G4endl;
687 }
688 
689 
692 {
693  fDisplayThreshold = flag;
694 }
695 
698 {
699  if(fDisplayThreshold==0) return;
701  fDisplayThreshold = 0;
702 }
703 
704 
707 {
708  G4bool ascii = fStoredInAscii;
709  G4String dir = directory;
710  if (dir.isNull()) dir = directoryPhysicsTable;
711  else directoryPhysicsTable = dir;
712 
713  // store CutsTable info
714  if (!fCutsTable->StoreCutsTable(dir, ascii)) {
715  G4Exception("G4VUserPhysicsList::StorePhysicsTable",
716  "Run0281", JustWarning,
717  "Fail to store Cut Table");
718  return false;
719  }
720 #ifdef G4VERBOSE
721  if (verboseLevel>2){
722  G4cout << "G4VUserPhysicsList::StorePhysicsTable "
723  << " Store material and cut values successfully" << G4endl;
724  }
725 #endif
726 
727  G4bool success= true;
728 
729  // loop over all particles in G4ParticleTable
731  while( (*theParticleIterator)() ){
733  // Store physics tables for every process for this particle type
734  G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
735  G4int j;
736  for ( j=0; j < pVector->size(); ++j) {
737  if (!(*pVector)[j]->StorePhysicsTable(particle,dir,ascii)){
738  G4String comment = "Fail to store physics table for ";
739  comment += (*pVector)[j]->GetProcessName();
740  comment += "(" + particle->GetParticleName() + ")";
741  G4Exception("G4VUserPhysicsList::StorePhysicsTable",
742  "Run0282", JustWarning,
743  comment);
744  success = false;
745  }
746  }
747  // end loop over processes
748  }
749  // end loop over particles
750  return success;
751 }
752 
753 
754 
757 {
758  fRetrievePhysicsTable = true;
759  if(!directory.isNull()) {
760  directoryPhysicsTable = directory;
761  }
763  fIsRestoredCutValues = false;
764 }
765 
768  const G4String& directory,
769  G4bool ascii)
770 {
771  G4int j;
772  G4bool success[100];
773  // Retrieve physics tables for every process for this particle type
774  G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
775  for ( j=0; j < pVector->size(); ++j) {
776  success[j] =
777  (*pVector)[j]->RetrievePhysicsTable(particle,directory,ascii);
778 
779  if (!success[j]) {
780 #ifdef G4VERBOSE
781  if (verboseLevel>2){
782  G4cout << "G4VUserPhysicsList::RetrievePhysicsTable "
783  << " Fail to retrieve Physics Table for "
784  << (*pVector)[j]->GetProcessName() << G4endl;
785  G4cout << "Calculate Physics Table for "
786  << particle->GetParticleName() << G4endl;
787  }
788 #endif
789  (*pVector)[j]->BuildPhysicsTable(*particle);
790  }
791  }
792  for ( j=0; j < pVector->size(); ++j) {
793  // temporary addition to make the integral schema
794  if (!success[j]) BuildIntegralPhysicsTable((*pVector)[j], particle);
795  }
796 }
797 
798 
801 {
802 #ifdef G4VERBOSE
803  if (verboseLevel>2){
804  G4cout << "G4VUserPhysicsList::SetApplyCuts for " << name << G4endl;
805  }
806 #endif
807  if(name=="all") {
812  } else {
814  }
815 }
816 
819 {
821 }
822 
823 
826 {
829  }
830 }
831 
834 {
836 }
837 
840 {
842 }
843 
846  G4ParticleDefinition* particle)
847 {
848  return thePLHelper->RegisterProcess(process, particle);
849 }
850 
853 {
855  // set verboseLevel for G4ProductionCutsTable same as one for G4VUserPhysicsList:
857 
859 
860 #ifdef G4VERBOSE
861  if (verboseLevel >1){
862  G4cout << "G4VUserPhysicsList::SetVerboseLevel :"
863  << " Verbose level is set to " << verboseLevel << G4endl;
864  }
865 #endif
866 }
867 
868 
871 
874 {
875 #ifdef G4VERBOSE
876  if (verboseLevel>0){
877  G4cout << "G4VUserPhysicsList::ResetCuts() is obsolete."
878  << " This method gives no effect and you can remove it. "<< G4endl;
879  }
880 #endif
881 }