CMS 3D CMS Logo

MonopoleSteppingAction.cc
Go to the documentation of this file.
2 
6 
10 
11 #include "G4Event.hh"
12 #include "G4ParticleTable.hh"
13 #include "G4Run.hh"
14 #include "G4Track.hh"
15 
16 #include "CLHEP/Units/GlobalPhysicalConstants.h"
17 #include "CLHEP/Units/GlobalSystemOfUnits.h"
18 
20  mode = p.getUntrackedParameter<bool>("ChangeFromFirstStep", true);
21  edm::LogVerbatim("SimG4CoreWatcher") << "MonopoleSeppingAction set mode for"
22  << " start at first step to " << mode;
23 }
24 
26 
29  edm::LogVerbatim("SimG4CoreWatcher") << "MonopoleSteppingAction::Initialize ESGetToken for MagneticField";
30 }
31 
34  const GlobalPoint p(0, 0, 0);
35  bZ = (bField->inTesla(p)).z();
36  edm::LogVerbatim("SimG4CoreWatcher") << "Magnetic Field (X): " << (bField->inTesla(p)).x()
37  << " Y: " << (bField->inTesla(p)).y() << " Z: " << bZ;
38 }
39 
41  G4ParticleTable *partTable = G4ParticleTable::GetParticleTable();
42  magCharge = CLHEP::e_SI / CLHEP::fine_structure_const * 0.5;
43  for (int ii = 0; ii < partTable->size(); ii++) {
44  G4ParticleDefinition *particle = partTable->GetParticle(ii);
45  std::string particleName = (particle->GetParticleName()).substr(0, 8);
46  if (strcmp(particleName.c_str(), "Monopole") == 0) {
47  magCharge = CLHEP::e_SI * ((Monopole *)(particle))->MagneticCharge();
48  pdgCode.push_back(particle->GetPDGEncoding());
49  }
50  }
51  edm::LogVerbatim("SimG4CoreWatcher") << "MonopoleSeppingAction Finds " << pdgCode.size() << " candidates";
52  for (unsigned int ii = 0; ii < pdgCode.size(); ++ii) {
53  edm::LogVerbatim("SimG4CoreWatcher") << "PDG Code[" << ii << "] = " << pdgCode[ii];
54  }
55  cMevToJ = CLHEP::e_SI / CLHEP::eV;
56  cMeVToKgMByS = CLHEP::e_SI * CLHEP::meter / (CLHEP::eV * CLHEP::c_light * CLHEP::second);
57  cInMByS = CLHEP::c_light * CLHEP::second / CLHEP::meter;
58  edm::LogVerbatim("SimG4CoreWatcher") << "MonopoleSeppingAction Constants:"
59  << " MevToJoules = " << cMevToJ << ", MevToKgm/s = " << cMeVToKgMByS
60  << ", c in m/s = " << cInMByS << ", mag. charge = " << magCharge;
61 }
62 
64  actOnTrack = false;
65  if (!pdgCode.empty()) {
66  const G4Track *aTrack = (*trk)();
67  int code = aTrack->GetDefinition()->GetPDGEncoding();
68  if (std::count(pdgCode.begin(), pdgCode.end(), code) > 0) {
69  actOnTrack = true;
70  eStart = aTrack->GetTotalEnergy();
71  pxStart = aTrack->GetMomentum().x();
72  pyStart = aTrack->GetMomentum().y();
73  pzStart = aTrack->GetMomentum().z();
74  dirxStart = aTrack->GetMomentumDirection().x();
75  diryStart = aTrack->GetMomentumDirection().y();
76  dirzStart = aTrack->GetMomentumDirection().z();
77  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction Track " << code << " Flag " << actOnTrack
78  << " (px,py,pz,E) = (" << pxStart / GeV << ", " << pyStart / GeV << ", "
79  << pzStart / GeV << ", " << eStart / GeV << ")";
80  }
81  }
82 }
83 
84 void MonopoleSteppingAction::update(const G4Step *aStep) {
85  if (actOnTrack) {
86  double eT, pT, pZ, tStep;
87  G4Track *aTrack = aStep->GetTrack();
88  G4ThreeVector initialPosition(0, 0, 0);
89  if (mode) {
90  tStep = aTrack->GetGlobalTime();
93  pZ = pzStart;
94  } else {
95  const G4ThreeVector &dirStep = aTrack->GetMomentumDirection();
96  double lStep = aTrack->GetStepLength();
97  double xStep = aTrack->GetPosition().x() - lStep * dirStep.x();
98  double yStep = aTrack->GetPosition().y() - lStep * dirStep.y();
99  double zStep = aTrack->GetPosition().z() - lStep * dirStep.z();
100  double vStep = aTrack->GetVelocity();
101  initialPosition = G4ThreeVector(xStep, yStep, zStep);
102  tStep = (vStep > 0. ? (lStep / vStep) : 0.);
103  double eStep = aTrack->GetTotalEnergy();
104  eT = eStep * dirStep.perp();
105  pT = aTrack->GetMomentum().perp();
106  pZ = aTrack->GetMomentum().z();
107  }
108  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: tStep " << tStep << " eT " << eT << " pT " << pT << " pZ "
109  << pZ;
110  double eT1 = eT * cMevToJ;
111  double pT1 = pT * cMeVToKgMByS;
112  double pZ1 = pZ * cMeVToKgMByS;
113  double ts1 = tStep / CLHEP::second;
114  double eT2 = (eT1 > 0. ? (1. / eT1) : 0.);
115  double fac0 = magCharge * bZ * cInMByS;
116  double fac2 = pZ1 * cInMByS * eT2;
117  double fac1 = fac0 * cInMByS * ts1 * eT2 + fac2;
118  double fact1 = (eT1 / fac0) * (std::sqrt(1 + fac1 * fac1) - std::sqrt(1 + fac2 * fac2));
119  double fact2 = (pT1 / (magCharge * bZ) * (asinh(fac1) - asinh(fac2)));
120  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: Factor eT = " << eT << " " << eT1 << " " << eT2
121  << ", pT = " << pT << " " << pT1 << ", pZ = " << pZ << " " << pZ1
122  << ", time = " << tStep << " " << ts1 << ", Factors " << fac0 << " " << fac1 << " "
123  << fac2 << ", Final ones " << fact1 << " " << fact2;
124  G4ThreeVector ez(0, 0, 1), et(std::sqrt(0.5), std::sqrt(0.5), 0);
125  G4ThreeVector displacement = (fact1 * ez + fact2 * et) * CLHEP::m;
126 
127  LogDebug("SimG4CoreWatcher") << "MonopoleSeppingAction: Initial " << initialPosition << " Displacement "
128  << displacement << " Final " << (initialPosition + displacement);
129  aTrack->SetPosition(initialPosition + displacement);
130  }
131 }
Log< level::Info, true > LogVerbatim
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > tok_bFieldH_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void registerConsumes(edm::ConsumesCollector) override
U second(std::pair< T, U > const &p)
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
def fac0(jeta)
marinas correction form
T sqrt(T t)
Definition: SSEVec.h:19
MonopoleSteppingAction(edm::ParameterSet const &p)
ii
Definition: cuy.py:589
void beginRun(edm::EventSetup const &) override
#define e_SI
#define LogDebug(id)