4 #include "G4ParticleDefinition.hh"
5 #include "G4ProcessManager.hh"
7 #include "G4StepLimiter.hh"
8 #include "G4Transportation.hh"
9 #include "G4MonopoleTransportation.hh"
10 #include "G4mplIonisation.hh"
11 #include "G4mplIonisationWithDeltaModel.hh"
12 #include "G4hMultipleScattering.hh"
13 #include "G4hIonisation.hh"
14 #include "G4hhIonisation.hh"
16 #include "CLHEP/Units/GlobalSystemOfUnits.h"
21 G4VPhysicsConstructor(
"Monopole Physics"), fieldBuilder(fB_) {
29 if (pdt && mass > 0.0) {
31 for (HepPDT::ParticleDataTable::const_iterator p=pdt->begin();
32 p != pdt->end(); ++
p,++
ii) {
34 std::string particleName = (particle.name()).substr(0,8);
35 if (strcmp(particleName.c_str(),
"Monopole") == 0) {
36 names.push_back(particle.name());
37 masses.push_back(mass*CLHEP::GeV);
38 elCharges.push_back((
int)(particle.charge()));
42 if (
verbose > 0) G4cout <<
"CMSMonopolePhysics: Monopole[" << ii
43 <<
"] " << particleName <<
" Mass "
44 << particle.mass() <<
" GeV, Magnetic Charge "
45 << magCharge <<
", Electric Charge "
46 << particle.charge() << G4endl;
47 }
else if(strcmp(particleName.c_str(),
"AntiMono") == 0) {
48 names.push_back(particle.name());
49 masses.push_back(mass*CLHEP::GeV);
50 elCharges.push_back((
int)(particle.charge()));
54 if (
verbose > 0) G4cout <<
"CMSMonopolePhysics: Monopole[" << ii
55 <<
"] " << particleName <<
" Mass "
56 << particle.mass() <<
" GeV, Magnetic Charge "
57 << magCharge <<
", Electric Charge "
58 << particle.charge() << G4endl;
63 if (
verbose > 0) G4cout <<
"CMSMonopolePhysics has " <<
names.size()
64 <<
" monopole candidates and delta Ray option "
65 << deltaRay << G4endl;
79 <<
" of mass " <<
masses[
ii]/CLHEP::GeV
92 G4cout <<
"### CMSMonopolePhysics ConstructProcess()" << G4endl;
96 G4ProcessManager* pmanager =
new G4ProcessManager(mpl);
97 mpl->SetProcessManager(pmanager);
99 G4String particleName = mpl->GetParticleName();
100 G4double magn = mpl->MagneticCharge();
101 G4double mass = mpl->GetPDGMass();
103 G4cout <<
"### CMSMonopolePhysics instantiates for " << particleName
104 <<
" at " << mpl <<
" Mass " << mass/CLHEP::GeV
105 <<
" GeV Mag " << magn <<
" Process manager " << pmanager
110 G4double emin = mass/20000.;
111 if(emin < CLHEP::keV) emin = CLHEP::keV;
112 G4double emax =
std::max(10.*CLHEP::TeV, mass*100);
113 G4int nbin = G4int(std::log10(emax/emin));
114 if (nbin < 1) nbin = 1;
119 G4cout <<
"### Magnetic charge " << magn <<
" and electric charge "
120 << mpl->GetPDGCharge() <<
"\n # of bins in dE/dx table = "
121 << nbin <<
" in the range " << emin <<
":" << emax << G4endl;
124 pmanager->AddProcess(
new G4Transportation(
verbose), -1, 0, 0);
129 if (mpl->GetPDGCharge() != 0.0) {
131 G4hMultipleScattering* hmsc =
new G4hMultipleScattering();
132 pmanager->AddProcess(hmsc, -1, idx, idx);
136 G4hIonisation* hhioni =
new G4hIonisation();
137 hhioni->SetDEDXBinning(nbin);
138 hhioni->SetMinKinEnergy(emin);
139 hhioni->SetMaxKinEnergy(emax);
140 pmanager->AddProcess(hhioni, -1, idx, idx);
142 G4hhIonisation* hhioni =
new G4hhIonisation();
143 hhioni->SetDEDXBinning(nbin);
144 hhioni->SetMinKinEnergy(emin);
145 hhioni->SetMaxKinEnergy(emax);
146 pmanager->AddProcess(hhioni, -1, idx, idx);
151 G4mplIonisation* mplioni =
new G4mplIonisation(magn);
152 mplioni->SetDEDXBinning(nbin);
153 mplioni->SetMinKinEnergy(emin);
154 mplioni->SetMaxKinEnergy(emax);
156 G4mplIonisationWithDeltaModel*
mod =
157 new G4mplIonisationWithDeltaModel(magn,
"PAI");
158 mplioni->AddEmModel(0,mod,mod);
160 pmanager->AddProcess(mplioni, -1, idx, idx);
163 pmanager->AddProcess(
new G4StepLimiter(), -1, -1, idx);
164 if (
verbose > 1) pmanager->DumpInfo();
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > pdgEncodings
HepPDT::ParticleDataTable ParticleDataTable
virtual ~CMSMonopolePhysics()
const T & max(const T &a, const T &b)
std::vector< G4Monopole * > monopoles
std::vector< int > elCharges
HepPDT::ParticleData ParticleData
std::vector< std::string > names
tuple idx
DEBUGGING if hasattr(process,"trackMonIterativeTracking2012"): print "trackMonIterativeTracking2012 D...
std::vector< double > masses
sim::FieldBuilder * fieldBuilder
T mod(const T &a, const T &b)
CMSMonopolePhysics(const HepPDT::ParticleDataTable *table_, sim::FieldBuilder *fB_, const edm::ParameterSet &p)