00001 #include "SimG4Core/PhysicsLists/interface/CMSMonopolePhysics.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003
00004 #include "G4ParticleDefinition.hh"
00005 #include "G4ProcessManager.hh"
00006
00007 #include "G4StepLimiter.hh"
00008 #include "G4Transportation.hh"
00009 #include "G4MonopoleTransportation.hh"
00010 #include "CMSG4mplIonisation.hh"
00011 #include "G4mplIonisationWithDeltaModel.hh"
00012 #include "G4hMultipleScattering.hh"
00013 #include "G4hIonisation.hh"
00014 #include "G4hhIonisation.hh"
00015
00016 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00017
00018 CMSMonopolePhysics::CMSMonopolePhysics(const HepPDT::ParticleDataTable * pdt,
00019 sim::FieldBuilder * fB_,
00020 const edm::ParameterSet & p) :
00021 G4VPhysicsConstructor("Monopole Physics"), fieldBuilder(fB_) {
00022
00023 verbose = p.getUntrackedParameter<int>("Verbosity",0);
00024 magCharge = p.getUntrackedParameter<int>("MonopoleCharge",1);
00025 deltaRay = p.getUntrackedParameter<bool>("MonopoleDeltaRay",true);
00026 multiSc = p.getUntrackedParameter<bool>("MonopoleMultiScatter",false);
00027 transport = p.getUntrackedParameter<bool>("MonopoleTransport",true);
00028 if (pdt) {
00029 int ii=0;
00030 for (HepPDT::ParticleDataTable::const_iterator p=pdt->begin();
00031 p != pdt->end(); ++p,++ii) {
00032 HepPDT::ParticleData particle = (p->second);
00033 std::string particleName = (particle.name()).substr(0,8);
00034 if (strcmp(particleName.c_str(),"Monopole") == 0) {
00035 names.push_back(particle.name());
00036 masses.push_back((particle.mass())*CLHEP::GeV);
00037 elCharges.push_back((int)(particle.charge()));
00038 pdgEncodings.push_back(particle.pid());
00039 monopoles.push_back(0);
00040 if (verbose > 0) G4cout << "CMSMonopolePhysics: Monopole[" << ii
00041 << "] " << particleName << " Mass "
00042 << particle.mass() << " GeV, Magnetic Charge "
00043 << magCharge << ", Electric Charge "
00044 << particle.charge() << G4endl;
00045 }
00046 }
00047 }
00048 if (verbose > 0) G4cout << "CMSMonopolePhysics has " << names.size()
00049 << " monopole candidates and delta Ray option "
00050 << deltaRay << G4endl;
00051 }
00052
00053 CMSMonopolePhysics::~CMSMonopolePhysics() {}
00054
00055 void CMSMonopolePhysics::ConstructParticle() {
00056
00057 for (unsigned int ii=0; ii<names.size(); ++ii) {
00058 if (!monopoles[ii]) {
00059 G4Monopole* mpl = new G4Monopole(names[ii], pdgEncodings[ii],
00060 masses[ii], ((pdgEncodings[ii] > 0 ) ? magCharge : -magCharge), elCharges[ii]);;
00061 monopoles[ii] = mpl;
00062 if (verbose > 0) G4cout << "Create G4Monopole " << names[ii]
00063 << " of mass " << masses[ii]/CLHEP::GeV
00064 << " GeV, magnetic charge " << ((pdgEncodings[ii] > 0) ? magCharge : -magCharge)
00065 << ", electric charge " << elCharges[ii]
00066 << " and PDG encoding " << pdgEncodings[ii]
00067 << " at " << monopoles[ii] << G4endl;
00068 }
00069 }
00070 }
00071
00072 void CMSMonopolePhysics::ConstructProcess() {
00073
00074
00075 if (verbose > 0)
00076 G4cout << "### CMSMonopolePhysics ConstructProcess()" << G4endl;
00077 for (unsigned int ii=0; ii<monopoles.size(); ++ii) {
00078 if (monopoles[ii]) {
00079 G4Monopole* mpl = monopoles[ii];
00080 G4ProcessManager* pmanager = new G4ProcessManager(mpl);
00081 mpl->SetProcessManager(pmanager);
00082
00083 G4String particleName = mpl->GetParticleName();
00084 G4double magn = mpl->MagneticCharge();
00085 if (verbose > 1)
00086 G4cout << "### CMSMonopolePhysics instantiates for " << particleName
00087 << " at " << mpl << " Mass " << mpl->GetPDGMass()/CLHEP::GeV
00088 << " GeV Mag " << magn << " Process manager " << pmanager
00089 << G4endl;
00090
00091
00092
00093 G4double emin = mpl->GetPDGMass()/20000.;
00094 if(emin < CLHEP::keV) emin = CLHEP::keV;
00095 G4double emax = 10.*CLHEP::TeV;
00096 G4int nbin = G4int(std::log10(emin/CLHEP::eV));
00097 emin = std::pow(10.,G4double(nbin))*CLHEP::eV;
00098
00099 nbin = G4int(std::log10(emax/emin));
00100 if (nbin < 1) nbin = 1;
00101 nbin *= 10;
00102
00103 G4int idx = 1;
00104 if (verbose > 1)
00105 G4cout << "### Magnetic charge " << magn << " and electric charge "
00106 << mpl->GetPDGCharge() <<"\n # of bins in dE/dx table = "
00107 << nbin << " in the range " << emin << ":" << emax << G4endl;
00108
00109 if (magn == 0.0 || (!transport)) {
00110 pmanager->AddProcess( new G4Transportation(verbose), -1, 0, 0);
00111 } else {
00112 pmanager->AddProcess( new G4MonopoleTransportation(mpl,fieldBuilder,verbose), -1, 0, 0);
00113 }
00114
00115 if (mpl->GetPDGCharge() != 0.0) {
00116 if (multiSc) {
00117 G4hMultipleScattering* hmsc = new G4hMultipleScattering();
00118 pmanager->AddProcess(hmsc, -1, idx, idx);
00119 ++idx;
00120 }
00121 if (deltaRay) {
00122 G4hIonisation* hhioni = new G4hIonisation();
00123 hhioni->SetDEDXBinning(nbin);
00124 hhioni->SetMinKinEnergy(emin);
00125 hhioni->SetMaxKinEnergy(emax);
00126 pmanager->AddProcess(hhioni, -1, idx, idx);
00127 } else {
00128 G4hhIonisation* hhioni = new G4hhIonisation();
00129 hhioni->SetDEDXBinning(nbin);
00130 hhioni->SetMinKinEnergy(emin);
00131 hhioni->SetMaxKinEnergy(emax);
00132 pmanager->AddProcess(hhioni, -1, idx, idx);
00133 }
00134 ++idx;
00135 }
00136 if(magn != 0.0) {
00137 CMSG4mplIonisation* mplioni = new CMSG4mplIonisation(magn);
00138 mplioni->SetDEDXBinning(nbin);
00139 mplioni->SetMinKinEnergy(emin);
00140 mplioni->SetMaxKinEnergy(emax);
00141 if (deltaRay) {
00142 G4mplIonisationWithDeltaModel* mod = new G4mplIonisationWithDeltaModel(magn,"PAI");
00143 mplioni->AddEmModel(0,mod,mod);
00144 }
00145 pmanager->AddProcess(mplioni, -1, idx, idx);
00146 ++idx;
00147 }
00148 pmanager->AddProcess( new G4StepLimiter(), -1, -1, idx);
00149 if (verbose > 1) pmanager->DumpInfo();
00150 }
00151 }
00152 }