CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MonopoleSteppingAction.cc
Go to the documentation of this file.
2 
6 #include "SimG4Core/Physics/interface/G4Monopole.hh"
7 
11 
14 
15 #include "G4Track.hh"
16 #include "G4Run.hh"
17 #include "G4Event.hh"
18 #include "G4ParticleTable.hh"
19 
20 #include "CLHEP/Units/GlobalSystemOfUnits.h"
21 #include "CLHEP/Units/GlobalPhysicalConstants.h"
22 
24  actOnTrack (false), bZ(0) {
25  mode = p.getUntrackedParameter<bool>("ChangeFromFirstStep",true);
26  edm::LogInfo("SimG4CoreWatcher") << "MonopoleSeppingAction set mode for"
27  << " start at first step to " << mode;
28 }
29 
31 
33 
34  const edm::EventSetup* iSetup = (*job)();
36  iSetup->get<IdealMagneticFieldRecord>().get(bFieldH);
37  const MagneticField *bField = bFieldH.product();
38  const GlobalPoint p(0,0,0);
39  bZ = (bField->inTesla(p)).z();
40  edm::LogInfo("SimG4CoreWatcher") << "Magnetic Field (X): "
41  << (bField->inTesla(p)).x() << " Y: "
42  << (bField->inTesla(p)).y() << " Z: "
43  << bZ;
44 }
45 
47 
48  G4ParticleTable * partTable = G4ParticleTable::GetParticleTable();
49  magCharge = CLHEP::e_SI/CLHEP::fine_structure_const * 0.5;
50  for (int ii= 0; ii < partTable->size(); ii++) {
51  G4ParticleDefinition * particle = partTable->GetParticle(ii);
52  std::string particleName = (particle->GetParticleName()).substr(0,8);
53  if (strcmp(particleName.c_str(),"Monopole") == 0) {
54  magCharge = CLHEP::e_SI*((G4Monopole*)(particle))->MagneticCharge();
55  pdgCode.push_back(particle->GetPDGEncoding());
56  }
57  }
58  edm::LogInfo("SimG4CoreWatcher") << "MonopoleSeppingAction Finds "
59  << pdgCode.size() << " candidates";
60  for (unsigned int ii=0; ii<pdgCode.size(); ++ii) {
61  edm::LogInfo("SimG4CoreWatcher") << "PDG Code[" << ii << "] = "
62  << pdgCode[ii];
63  }
64  cMevToJ = CLHEP::e_SI/CLHEP::eV;
65  cMeVToKgMByS = CLHEP::e_SI*CLHEP::meter/(CLHEP::eV*CLHEP::c_light*CLHEP::second);
66  cInMByS = CLHEP::c_light*CLHEP::second/CLHEP::meter;
67  edm::LogInfo("SimG4CoreWatcher") << "MonopoleSeppingAction Constants:"
68  << " MevToJoules = " << cMevToJ
69  << ", MevToKgm/s = " << cMeVToKgMByS
70  << ", c in m/s = " << cInMByS
71  << ", mag. charge = " << magCharge;
72 }
73 
75 
76  actOnTrack = false;
77  if (pdgCode.size() > 0) {
78  const G4Track * aTrack = (*trk)();
79  int code = aTrack->GetDefinition()->GetPDGEncoding();
80  if (std::count(pdgCode.begin(),pdgCode.end(),code) > 0) {
81  actOnTrack = true;
82  eStart = aTrack->GetTotalEnergy();
83  pxStart = aTrack->GetMomentum().x();
84  pyStart = aTrack->GetMomentum().y();
85  pzStart = aTrack->GetMomentum().z();
86  dirxStart = aTrack->GetMomentumDirection().x();
87  diryStart = aTrack->GetMomentumDirection().y();
88  dirzStart = aTrack->GetMomentumDirection().z();
89  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction Track "
90  << code << " Flag " << actOnTrack
91  << " (px,py,pz,E) = (" << pxStart/GeV
92  << ", " << pyStart/GeV << ", "
93  <<pzStart/GeV <<", " <<eStart/GeV <<")";
94  }
95  }
96 }
97 
98 void MonopoleSteppingAction::update(const G4Step* aStep) {
99 
100  if (actOnTrack) {
101  double eT, pT, pZ, tStep;
102  G4Track* aTrack = aStep->GetTrack();
103  G4ThreeVector initialPosition(0,0,0);
104  if (mode) {
105  tStep = aTrack->GetGlobalTime();
108  pZ = pzStart;
109  } else {
110  G4ThreeVector dirStep = aTrack->GetMomentumDirection();
111  double lStep = aTrack->GetStepLength();
112  double xStep = aTrack->GetPosition().x()-lStep*dirStep.x();
113  double yStep = aTrack->GetPosition().y()-lStep*dirStep.y();
114  double zStep = aTrack->GetPosition().z()-lStep*dirStep.z();
115  double vStep = aTrack->GetVelocity();
116  initialPosition = G4ThreeVector(xStep,yStep,zStep);
117  tStep = (vStep > 0. ? (lStep/vStep) : 0.);
118  double eStep = aTrack->GetTotalEnergy();
119  eT = eStep*dirStep.perp();
120  pT = aTrack->GetMomentum().perp();
121  pZ = aTrack->GetMomentum().z();
122  }
123  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: tStep " <<tStep
124  << " eT " << eT << " pT " << pT << " pZ "
125  << pZ;
126  double eT1 = eT*cMevToJ;
127  double pT1 = pT*cMeVToKgMByS;
128  double pZ1 = pZ*cMeVToKgMByS;
129  double ts1 = tStep/CLHEP::second;
130  double eT2 = (eT1 > 0. ? (1./eT1) : 0.);
131  double fac0 = magCharge*bZ*cInMByS;
132  double fac2 = pZ1*cInMByS*eT2;
133  double fac1 = fac0*cInMByS*ts1*eT2 + fac2;
134  double fact1 = (eT1/fac0)*(std::sqrt(1+fac1*fac1)-std::sqrt(1+fac2*fac2));
135  double fact2 = (pT1/(magCharge*bZ)*(asinh(fac1)-asinh(fac2)));
136  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: Factor eT = "
137  << eT << " " << eT1 << " " << eT2
138  << ", pT = " << pT << " " << pT1
139  << ", pZ = " << pZ << " " << pZ1
140  << ", time = " << tStep << " " << ts1
141  << ", Factors " << fac0 << " " << fac1
142  << " " << fac2 << ", Final ones "
143  << fact1 << " " << fact2;
144  G4ThreeVector ez(0,0,1), et(std::sqrt(0.5),std::sqrt(0.5),0);
145  G4ThreeVector displacement = (fact1*ez + fact2*et)*CLHEP::m;
146 
147  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: Initial "
148  << initialPosition << " Displacement "
149  << displacement << " Final "
150  << (initialPosition+displacement);
151  aTrack->SetPosition(initialPosition+displacement);
152  }
153 }
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
int ii
Definition: cuy.py:588
float float float z
U second(std::pair< T, U > const &p)
#define e_SI
T sqrt(T t)
Definition: SSEVec.h:48
MonopoleSteppingAction(edm::ParameterSet const &p)
void update(const BeginOfJob *)
This routine will be called when the appropriate signal arrives.
const T & get() const
Definition: EventSetup.h:55
Definition: DDAxes.h:10