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

Last change on this file since 834 was 828, checked in by garnier, 17 years ago

import all except CVS

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