| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

DsPhysConsOptical.cc

Go to the documentation of this file.
00001 #define USE_CUSTOM_CERENKOV
00002 #define USE_CUSTOM_SCINTILLATION
00003 
00004 #include "DsPhysConsOptical.h"
00005 #include "DsG4OpRayleigh.h"
00006 
00007 #ifdef USE_CUSTOM_CERENKOV
00008 #include "DsG4Cerenkov.h"
00009 #else
00010 #include "G4Cerenkov.hh"
00011 #endif
00012 
00013 #ifdef USE_CUSTOM_SCINTILLATION
00014 #include "DsG4Scintillation.h"
00015 #else
00016 #include "G4Scintillation.hh"
00017 #endif
00018 
00019 #include "GaudiKernel/DeclareFactoryEntries.h" 
00020 #include "GaudiKernel/IProperty.h" 
00021 
00022 #include "G4OpAbsorption.hh"
00023 #include "G4OpRayleigh.hh"
00024 #include "G4OpBoundaryProcess.hh"
00025 #include "G4ProcessManager.hh"
00026 #include "G4FastSimulationManagerProcess.hh"
00027 
00028 DECLARE_TOOL_FACTORY(DsPhysConsOptical);
00029 
00030 
00031 DsPhysConsOptical::DsPhysConsOptical(const std::string& type,
00032                                      const std::string& name, 
00033                                      const IInterface* parent)
00034     : GiGaPhysConstructorBase(type,name,parent)
00035 {
00036     declareProperty("CerenMaxPhotonsPerStep",m_cerenMaxPhotonPerStep = 300,
00037                     "Limit step to at most this many (unscaled) Cerenkov photons.");
00038 
00039     declareProperty("ScintDoReemission",m_doReemission = true,
00040                     "Do reemission in scintilator.");
00041     declareProperty("ScintDoScintAndCeren",m_doScintAndCeren = true,
00042                     "Do both scintillation and Cerenkov in scintilator.");
00043 
00044     declareProperty("UseCerenkov", m_useCerenkov=true, 
00045                     "Use the Cerenkov process?");
00046 
00047     declareProperty("UseScintillation",m_useScintillation=true,
00048                     "Use the Scintillation process?");
00049     declareProperty("UseRayleigh", m_useRayleigh=true,
00050                     "Use the Rayleigh scattering process?");
00051     declareProperty("UseAbsorption", m_useAbsorption=true,
00052                     "Use light absorption process?");
00053     declareProperty("UseFastMu300nsTrick", m_useFastMu300nsTrick=false,
00054                     "Use Fast muon simulation?");
00055     declareProperty("ScintillationYieldFactor",m_ScintillationYieldFactor = 1.0,
00056                     "Scale the number of scintillation photons per MeV by this much.");
00057    
00058     declareProperty("BirksConstant1", m_birksConstant1 = 6.5e-3*g/cm2/MeV, 
00059                     "Birks constant C1");
00060     declareProperty("BirksConstant2", m_birksConstant2 = 2.1e-6*(g/cm2/MeV)*(g/cm2/MeV), 
00061                    "Birks constant C2");
00062 
00063     declareProperty("PreQE", m_preQE=0.32,
00064                     "Miximal possible PMT QE value to be applied on the stage of "
00065                     "generating photons in scintillation/Cerenkov processes");
00066 
00067     declareProperty("ApplyPreQE",m_applyPreQE=true,
00068                     "Apply preQE on stage of generating photons of not.");
00069 }
00070 
00071 StatusCode DsPhysConsOptical::initialize()
00072 {
00073     info()<<"Photons prescaling is "<<( m_applyPreQE?"on":"off" )
00074           <<". Preliminary applied efficiency is "<<m_preQE<<endreq;
00075     info()<<"Make sure the same values are set for the DsPmtSensDet"<<endreq;
00076     return this->GiGaPhysConstructorBase::initialize();
00077 }
00078 
00079 DsPhysConsOptical::~DsPhysConsOptical()
00080 {
00081 }
00082 
00083 void DsPhysConsOptical::ConstructParticle()
00084 {
00085 }
00086 
00087 void DsPhysConsOptical::ConstructProcess()
00088 {
00089 #ifdef USE_CUSTOM_CERENKOV
00090     
00091     info () << "Using customized DsG4Cerenkov." << endreq;
00092     DsG4Cerenkov* cerenkov = 0;
00093     if (m_useCerenkov) {
00094         cerenkov = new DsG4Cerenkov();
00095         cerenkov->SetMaxNumPhotonsPerStep(m_cerenMaxPhotonPerStep);
00096         cerenkov->SetApplyPreQE(m_applyPreQE);
00097         cerenkov->SetPreQE(m_preQE);
00098         // cerenkov->SetPhotonWeight(m_cerenPhotonScaleWeight);
00099         cerenkov->SetTrackSecondariesFirst(true);
00100     }
00101 #else
00102     info () << "Using standard G4Cerenkov." << endreq;
00103     G4Cerenkov* cerenkov = 0;
00104     if (m_useCerenkov) {
00105         cerenkov = new G4Cerenkov();
00106         cerenkov->SetMaxNumPhotonsPerStep(m_cerenMaxPhotonPerStep);
00107         cerenkov->SetTrackSecondariesFirst(true);
00108     }
00109 #endif
00110 
00111 #ifdef USE_CUSTOM_SCINTILLATION
00112     DsG4Scintillation* scint = 0;
00113     info() << "Using customized DsG4Scintillation." << endreq;
00114     scint = new DsG4Scintillation();
00115     scint->SetBirksConstant1(m_birksConstant1);
00116     scint->SetBirksConstant2(m_birksConstant2);
00117     scint->SetDoReemission(m_doReemission);
00118     scint->SetDoBothProcess(m_doScintAndCeren);
00119     // scint->SetPhotonWeight(m_scintPhotonScaleWeight);
00120     scint->SetApplyPreQE(m_applyPreQE);
00121     scint->SetPreQE(m_preQE);
00122     scint->SetScintillationYieldFactor(m_ScintillationYieldFactor); //1.);
00123     scint->SetUseFastMu300nsTrick(m_useFastMu300nsTrick);
00124     scint->SetTrackSecondariesFirst(true);
00125     if (!m_useScintillation) {
00126         scint->SetNoOp();
00127     }
00128 #else  // standard G4 scint
00129     G4Scintillation* scint = 0;
00130     if (m_useScintillation) {
00131         info() << "Using standard G4Scintillation." << endreq;
00132         scint = new G4Scintillation();
00133         scint->SetScintillationYieldFactor(m_ScintillationYieldFactor); // 1.);
00134         scint->SetTrackSecondariesFirst(true);
00135     }
00136 #endif
00137 
00138     G4OpAbsorption* absorb = 0;
00139     if (m_useAbsorption) {
00140         absorb = new G4OpAbsorption();
00141     }
00142 
00143     DsG4OpRayleigh* rayleigh = 0;
00144     if (m_useRayleigh) {
00145         rayleigh = new DsG4OpRayleigh();
00146         //        rayleigh->SetVerboseLevel(2);
00147     }
00148 
00149     G4OpBoundaryProcess* boundproc = new G4OpBoundaryProcess();
00150     boundproc->SetModel(unified);
00151 
00152     G4FastSimulationManagerProcess* fast_sim_man
00153         = new G4FastSimulationManagerProcess("fast_sim_man");
00154     
00155     theParticleIterator->reset();
00156     while( (*theParticleIterator)() ) {
00157 
00158         G4ParticleDefinition* particle = theParticleIterator->value();
00159         G4ProcessManager* pmanager = particle->GetProcessManager();
00160     
00161         // Caution: as of G4.9, Cerenkov becomes a Discrete Process.
00162         // This code assumes a version of G4Cerenkov from before this version.
00163 
00164         if(cerenkov && cerenkov->IsApplicable(*particle)) {
00165             pmanager->AddProcess(cerenkov);
00166             pmanager->SetProcessOrdering(cerenkov, idxPostStep);
00167             debug() << "Process: adding Cherenkov to " 
00168                     << particle->GetParticleName() << endreq;
00169         }
00170 
00171         if(scint && scint->IsApplicable(*particle)) {
00172             pmanager->AddProcess(scint);
00173             pmanager->SetProcessOrderingToLast(scint, idxAtRest);
00174             pmanager->SetProcessOrderingToLast(scint, idxPostStep);
00175             debug() << "Process: adding Scintillation to "
00176                     << particle->GetParticleName() << endreq;
00177         }
00178 
00179         if (particle == G4OpticalPhoton::Definition()) {
00180             if (absorb)
00181                 pmanager->AddDiscreteProcess(absorb);
00182             if (rayleigh)
00183                 pmanager->AddDiscreteProcess(rayleigh);
00184             pmanager->AddDiscreteProcess(boundproc);
00185             //pmanager->AddDiscreteProcess(pee);
00186             pmanager->AddDiscreteProcess(fast_sim_man);
00187         }
00188     }
00189 }
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Mon Apr 11 20:53:24 2011 for DetSim by doxygen 1.4.7