// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // /// \file OpNovice/src/OpNovicePhysicsList.cc /// \brief Implementation of the OpNovicePhysicsList class // // // //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... #include "globals.hh" #include "OpNovicePhysicsList.hh" #include "OpNovicePhysicsListMessenger.hh" #include "G4ParticleDefinition.hh" #include "G4ParticleTypes.hh" #include "G4ParticleTable.hh" #include "G4BosonConstructor.hh" #include "G4LeptonConstructor.hh" #include "G4MesonConstructor.hh" #include "G4BaryonConstructor.hh" #include "G4IonConstructor.hh" #include "G4ShortLivedConstructor.hh" #include "G4ProcessManager.hh" #include "G4Cerenkov.hh" #include "G4Scintillation.hh" #include "G4OpAbsorption.hh" #include "G4OpRayleigh.hh" #include "G4OpMieHG.hh" #include "G4OpBoundaryProcess.hh" #include "G4LossTableManager.hh" #include "G4EmSaturation.hh" G4ThreadLocal G4int OpNovicePhysicsList::fVerboseLevel = 1; G4ThreadLocal G4int OpNovicePhysicsList::fMaxNumPhotonStep = 20; G4ThreadLocal G4Cerenkov* OpNovicePhysicsList::fCerenkovProcess = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessDef = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessAlpha = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessProton = 0; //G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessNeutron = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessDeuteron = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessTriton = 0; G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcessNuc = 0; G4ThreadLocal G4OpAbsorption* OpNovicePhysicsList::fAbsorptionProcess = 0; G4ThreadLocal G4OpRayleigh* OpNovicePhysicsList::fRayleighScatteringProcess = 0; G4ThreadLocal G4OpMieHG* OpNovicePhysicsList::fMieHGScatteringProcess = 0; G4ThreadLocal G4OpBoundaryProcess* OpNovicePhysicsList::fBoundaryProcess = 0; //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... OpNovicePhysicsList::OpNovicePhysicsList() : G4VUserPhysicsList() { fMessenger = new OpNovicePhysicsListMessenger(this); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... OpNovicePhysicsList::~OpNovicePhysicsList() { delete fMessenger; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::ConstructParticle() { // In this method, static member functions should be called // for all particles which you want to use. // This ensures that objects of these particle types will be // created in the program. G4BosonConstructor bConstructor; bConstructor.ConstructParticle(); G4LeptonConstructor lConstructor; lConstructor.ConstructParticle(); G4MesonConstructor mConstructor; mConstructor.ConstructParticle(); G4BaryonConstructor rConstructor; rConstructor.ConstructParticle(); G4IonConstructor iConstructor; iConstructor.ConstructParticle(); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::ConstructProcess() { AddTransportation(); ConstructDecay(); ConstructEM(); ConstructOp(); ConstructHad(); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... #include "G4Decay.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::ConstructDecay() { // Add Decay Process G4Decay* theDecayProcess = new G4Decay(); auto particleIterator=GetParticleIterator(); particleIterator->reset(); while( (*particleIterator)() ){ G4ParticleDefinition* particle = particleIterator->value(); G4ProcessManager* pmanager = particle->GetProcessManager(); if (theDecayProcess->IsApplicable(*particle)) { pmanager ->AddProcess(theDecayProcess); // set ordering for PostStepDoIt and AtRestDoIt pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep); pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest); } } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... #include "G4ComptonScattering.hh" #include "G4GammaConversion.hh" #include "G4PhotoElectricEffect.hh" #include "G4eMultipleScattering.hh" #include "G4MuMultipleScattering.hh" #include "G4hMultipleScattering.hh" #include "G4eIonisation.hh" #include "G4eBremsstrahlung.hh" #include "G4eplusAnnihilation.hh" #include "G4MuIonisation.hh" #include "G4MuBremsstrahlung.hh" #include "G4MuPairProduction.hh" #include "G4hIonisation.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::ConstructEM() { auto particleIterator=GetParticleIterator(); particleIterator->reset(); while( (*particleIterator)() ){ G4ParticleDefinition* particle = particleIterator->value(); G4ProcessManager* pmanager = particle->GetProcessManager(); G4String particleName = particle->GetParticleName(); if (particleName == "gamma") { // gamma // Construct processes for gamma pmanager->AddDiscreteProcess(new G4GammaConversion()); pmanager->AddDiscreteProcess(new G4ComptonScattering()); pmanager->AddDiscreteProcess(new G4PhotoElectricEffect()); } else if (particleName == "e-") { //electron // Construct processes for electron pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1); pmanager->AddProcess(new G4eIonisation(), -1, 2, 2); pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3); } else if (particleName == "e+") { //positron // Construct processes for positron pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1); pmanager->AddProcess(new G4eIonisation(), -1, 2, 2); pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3); pmanager->AddProcess(new G4eplusAnnihilation(), 0,-1, 4); } else if( particleName == "mu+" || particleName == "mu-" ) { //muon // Construct processes for muon pmanager->AddProcess(new G4MuMultipleScattering(),-1, 1, 1); pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2); pmanager->AddProcess(new G4MuBremsstrahlung(), -1, 3, 3); pmanager->AddProcess(new G4MuPairProduction(), -1, 4, 4); } else { if ((particle->GetPDGCharge() != 0.0) && (particle->GetParticleName() != "chargedgeantino") && !particle->IsShortLived()) { // all others charged particles except geantino pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1); pmanager->AddProcess(new G4hIonisation(), -1,2,2); } } } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... #include "G4Threading.hh" void OpNovicePhysicsList::ConstructOp() { fCerenkovProcess = new G4Cerenkov("Cerenkov"); fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep); fCerenkovProcess->SetMaxBetaChangePerStep(10.0); fCerenkovProcess->SetTrackSecondariesFirst(true); //fScintillationProcess->SetScintillationByParticleType(false); // gammas (electrons) and neutrons (protons) the same // needs line in detector construction like: // LaCl3_mt->AddConstProperty("SCINTILLATIONYIELD",1./keV); //fScintillationProcess->SetScintillationByParticleType(true); // ionizing particles different // needs lines in detector construction like: //G4double energies[2] = {1.0*keV, 1.0*MeV}; //G4double yield[2] = {1., 1000.}; //LaCl3_mt->AddProperty("ELECTRONSCINTILLATIONYIELD", energies, yield, 2); //LaCl3_mt->AddProperty("CARBONSCINTILLATIONYIELD", energies, yield, 2); //LaCl3_mt->AddProperty("PROTONSCINTILLATIONYIELD", energies, yield, 2); // default scintillation process fScintillationProcessDef = new G4Scintillation("Scintillation"); fScintillationProcessDef->SetScintillationByParticleType(true); // ionizing particles different // fScintillationProcessDef->DumpPhysicsTable(); fScintillationProcessDef->SetTrackSecondariesFirst(true); fScintillationProcessDef->SetScintillationYieldFactor(1.0); // particle-dependent scale factor for SCINTILLATIONYIELD in G4MaterialPropertiesTable fScintillationProcessDef->SetScintillationExcitationRatio(0.80); // particle-dependent fast yield/total yield - overwrites the YIELDRATIO obtained from G4MaterialPropertiesTable of sintillator material // -- values for LaCl3 (e-, alpha from McFee NIMA 2016 Vol. 828, pp 105-115, proton is a guess) fAbsorptionProcess = new G4OpAbsorption(); fRayleighScatteringProcess = new G4OpRayleigh(); fMieHGScatteringProcess = new G4OpMieHG(); fBoundaryProcess = new G4OpBoundaryProcess(); // scintillation process for alpha: fScintillationProcessAlpha = new G4Scintillation("Scintillation"); fScintillationProcessAlpha->SetScintillationByParticleType(true); // fScintillationProcessAlpha->DumpPhysicsTable(); fScintillationProcessAlpha->SetTrackSecondariesFirst(true); fScintillationProcessAlpha->SetScintillationYieldFactor(1.0); fScintillationProcessAlpha->SetScintillationExcitationRatio(0.72); // for LaCl3 (e-, alpha from McFee NIMA 2016 Vol. 828, pp 105-115, proton is a guess) // scintillation process for proton: fScintillationProcessProton = new G4Scintillation("Scintillation"); fScintillationProcessProton->SetScintillationByParticleType(true); // fScintillationProcessProton->DumpPhysicsTable(); fScintillationProcessProton->SetTrackSecondariesFirst(true); fScintillationProcessProton->SetScintillationYieldFactor(1.0); fScintillationProcessProton->SetScintillationExcitationRatio(0.76); // for LaCl3 (e-, alpha from McFee NIMA 2016 Vol. 828, pp 105-115. proton is a guess) // scintillation process for neutron: //fScintillationProcessNeutron = new G4Scintillation("Scintillation"); //fScintillationProcessNeutron->SetScintillationByParticleType(true); // fScintillationProcessNeutron->DumpPhysicsTable(); //fScintillationProcessNeutron->SetTrackSecondariesFirst(true); //fScintillationProcessNeutron->SetScintillationYieldFactor(1.0); //fScintillationProcessNeutron->SetScintillationExcitationRatio(1.0); // scintillation process for deuteron: fScintillationProcessDeuteron = new G4Scintillation("Scintillation"); fScintillationProcessDeuteron->SetScintillationByParticleType(true); // fScintillationProcessDeuteron->DumpPhysicsTable(); fScintillationProcessDeuteron->SetTrackSecondariesFirst(true); fScintillationProcessDeuteron->SetScintillationYieldFactor(1.0); fScintillationProcessDeuteron->SetScintillationExcitationRatio(0.76); //using proton value // scintillation process for triton: fScintillationProcessTriton = new G4Scintillation("Scintillation"); fScintillationProcessTriton->SetScintillationByParticleType(true); // fScintillationProcessTriton->DumpPhysicsTable(); fScintillationProcessTriton->SetTrackSecondariesFirst(true); fScintillationProcessTriton->SetScintillationYieldFactor(1.0); fScintillationProcessTriton->SetScintillationExcitationRatio(0.76); //using proton value // scintillation process for heavy nuclei fScintillationProcessNuc = new G4Scintillation("Scintillation"); fScintillationProcessNuc->SetScintillationByParticleType(true); // fScintillationProcessNuc->DumpPhysicsTable(); fScintillationProcessNuc->SetTrackSecondariesFirst(true); fScintillationProcessNuc->SetScintillationYieldFactor(1.0); fScintillationProcessNuc->SetScintillationExcitationRatio(0.72); // using alpha value fCerenkovProcess->SetVerboseLevel(fVerboseLevel); fScintillationProcessDef->SetVerboseLevel(fVerboseLevel); fScintillationProcessAlpha->SetVerboseLevel(fVerboseLevel); fScintillationProcessProton->SetVerboseLevel(fVerboseLevel); //fScintillationProcessNeutron->SetVerboseLevel(fVerboseLevel); fScintillationProcessDeuteron->SetVerboseLevel(fVerboseLevel); fScintillationProcessTriton->SetVerboseLevel(fVerboseLevel); fScintillationProcessNuc->SetVerboseLevel(fVerboseLevel); fAbsorptionProcess->SetVerboseLevel(fVerboseLevel); fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel); fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel); fBoundaryProcess->SetVerboseLevel(fVerboseLevel); // Note:Birk's correction via fEmSaturation and specifying scintillation by particle type are physically mutually exclusive. So following is commented out. // Use Birks Correction in the Scintillation process //if(G4Threading::IsMasterThread()) //{ // G4EmSaturation* emSaturation = // G4LossTableManager::Instance()->EmSaturation(); // fScintillationProcessDef->AddSaturation(emSaturation); //} auto particleIterator=GetParticleIterator(); particleIterator->reset(); while( (*particleIterator)() ){ G4ParticleDefinition* particle = particleIterator->value(); G4ProcessManager* pmanager = particle->GetProcessManager(); G4String particleName = particle->GetParticleName(); if (fCerenkovProcess->IsApplicable(*particle)) { pmanager->AddProcess(fCerenkovProcess); pmanager->SetProcessOrdering(fCerenkovProcess,idxPostStep); } //if (fScintillationProcess->IsApplicable(*particle)) { // pmanager->AddProcess(fScintillationProcess); // pmanager->SetProcessOrderingToLast(fScintillationProcess, idxAtRest); // pmanager->SetProcessOrderingToLast(fScintillationProcess, idxPostStep); //} if (fScintillationProcessDef->IsApplicable(*particle)) { if(particle->GetParticleName() == "GenericIon") { pmanager->AddProcess(fScintillationProcessNuc); // AtRestDiscrete pmanager->SetProcessOrderingToLast(fScintillationProcessNuc,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessNuc,idxPostStep); } else if(particle->GetParticleName() == "alpha") { pmanager->AddProcess(fScintillationProcessAlpha); pmanager->SetProcessOrderingToLast(fScintillationProcessAlpha,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessAlpha,idxPostStep); } else if(particle->GetParticleName() == "proton") { pmanager->AddProcess(fScintillationProcessProton); pmanager->SetProcessOrderingToLast(fScintillationProcessProton,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessProton,idxPostStep); } //else if(particle->GetParticleName() == "neutron") // { // pmanager->AddProcess(fScintillationProcessNeutron); // pmanager->SetProcessOrderingToLast(fScintillationProcessNeutron,idxAtRest); // pmanager->SetProcessOrderingToLast(fScintillationProcessNeutron,idxPostStep); // } else if(particle->GetParticleName() == "deuteron") { pmanager->AddProcess(fScintillationProcessDeuteron); pmanager->SetProcessOrderingToLast(fScintillationProcessDeuteron,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessDeuteron,idxPostStep); } else if(particle->GetParticleName() == "triton") { pmanager->AddProcess(fScintillationProcessTriton); pmanager->SetProcessOrderingToLast(fScintillationProcessTriton,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessTriton,idxPostStep); } else { pmanager->AddProcess(fScintillationProcessDef); pmanager->SetProcessOrderingToLast(fScintillationProcessDef,idxAtRest); pmanager->SetProcessOrderingToLast(fScintillationProcessDef,idxPostStep); } } if (particleName == "opticalphoton") { G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl; pmanager->AddDiscreteProcess(fAbsorptionProcess); pmanager->AddDiscreteProcess(fRayleighScatteringProcess); pmanager->AddDiscreteProcess(fMieHGScatteringProcess); pmanager->AddDiscreteProcess(fBoundaryProcess); } } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... // Hadronic Processes // Elastic processes #include "G4HadronElasticProcess.hh" // Neutron high-precision models: <20 MeV #include "G4NeutronHPElastic.hh" #include "G4NeutronHPElasticData.hh" #include "G4NeutronHPCapture.hh" #include "G4NeutronHPCaptureData.hh" #include "G4NeutronHPInelastic.hh" #include "G4NeutronHPInelasticData.hh" // Capture processes #include "G4NeutronRadCapture.hh" //G4LCapture.hh is obsolete #include "G4HadronCaptureProcess.hh" // Inelastic processes #include "G4ProtonInelasticProcess.hh" #include "G4BinaryCascade.hh" #include "G4ProtonInelasticCrossSection.hh" #include "G4NeutronInelasticProcess.hh" #include "G4DeuteronInelasticProcess.hh" #include "G4QMDReaction.hh" #include "G4IonsShenCrossSection.hh" #include "G4TritonInelasticProcess.hh" #include "G4AlphaInelasticProcess.hh" #include "G4HadronElastic.hh" void OpNovicePhysicsList::ConstructHad() { G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess; G4HadronElastic* theElasticModel = new G4HadronElastic; theElasticProcess->RegisterMe(theElasticModel); auto theParticleIterator=GetParticleIterator(); theParticleIterator->reset(); while ((*theParticleIterator)()) { G4ParticleDefinition* particle = theParticleIterator->value(); G4ProcessManager* pmanager = particle->GetProcessManager(); G4String particleName = particle->GetParticleName(); if (particleName == "proton") { pmanager->AddDiscreteProcess(theElasticProcess); G4ProtonInelasticProcess* theInelasticProcess = new G4ProtonInelasticProcess("inelastic"); G4BinaryCascade * theProtonInelasticModel =new G4BinaryCascade; theInelasticProcess->RegisterMe(theProtonInelasticModel); G4ProtonInelasticCrossSection * theProtonData = new G4ProtonInelasticCrossSection; theInelasticProcess->AddDataSet(theProtonData); pmanager->AddDiscreteProcess(theInelasticProcess); } else if(particleName == "neutron") { // elastic scattering G4HadronElasticProcess* theNeutronElasticProcess = new G4HadronElasticProcess; G4HadronElastic* theElasticModel1 = new G4HadronElastic; G4NeutronHPElastic * theElasticNeutron = new G4NeutronHPElastic; theNeutronElasticProcess->RegisterMe(theElasticModel1); theElasticModel1->SetMinEnergy(19*MeV); theNeutronElasticProcess->RegisterMe(theElasticNeutron); G4NeutronHPElasticData * theNeutronData = new G4NeutronHPElasticData; theNeutronElasticProcess->AddDataSet(theNeutronData); pmanager->AddDiscreteProcess(theNeutronElasticProcess); // inelastic scattering G4NeutronInelasticProcess* theInelasticProcess = new G4NeutronInelasticProcess("inelastic"); G4NeutronHPInelastic * theLENeutronInelasticModel =new G4NeutronHPInelastic; theInelasticProcess->RegisterMe(theLENeutronInelasticModel); G4NeutronHPInelasticData * theNeutronData1 = new G4NeutronHPInelasticData; theInelasticProcess->AddDataSet(theNeutronData1); pmanager->AddDiscreteProcess(theInelasticProcess); // capture G4HadronCaptureProcess* theCaptureProcess =new G4HadronCaptureProcess; G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture; theCaptureModel->SetMinEnergy(19*MeV); theCaptureProcess->RegisterMe(theCaptureModel); G4NeutronHPCapture * theLENeutronCaptureModel = new G4NeutronHPCapture; theCaptureProcess->RegisterMe(theLENeutronCaptureModel); G4NeutronHPCaptureData * theNeutronData3 = new G4NeutronHPCaptureData; theCaptureProcess->AddDataSet(theNeutronData3); pmanager->AddDiscreteProcess(theCaptureProcess); } else if(particleName == "deuteron") { pmanager->AddDiscreteProcess(theElasticProcess); G4DeuteronInelasticProcess* theInelasticProcess = new G4DeuteronInelasticProcess("inelastic"); G4QMDReaction * theDeuteronInelasticModel =new G4QMDReaction; theInelasticProcess->RegisterMe(theDeuteronInelasticModel); G4IonsShenCrossSection * theDeuteronData = new G4IonsShenCrossSection; theInelasticProcess->AddDataSet(theDeuteronData); pmanager->AddDiscreteProcess(theInelasticProcess); } else if(particleName == "triton") { pmanager->AddDiscreteProcess(theElasticProcess); G4TritonInelasticProcess* theInelasticProcess = new G4TritonInelasticProcess("inelastic"); G4QMDReaction * theTritonInelasticModel =new G4QMDReaction; theInelasticProcess->RegisterMe(theTritonInelasticModel); G4IonsShenCrossSection * theTritonData = new G4IonsShenCrossSection; theInelasticProcess->AddDataSet(theTritonData); pmanager->AddDiscreteProcess(theInelasticProcess); } else if(particleName == "alpha") { pmanager->AddDiscreteProcess(theElasticProcess); G4AlphaInelasticProcess* theInelasticProcess = new G4AlphaInelasticProcess("inelastic"); G4QMDReaction * theAlphaInelasticModel =new G4QMDReaction; theInelasticProcess->RegisterMe(theAlphaInelasticModel); G4IonsShenCrossSection * theAlphaData = new G4IonsShenCrossSection; theInelasticProcess->AddDataSet(theAlphaData); pmanager->AddDiscreteProcess(theInelasticProcess); } } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::SetVerbose(G4int verbose) { fVerboseLevel = verbose; fCerenkovProcess->SetVerboseLevel(fVerboseLevel); fScintillationProcessDef->SetVerboseLevel(fVerboseLevel); fScintillationProcessAlpha->SetVerboseLevel(fVerboseLevel); fScintillationProcessProton->SetVerboseLevel(fVerboseLevel); //fScintillationProcessNeutron->SetVerboseLevel(fVerboseLevel); fScintillationProcessNuc->SetVerboseLevel(fVerboseLevel); fAbsorptionProcess->SetVerboseLevel(fVerboseLevel); fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel); fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel); fBoundaryProcess->SetVerboseLevel(fVerboseLevel); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::SetNbOfPhotonsCerenkov(G4int MaxNumber) { fMaxNumPhotonStep = MaxNumber; fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... void OpNovicePhysicsList::SetCuts() { // " G4VUserPhysicsList::SetCutsWithDefault" method sets // the default cut value for all particle types // SetCutsWithDefault(); if (verboseLevel>0) DumpCutValuesTable(); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......