Public Member Functions | |
BLCMDphysics () | |
Default Constructor. | |
~BLCMDphysics () | |
Destructor. | |
void | setDoStochastics (BLSetValue value, G4int warn=1) |
setDoStochastics() sets whether or not stochastic processes are to be enabled. | |
G4VUserPhysicsList * | getPhysicsList () |
getPhysicsList() returns the G4PhysicsList. | |
G4String | commandName () |
commandName() returns "physics". | |
int | command (BLArgumentVector &argv, BLArgumentMap &namedArgs) |
command() implements the command. | |
void | defineNamedArgs () |
defineNamedArgs() defines the named arguments. | |
void | help (bool detailed) |
help() provides special help text for the physics command. | |
void | callback (int type) |
callback() is from BLCallback, and adds a StepLimiter to the physics list for every particle. | |
bool | isLocalPhysicsList (G4String name) |
isLocalPhysicsList() returns true if thsi is a local physics list. | |
G4VUserPhysicsList * | getLocalPhysicsList (G4String name) |
getLocalPhysicsList() returns a local physics list. | |
void | listPhysicsLists () |
listPhysicsLists() will list all known physics lists, with synopsis. | |
void | printCase (G4String name, G4String description) |
printCase() prints the description of a physics list. | |
G4String | getSynopsis (G4String listName) |
getSynopsis() returns the synopsis for a physics list, by name. | |
Private Attributes | |
G4VUserPhysicsList * | physicsList |
G4int | doStochastics |
G4String | disable |
G4int | list |
G4int | synchrotronRadiation |
G4int | synchrotronRadiationMuon |
G4int | gammaToMuPair |
std::vector< G4String > | disableList |
G4PhysListFactory * | physListFactory |
BLCMDphysics::BLCMDphysics | ( | ) |
Default Constructor.
References BLCMDTYPE_PHYSICS, disable, doStochastics, gammaToMuPair, list, BLPhysics::minRangeCut, physicsList, physListFactory, BLCommand::registerCommand(), BLCommand::setDescription(), BLCommand::setSynopsis(), synchrotronRadiation, and synchrotronRadiationMuon.
00102 : BLPhysics(), BLCommand(), disableList() 00103 { 00104 registerCommand(BLCMDTYPE_PHYSICS); 00105 setSynopsis("Defines the physics processes and controls them."); 00106 setDescription("Exactly one physics command must be present.\n" 00107 "This command implements the geant4 physics lists\n" 00108 "of physics processes.\n" 00109 "The command is 'physics QGSP' for the QGSP set, and\n" 00110 "similarly for the other sets. With no argument it prints\n" 00111 "the available physics lists.\n" 00112 "Note that stochastic processes are always disabled while\n" 00113 "tracking the tune and reference particles.\n" 00114 "The only non-stochastic processes are Transportation and\n" 00115 "ionization energy loss (with fluctuations disabled).\n" 00116 "For muon beam studies you may want to disable the Decay " 00117 "process.\n\n" 00118 "NOTE: this command defines the particles used throughout " 00119 "the simulation, so this command must come before others " 00120 "that use particle names.\n\n" 00121 "NOTE: the rare decy mode for pi+/pi- to e nu is added.\n\n" 00122 "The default all-around physics list for HEP is called " 00123 "'QGSP_BERT'."); 00124 00125 physicsList = 0; 00126 doStochastics = 1; 00127 list = 0; 00128 synchrotronRadiation = 0; 00129 synchrotronRadiationMuon = 0; 00130 gammaToMuPair = 0; 00131 disable = ""; 00132 minRangeCut = 1.0 * mm; 00133 physListFactory = 0; 00134 }
void BLCMDphysics::setDoStochastics | ( | BLSetValue | value, | |
G4int | warn = 1 | |||
) | [virtual] |
setDoStochastics() sets whether or not stochastic processes are to be enabled.
Implements BLPhysics.
References disableList, doStochastics, FORCE_OFF, FORCE_ON, BLPhysics::isStochasticProcess(), list, BLPhysics::minRangeCut, NORMAL, and BLPhysics::stochasticsEnabled.
00224 { 00225 if(value == NORMAL) { 00226 if(doStochastics) 00227 value = FORCE_ON; 00228 else 00229 value = FORCE_OFF; 00230 } 00231 00232 G4UImanager* UI = G4UImanager::GetUIpointer(); 00233 G4String cmd; 00234 char rangeCmd[64]; 00235 switch(value) { 00236 case FORCE_OFF: 00237 // turn off delta rays (and other particle production) by 00238 // setting production cut very high in energy by setting 00239 // the required range large. 00240 UI->ApplyCommand("/run/setCut 1 km"); 00241 // Turn off enegry-loss straggling 00242 UI->ApplyCommand("/process/eLoss/fluct false"); 00243 cmd = "/process/inactivate "; 00244 switch(warn) { 00245 case 0: break; 00246 case 1: printf("*** NOTE: All stochastic processes are disabled!\n"); 00247 break; 00248 case 2: G4Exception("physics", 00249 "All stochastic processes disabled", 00250 JustWarning,""); 00251 break; 00252 } 00253 stochasticsEnabled = false; 00254 break; 00255 case FORCE_ON: 00256 // turn on delta rays (and other particle production) 00257 sprintf(rangeCmd,"/run/setCut %.6g mm",minRangeCut/mm); 00258 UI->ApplyCommand(rangeCmd); 00259 // Turn on enegry-loss straggling 00260 UI->ApplyCommand("/process/eLoss/fluct true"); 00261 cmd = "/process/activate "; 00262 if(warn != 0) printf("Stochastic processes are enabled.\n"); 00263 stochasticsEnabled = true; 00264 break; 00265 case NORMAL: 00266 return; 00267 } 00268 00269 // run through list of processes, applying cmd to stochastic ones 00270 G4ProcessTable *pt = G4ProcessTable::GetProcessTable(); 00271 G4ProcessTable::G4ProcNameVector *nv = pt->GetNameList(); 00272 if(list) printf("List of physics processes:\n"); 00273 for(unsigned int i=0; i<nv->size(); ++i) { 00274 G4String name = (*nv)[i]; 00275 if(list) printf("\t%s\n",name.c_str()); 00276 if(isStochasticProcess(name)) 00277 UI->ApplyCommand(cmd+name); 00278 } 00279 if(list) printf("\n"); 00280 list = 0; // print list only once 00281 00282 // now disable all proceses in the disableList 00283 cmd = "/process/inactivate "; 00284 for(unsigned i=0; i<disableList.size(); ++i) { 00285 switch(warn) { 00286 case 0: break; 00287 case 1: printf("*** NOTE: Disabling physics process '%s'\n", 00288 disableList[i].c_str()); 00289 break; 00290 case 2: G4Exception("physics","Physics process disabled", 00291 JustWarning,disableList[i].c_str()); 00292 break; 00293 } 00294 UI->ApplyCommand(cmd+disableList[i]); 00295 } 00296 }
G4VUserPhysicsList* BLCMDphysics::getPhysicsList | ( | ) | [inline, virtual] |
getPhysicsList() returns the G4PhysicsList.
Implements BLPhysics.
References physicsList.
00066 { return physicsList; }
G4String BLCMDphysics::commandName | ( | ) | [inline, virtual] |
int BLCMDphysics::command | ( | BLArgumentVector & | argv, | |
BLArgumentMap & | namedArgs | |||
) | [virtual] |
command() implements the command.
Implements BLCommand.
References disable, disableList, getLocalPhysicsList(), BLManager::getObject(), BLCommand::handleNamedArgs(), isLocalPhysicsList(), listPhysicsLists(), physicsList, physListFactory, BLCommand::print(), BLCommand::printError(), BLManager::registerCallback(), BLManager::registerPhysics(), and BLCommand::splitString().
00137 { 00138 if(physicsList != 0 && argv.size() > 0) { 00139 printError("physics: Multiple physics commands not allowed!"); 00140 return -1; 00141 } 00142 00143 if(argv.size() == 0) 00144 return -1; 00145 00146 // handle obsolete physics lists 00147 if(argv[0] == "QGSC") { 00148 argv[0] = "QGSC_BERT"; 00149 printf("Physics list QGSC is obsolete, replaced by QGSC_BERT\n"); 00150 } 00151 if(argv[0] == "QGSP") { 00152 argv[0] = "QGSP_BERT"; 00153 printf("Physics list QGSP is obsolete, replaced by QGSP_BERT\n"); 00154 } 00155 00156 if(!physListFactory) physListFactory = new G4PhysListFactory(); 00157 if(isLocalPhysicsList(argv[0])) 00158 physicsList = getLocalPhysicsList(argv[0]); 00159 if(!physicsList) { 00160 if(!physListFactory->IsReferencePhysList(argv[0])) { 00161 printError("physics: Unknown physics list '%s'", 00162 argv[0].c_str()); 00163 printError("The choices are:"); 00164 listPhysicsLists(); 00165 return -1; 00166 } 00167 physicsList = physListFactory->GetReferencePhysList(argv[0]); 00168 } 00169 00170 handleNamedArgs(namedArgs); 00171 00172 disableList = splitString(disable,',',true); 00173 00174 BLManager::getObject()->registerPhysics(this); 00175 00176 // now is too early to manage the physics processes; do it in a callback 00177 BLManager::getObject()->registerCallback(this,0); 00178 00179 print(argv[0]); 00180 00181 return 0; 00182 }
void BLCMDphysics::defineNamedArgs | ( | ) | [virtual] |
defineNamedArgs() defines the named arguments.
Reimplemented from BLCommand.
References BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), disable, doStochastics, gammaToMuPair, list, BLPhysics::minRangeCut, synchrotronRadiation, and synchrotronRadiationMuon.
00185 { 00186 argString(disable,"disable","Comma-separated list of processes to disable (e.g. 'Decay,msc')."); 00187 argString(disable,"inactivate","Synonym of disable."); 00188 argString(disable,"deactivate","Synonym of disable."); 00189 argInt(doStochastics,"doStochastics","Set to zero to disable all stochastic processes."); 00190 argDouble(minRangeCut,"minRangeCut","Minimum range cut for particle production (1 mm)"); 00191 argInt(list,"list","Nonzero to list the processes (later in output)."); 00192 argInt(gammaToMuPair,"gammaToMuPair","Nonzero to add gamma->mu+mu- (0)."); 00193 argInt(synchrotronRadiation,"synchrotronRadiation","Nonzero to add synchrotron radiation to the physics list for e- and e+."); 00194 argInt(synchrotronRadiationMuon,"synchrotronRadiationMuon","Nonzero to add synchrotron radiation to the physics list for mu- and mu+ NOTE: This is experimental, and may not work correctly."); 00195 }
void BLCMDphysics::help | ( | bool | detailed | ) | [virtual] |
help() provides special help text for the physics command.
Reimplemented from BLCommand.
References BLCommand::help(), and listPhysicsLists().
00198 { 00199 BLCommand::help(detailed); 00200 if(detailed) { 00201 printf("\n ----- PHYSICS LISTS -----\n\n"); 00202 printf(" Further guidance in selecting physics lists is available at:\n" 00203 "http://geant4.web.cern.ch/geant4/support/physicsLists/referencePL/index.shtml\n\n"); 00204 printf(" The default all-around physics list for HEP is called 'QGSP_BERT'.\n\n"); 00205 printf(" LHEP uses exclusively parameterized modeling.\n"); 00206 printf(" FTF lists use the FRITIOF description of string excitation\n" 00207 " and fragmentation.\n"); 00208 printf(" QGSP lists use the quark gluon string model\n"); 00209 printf(" QGSC are as QGSP except applying CHIPS modeling for the nuclear\n" 00210 " de-excitation.\n"); 00211 printf(" _BERT uses Geant4 Bertini cascade below ~ 10 GeV.\n"); 00212 printf(" _BIC uses Geant4 Binary cascade below ~ 10 GeV.\n"); 00213 printf(" _EMV suffix indicates a faster but less accurate EM modeling.\n"); 00214 printf(" _EMX suffix indicates the low-energy EM processes.\n"); 00215 printf(" _HP suffix uses the data driven high precision neutron package\n" 00216 " (thermal to 20 MeV).\n"); 00217 printf(" _NQE suffix indicates a list for comparison with earlier release.\n"); 00218 printf("\n List of available physics lists:\n"); 00219 listPhysicsLists(); 00220 } 00221 }
void BLCMDphysics::callback | ( | int | type | ) | [virtual] |
callback() is from BLCallback, and adds a StepLimiter to the physics list for every particle.
ensure only mu-nu is present
ensure only mu-nu is present
Reimplemented from BLCallback.
References BLAssert, gammaToMuPair, synchrotronRadiation, and synchrotronRadiationMuon.
00299 { 00300 G4ParticleTable::G4PTblDicIterator *theParticleIterator = 00301 G4ParticleTable::GetParticleTable()->GetIterator(); 00302 theParticleIterator -> reset(); 00303 while( (*theParticleIterator)() ) { 00304 G4ParticleDefinition* particle = theParticleIterator->value(); 00305 G4ProcessManager* pmanager = particle -> GetProcessManager(); 00306 if(!pmanager) { 00307 printf("BLCMDphysics: particle '%s' has no ProcessManager!\n", 00308 particle->GetParticleName().c_str()); 00309 continue; 00310 } 00311 pmanager -> AddProcess(new G4StepLimiter(),-1,-1,ordDefault); 00312 } 00313 00314 // add pi+/pi- decay to e nu 00315 double br=1.230e-4; 00316 G4ParticleTable *table = G4ParticleTable::GetParticleTable(); 00317 { G4ParticleDefinition *pi=table->FindParticle("pi+"); 00318 G4DecayTable* table = pi->GetDecayTable(); 00319 BLAssert(table->entries() == 1); /// ensure only mu-nu is present 00320 G4VDecayChannel *munu = table->GetDecayChannel(0); 00321 munu->SetBR(1.0-br); 00322 table->Insert(new G4PhaseSpaceDecayChannel("pi+",br,2,"e+","nu_e")); 00323 } 00324 { G4ParticleDefinition *pi=table->FindParticle("pi-"); 00325 G4DecayTable* table = pi->GetDecayTable(); 00326 BLAssert(table->entries() == 1); /// ensure only mu-nu is present 00327 G4VDecayChannel *munu = table->GetDecayChannel(0); 00328 munu->SetBR(1.0-br); 00329 table->Insert(new G4PhaseSpaceDecayChannel("pi-",br,2,"e-","anti_nu_e")); 00330 } 00331 00332 if(synchrotronRadiation) { 00333 table->FindParticle("e-")->GetProcessManager()-> 00334 AddDiscreteProcess(new G4SynchrotronRadiation); 00335 table->FindParticle("e+")->GetProcessManager()-> 00336 AddDiscreteProcess(new G4SynchrotronRadiation); 00337 } 00338 00339 if(synchrotronRadiationMuon) { 00340 table->FindParticle("mu-")->GetProcessManager()-> 00341 AddDiscreteProcess(new G4SynchrotronRadiation); 00342 table->FindParticle("mu+")->GetProcessManager()-> 00343 AddDiscreteProcess(new G4SynchrotronRadiation); 00344 } 00345 00346 if(gammaToMuPair) { 00347 table->FindParticle("gamma")->GetProcessManager()-> 00348 AddDiscreteProcess(new G4GammaConversionToMuons); 00349 } 00350 }
bool BLCMDphysics::isLocalPhysicsList | ( | G4String | name | ) |
G4VUserPhysicsList * BLCMDphysics::getLocalPhysicsList | ( | G4String | name | ) |
void BLCMDphysics::listPhysicsLists | ( | ) |
listPhysicsLists() will list all known physics lists, with synopsis.
References getSynopsis(), physListFactory, and printCase().
Referenced by command(), and help().
00364 { 00365 if(!physListFactory) physListFactory = new G4PhysListFactory(); 00366 const std::vector<G4String> v = physListFactory->AvailablePhysLists(); 00367 for(unsigned i=0; i<v.size(); ++i) { 00368 G4String n = v[i]; 00369 G4String s = getSynopsis(n); 00370 printCase(n,s); 00371 } 00372 }
void BLCMDphysics::printCase | ( | G4String | name, | |
G4String | description | |||
) |
printCase() prints the description of a physics list.
References BLCommand::wrapWords().
Referenced by listPhysicsLists().
00375 { 00376 G4String indent1(" "), indent(" "); 00377 indent += indent1; 00378 indent1 += name; 00379 do { indent1 += " "; } while(indent1.size() < indent.size()); 00380 printf("%s",wrapWords(description,indent1,indent).c_str()); 00381 }
G4String BLCMDphysics::getSynopsis | ( | G4String | listName | ) |
getSynopsis() returns the synopsis for a physics list, by name.
References CASE.
Referenced by listPhysicsLists().
00384 { 00385 #define CASE(NAME,SYNOPSIS) \ 00386 if(listName == #NAME) return SYNOPSIS; 00387 00388 CASE(CHIPS,"No synopsis available."); 00389 CASE(FTFP_BERT,"For calorimetry. The FTF model is based on the FRITIOF description of string excitation and fragmentation. Uses Geant4 Bertini cascade for primary protons, neutrons, pions and Kaons below ~10GeV."); 00390 CASE(FTFP_BERT_EMV,"Like FTFP_BERT but with faster EM modeling."); 00391 CASE(FTFP_BERT_EMX,"Like FTFP_BERT but with low-energy EM modeling."); 00392 CASE(FTFP_BERT_TRV,"A variant of QGSP_BERT where the Geant4 Bertini cascade is only used for particles below ~5.5 GeV."); 00393 CASE(FTFP_BIC,"For calorimetry. The FTF model is based on the FRITIOF description of string excitation and fragmentation. Uses Geant4 binary cascade for primary protons, neutrons, pions and Kaons below ~10GeV."); 00394 CASE(LBE,"For low background experiments (e.g. underground)"); 00395 CASE(LHEP,"For calorimetry -- is the fastest, when it comes to CPU. It uses the LEP and HEP parametrized models for inelastic scattering. The modeling parametrizes the final states individual inelastic reactions, so you will not see resonances, and the detailed secondary angular distributions for O(100MeV) reactions may not be described perfectly. The average quantities will be well described."); 00396 CASE(LHEP_BERT,"For 'typical' HEP collider detectors and low-energy penetration shielding (< 5 GeV) -- is as LHEP for energies above 3GeV. Below 3 GeV the Bertini cascade is used for nucleon and pion induced reactions. The price to pay is reduced CPU performance. It does not include gamma and electro nuclear physics."); 00397 CASE(LHEP_EMV,"Like LHEP but with faster EM modeling."); 00398 CASE(QBBC,"No synopsis available."); 00399 CASE(QGSC,"For calorimetry and high energy physics trackers -- is as QGSP for the initial reaction, but uses chiral invariant phase-space decay (multi-quasmon fragmentation) to model the behavior of the system's fragmentation."); 00400 CASE(QGSC_BERT,"For calorimetry and high energy physics trackers -- is as QGSP for the initial reaction, but uses chiral invariant phase-space decay (multi-quasmon fragmentation) to model the behavior of the system's fragmentation. Uses Geant4 Bertini cascade for nucleon and pion induced reactions."); 00401 CASE(QGSP,"For calorimetry and high energy physics trackers and high-energy and medium-energy production targets -- uses theory driven modeling for the reactions of energetic pions, kaons, and nucleons. It employs quark gluon string model for the 'punch-through' interactions of the projectile with a nucleus, the string excitation cross-sections being calculated in quasi-eikonal approximation. A pre-equilibrium decay model with an extensive evaporation phase to model the behavior of the nucleus 'after the punch'. It uses current best pion cross-section."); 00402 CASE(QGSP_BERT,"Like QGSP, but using Geant4 Bertini cascade for primary protons, neutrons, pions and Kaons below ~10GeV. In comparison to experimental data we find improved agreement to data compared to QGSP which uses the low energy parameterised (LEP) model for all particles at these energies. The Bertini model produces more secondary neutrons and protons than the LEP model, yielding a better agreement to experimental data."); 00403 CASE(QGSP_BERT_EMV,"Like QGSP_BERT but with faster EM modeling."); 00404 CASE(QGSP_BERT_EMX,"Like QGSP_BERT but with low-energy EM modeling."); 00405 CASE(QGSP_BERT_HP,"Like QGSP_BERT but with _HP modeling for neutrons."); 00406 CASE(QGSP_BIC,"Like QGSP, but using Geant4 Binary cascade for primary protons and neutrons with energies below ~10GeV, thus replacing the use of the LEP model for protons and neutrons In comparison to teh LEP model, Binary cascade better describes production of secondary particles produced in interactions of protons and neutrons with nuclei."); 00407 CASE(QGSP_BIC_EMX,"Like QGSP_BIC but with low-energy EM modeling."); 00408 CASE(QGSP_BIC_HP,"Like QGSP_BIC but with _HP modeling for neutrons."); 00409 00410 return "No synopsis available."; 00411 }
G4VUserPhysicsList* BLCMDphysics::physicsList [private] |
Referenced by BLCMDphysics(), command(), and getPhysicsList().
G4int BLCMDphysics::doStochastics [private] |
Referenced by BLCMDphysics(), defineNamedArgs(), and setDoStochastics().
G4String BLCMDphysics::disable [private] |
Referenced by BLCMDphysics(), command(), and defineNamedArgs().
G4int BLCMDphysics::list [private] |
Referenced by BLCMDphysics(), defineNamedArgs(), and setDoStochastics().
G4int BLCMDphysics::synchrotronRadiation [private] |
Referenced by BLCMDphysics(), callback(), and defineNamedArgs().
G4int BLCMDphysics::synchrotronRadiationMuon [private] |
Referenced by BLCMDphysics(), callback(), and defineNamedArgs().
G4int BLCMDphysics::gammaToMuPair [private] |
Referenced by BLCMDphysics(), callback(), and defineNamedArgs().
std::vector<G4String> BLCMDphysics::disableList [private] |
Referenced by command(), and setDoStochastics().
G4PhysListFactory* BLCMDphysics::physListFactory [private] |
Referenced by BLCMDphysics(), command(), and listPhysicsLists().