CMS 3D CMS Logo

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