CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CMSMonopolePhysics.cc
Go to the documentation of this file.
3 
4 #include "G4ParticleDefinition.hh"
5 #include "G4ProcessManager.hh"
6 
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"
15 
16 #include "CLHEP/Units/GlobalSystemOfUnits.h"
17 
19  sim::FieldBuilder * fB_,
20  const edm::ParameterSet & p) :
21  G4VPhysicsConstructor("Monopole Physics"), fieldBuilder(fB_) {
22 
23  verbose = p.getUntrackedParameter<int>("Verbosity",0);
24  magCharge = p.getUntrackedParameter<int>("MonopoleCharge",1);
25  deltaRay = p.getUntrackedParameter<bool>("MonopoleDeltaRay",true);
26  multiSc = p.getUntrackedParameter<bool>("MonopoleMultiScatter",false);
27  transport = p.getUntrackedParameter<bool>("MonopoleTransport",true);
28  double mass = p.getUntrackedParameter<double>("MonopoleMass",200);
29  if (pdt && mass > 0.0) {
30  int ii=0;
31  for (HepPDT::ParticleDataTable::const_iterator p=pdt->begin();
32  p != pdt->end(); ++p,++ii) {
33  HepPDT::ParticleData particle = (p->second);
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()));
39  pdgEncodings.push_back(particle.pid());
40  monopoles.push_back(0);
41  //std::cout << "CMSMonopolePhysics: Monopole[" << ii
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()));
51  pdgEncodings.push_back(particle.pid());
52  monopoles.push_back(0);
53  //std::cout << "CMSMonopolePhysics: Monopole[" << ii
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;
59  }
60  }
61  }
62  // std::cout << "CMSMonopolePhysics has " << names.size()
63  if (verbose > 0) G4cout << "CMSMonopolePhysics has " << names.size()
64  << " monopole candidates and delta Ray option "
65  << deltaRay << G4endl;
66 }
67 
69 
71 
72  for (unsigned int ii=0; ii<names.size(); ++ii) {
73  if (!monopoles[ii]) {
74  G4Monopole* mpl = new G4Monopole(names[ii], pdgEncodings[ii],
75  masses[ii], ((pdgEncodings[ii] > 0 ) ? magCharge : -magCharge), elCharges[ii]);;
76  monopoles[ii] = mpl;
77  //std::cout << "Create G4Monopole " << names[ii]
78  if (verbose > 0) G4cout << "Create G4Monopole " << names[ii]
79  << " of mass " << masses[ii]/CLHEP::GeV
80  << " GeV, magnetic charge " << ((pdgEncodings[ii] > 0) ? magCharge : -magCharge)
81  << ", electric charge " << elCharges[ii]
82  << " and PDG encoding " << pdgEncodings[ii]
83  << " at " << monopoles[ii] << G4endl;
84  }
85  }
86 }
87 
89  // Add standard EM Processes
90 
91  if (verbose > 0)
92  G4cout << "### CMSMonopolePhysics ConstructProcess()" << G4endl;
93  for (unsigned int ii=0; ii<monopoles.size(); ++ii) {
94  if (monopoles[ii]) {
95  G4Monopole* mpl = monopoles[ii];
96  G4ProcessManager* pmanager = new G4ProcessManager(mpl);
97  mpl->SetProcessManager(pmanager);
98 
99  G4String particleName = mpl->GetParticleName();
100  G4double magn = mpl->MagneticCharge();
101  G4double mass = mpl->GetPDGMass();
102  if (verbose > 1)
103  G4cout << "### CMSMonopolePhysics instantiates for " << particleName
104  << " at " << mpl << " Mass " << mass/CLHEP::GeV
105  << " GeV Mag " << magn << " Process manager " << pmanager
106  << G4endl;
107 
108  // defined monopole parameters and binning
109 
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;
115  nbin *= 10;
116 
117  G4int idx = 1;
118  if (verbose > 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;
122 
123  if (magn == 0.0 || (!transport)) {
124  pmanager->AddProcess( new G4Transportation(verbose), -1, 0, 0);
125  } else {
126  pmanager->AddProcess( new G4MonopoleTransportation(mpl,fieldBuilder,verbose), -1, 0, 0);
127  }
128 
129  if (mpl->GetPDGCharge() != 0.0) {
130  if (multiSc) {
131  G4hMultipleScattering* hmsc = new G4hMultipleScattering();
132  pmanager->AddProcess(hmsc, -1, idx, idx);
133  ++idx;
134  }
135  if (deltaRay) {
136  G4hIonisation* hhioni = new G4hIonisation();
137  hhioni->SetDEDXBinning(nbin);
138  hhioni->SetMinKinEnergy(emin);
139  hhioni->SetMaxKinEnergy(emax);
140  pmanager->AddProcess(hhioni, -1, idx, idx);
141  } else {
142  G4hhIonisation* hhioni = new G4hhIonisation();
143  hhioni->SetDEDXBinning(nbin);
144  hhioni->SetMinKinEnergy(emin);
145  hhioni->SetMaxKinEnergy(emax);
146  pmanager->AddProcess(hhioni, -1, idx, idx);
147  }
148  ++idx;
149  }
150  if(magn != 0.0) {
151  G4mplIonisation* mplioni = new G4mplIonisation(magn);
152  mplioni->SetDEDXBinning(nbin);
153  mplioni->SetMinKinEnergy(emin);
154  mplioni->SetMaxKinEnergy(emax);
155  if (!deltaRay) {
156  G4mplIonisationWithDeltaModel* mod =
157  new G4mplIonisationWithDeltaModel(magn,"PAI");
158  mplioni->AddEmModel(0,mod,mod);
159  }
160  pmanager->AddProcess(mplioni, -1, idx, idx);
161  ++idx;
162  }
163  pmanager->AddProcess( new G4StepLimiter(), -1, -1, idx);
164  if (verbose > 1) pmanager->DumpInfo();
165  }
166  }
167 }
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
std::vector< int > pdgEncodings
HepPDT::ParticleDataTable ParticleDataTable
int ii
Definition: cuy.py:588
std::vector< G4Monopole * > monopoles
std::vector< int > elCharges
HepPDT::ParticleData ParticleData
std::vector< std::string > names
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< double > masses
sim::FieldBuilder * fieldBuilder
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
CMSMonopolePhysics(const HepPDT::ParticleDataTable *table_, sim::FieldBuilder *fB_, const edm::ParameterSet &p)