CMS 3D CMS Logo

HFGflash.cc
Go to the documentation of this file.
2 
3 #include "G4VPhysicalVolume.hh"
4 #include "G4Step.hh"
5 #include "G4Track.hh"
6 #include "G4Navigator.hh"
7 #include "G4NavigationHistory.hh"
8 #include "CLHEP/Units/PhysicalConstants.h"
9 #include "CLHEP/Units/SystemOfUnits.h"
10 
11 #include "Randomize.hh"
12 #include "G4TransportationManager.hh"
13 #include "G4VPhysicalVolume.hh"
14 #include "G4LogicalVolume.hh"
15 #include "G4VSensitiveDetector.hh"
16 #include "G4EventManager.hh"
17 #include "G4SteppingManager.hh"
18 #include "G4FastTrack.hh"
19 #include "G4ParticleTable.hh"
20 
21 #include "CLHEP/GenericFunctions/IncompleteGamma.hh"
22 
26 
29 
30 #include "CLHEP/Units/GlobalPhysicalConstants.h"
31 #include "CLHEP/Units/GlobalSystemOfUnits.h"
32 
33 #include <cmath>
34 
35 //#define EDM_ML_DEBUG
36 
38  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFGflash");
39  theBField = m_HF.getUntrackedParameter<double>("BField", 3.8);
40  theWatcherOn = m_HF.getUntrackedParameter<bool>("WatcherOn", true);
41  theFillHisto = m_HF.getUntrackedParameter<bool>("FillHisto", true);
42  edm::LogVerbatim("HFShower") << "HFGFlash:: Set B-Field to " << theBField << ", WatcherOn to " << theWatcherOn
43  << " and FillHisto to " << theFillHisto;
44 
46  theGflashStep = new G4Step();
47  theGflashNavigator = new G4Navigator();
48  // theGflashNavigator = 0;
49  theGflashTouchableHandle = new G4TouchableHistory();
50 
51 #ifdef EDM_ML_DEBUG
52  if (theFillHisto) {
54  if (tfile.isAvailable()) {
55  TFileDirectory showerDir = tfile->mkdir("GflashEMShowerProfile");
56 
57  em_incE = showerDir.make<TH1F>("em_incE", "Incoming energy (GeV)", 500, 0, 500.);
58  em_ssp_rho =
59  showerDir.make<TH1F>("em_ssp_rho", "Shower starting position;#rho (cm);Number of Events", 100, 100.0, 200.0);
60  em_ssp_z =
61  showerDir.make<TH1F>("em_ssp_z", "Shower starting position;z (cm);Number of Events", 2000, 0.0, 2000.0);
62  em_long =
63  showerDir.make<TH1F>("em_long", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
64  em_lateral = showerDir.make<TH1F>("em_lateral", "Lateral Profile;Radiation Length;Moliere Radius", 100, 0.0, 5.0);
65  em_2d = showerDir.make<TH2F>("em_2d",
66  "Lateral Profile vs. Shower Depth;Radiation Length;Moliere Radius",
67  800,
68  800.0,
69  1600.0,
70  100,
71  0.0,
72  5.0);
73  em_long_sd = showerDir.make<TH1F>("em_long_sd",
74  "Longitudinal Profile in Sensitive Detector;Radiation Length;Number of Spots",
75  800,
76  800.0,
77  1600.0);
79  showerDir.make<TH1F>("em_lateral_sd",
80  "Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",
81  100,
82  0.0,
83  5.0);
84  em_2d_sd =
85  showerDir.make<TH2F>("em_2d_sd",
86  "Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",
87  800,
88  800.0,
89  1600.0,
90  100,
91  0.0,
92  5.0);
93  em_ze = showerDir.make<TH2F>(
94  "em_ze", "Profile vs. Energy;Radiation Length;Moliere Radius", 800, 800.0, 1600.0, 1000, 0.0, 1.0);
95  em_ratio = showerDir.make<TH2F>(
96  "em_ratio", "Profile vs. Energy;Radiation Length;Moliere Radius", 800, 800.0, 1600.0, 1000, 0.0, 100.0);
97  em_ratio_selected = showerDir.make<TH2F>("em_ratio_selected",
98  "Profile vs. Energy;Radiation Length;Moliere Radius",
99  800,
100  800.0,
101  1600.0,
102  1000,
103  0.0,
104  100.0);
105  em_nSpots_sd =
106  showerDir.make<TH1F>("em_nSpots_sd",
107  "Number of Gflash Spots in Sensitive Detector;Number of Spots;Number of Events",
108  100,
109  0.0,
110  100);
111  em_ze_ratio = showerDir.make<TH1F>("em_ze_ratio", "Ratio of Energy and Z Position", 1000, 0.0, 0.001);
112  } else {
113  theFillHisto = false;
114  edm::LogVerbatim("HFShower") << "HFGFlash::No file is available for saving"
115  << " histos so the flag is set to false";
116  }
117  }
118 #endif
120 }
121 
123  delete theHelix;
124  delete theGflashStep;
125  delete theGflashNavigator;
126 }
127 
128 std::vector<HFGflash::Hit> HFGflash::gfParameterization(const G4Step* aStep, bool onlyLong) {
129  double tempZCalo = 26; // Gflash::Z[jCalorimeter]
130  double hfcriticalEnergy = 0.021; // Gflash::criticalEnergy
131 
132  std::vector<HFGflash::Hit> hit;
133  HFGflash::Hit oneHit;
134 
135  auto const preStepPoint = aStep->GetPreStepPoint();
136  auto const track = aStep->GetTrack();
137 
138  // This part of code is copied from the original GFlash Fortran code.
139  // reference : hep-ex/0001020v1
140 
141  const G4double energyCutoff = 1;
142  const G4int maxNumberOfSpots = 10000000;
143 
144  G4ThreeVector showerStartingPosition = track->GetPosition() / cm;
145  G4ThreeVector showerMomentum = preStepPoint->GetMomentum() / GeV;
147 
148  const G4double invgev = 1.0 / CLHEP::GeV;
149  G4double energy = preStepPoint->GetTotalEnergy() * invgev; // energy left in GeV
150  G4double logEinc = std::log(energy);
151 
152  G4double y = energy / hfcriticalEnergy; // y = E/Ec, criticalEnergy is in GeV
153  G4double logY = std::log(y);
154 
155  G4double nSpots = 93.0 * std::log(tempZCalo) * energy; // total number of spot due linearization
156  if (energy < 1.6)
157  nSpots = 140.4 * std::log(tempZCalo) * std::pow(energy, 0.876);
158 
159  // // implementing magnetic field effects
160  double charge = track->GetStep()->GetPreStepPoint()->GetCharge();
161  theHelix->initializeTrajectory(showerMomentum, showerStartingPosition, charge, theBField);
162 
163  G4double pathLength0 = theHelix->getPathLengthAtZ(showerStartingPosition.getZ());
164  G4double pathLength = pathLength0; // this will grow along the shower development
165 
166  //--- 2.2 Fix intrinsic properties of em. showers.
167 
168  G4double fluctuatedTmax = std::log(logY - 0.7157);
169  G4double fluctuatedAlpha = std::log(0.7996 + (0.4581 + 1.8628 / tempZCalo) * logY);
170 
171  G4double sigmaTmax = 1.0 / (-1.4 + 1.26 * logY);
172  G4double sigmaAlpha = 1.0 / (-0.58 + 0.86 * logY);
173  G4double rho = 0.705 - 0.023 * logY;
174  G4double sqrtPL = std::sqrt((1.0 + rho) / 2.0);
175  G4double sqrtLE = std::sqrt((1.0 - rho) / 2.0);
176 
177  G4double norm1 = G4RandGauss::shoot();
178  G4double norm2 = G4RandGauss::shoot();
179  G4double tempTmax = fluctuatedTmax + sigmaTmax * (sqrtPL * norm1 + sqrtLE * norm2);
180  G4double tempAlpha = fluctuatedAlpha + sigmaAlpha * (sqrtPL * norm1 - sqrtLE * norm2);
181 
182  // tmax, alpha, beta : parameters of gamma distribution
183  G4double tmax = std::exp(tempTmax);
184  G4double alpha = std::exp(tempAlpha);
185  G4double beta = (alpha - 1.0) / tmax;
186 
187  if (!alpha)
188  return hit;
189  if (!beta)
190  return hit;
191  if (alpha < 0.00001)
192  return hit;
193  if (beta < 0.00001)
194  return hit;
195 
196  // spot fluctuations are added to tmax, alpha, beta
197  G4double averageTmax = logY - 0.858;
198  G4double averageAlpha = 0.21 + (0.492 + 2.38 / tempZCalo) * logY;
199  G4double spotTmax = averageTmax * (0.698 + .00212 * tempZCalo);
200  G4double spotAlpha = averageAlpha * (0.639 + .00334 * tempZCalo);
201  G4double spotBeta = (spotAlpha - 1.0) / spotTmax;
202 
203  if (!spotAlpha)
204  return hit;
205  if (!spotBeta)
206  return hit;
207  if (spotAlpha < 0.00001)
208  return hit;
209  if (spotBeta < 0.00001)
210  return hit;
211 
212 #ifdef EDM_ML_DEBUG
213  edm::LogVerbatim("HFShower") << "HFGflash: Incoming energy = " << energy << " Position (rho,z) = ("
214  << showerStartingPosition.rho() << ", " << showerStartingPosition.z() << ")";
215 
216  if (theFillHisto) {
217  em_incE->Fill(energy);
218  em_ssp_rho->Fill(showerStartingPosition.rho());
219  em_ssp_z->Fill(std::abs(showerStartingPosition.z()));
220  }
221 #endif
222  // parameters for lateral distribution and fluctuation
223  G4double z1 = 0.0251 + 0.00319 * logEinc;
224  G4double z2 = 0.1162 - 0.000381 * tempZCalo;
225 
226  G4double k1 = 0.659 - 0.00309 * tempZCalo;
227  G4double k2 = 0.645;
228  G4double k3 = -2.59;
229  G4double k4 = 0.3585 + 0.0421 * logEinc;
230 
231  G4double p1 = 2.623 - 0.00094 * tempZCalo;
232  G4double p2 = 0.401 + 0.00187 * tempZCalo;
233  G4double p3 = 1.313 - 0.0686 * logEinc;
234 
235  // // @@@ dwjang, intial tuning by comparing 20-150GeV TB data
236  // // the width of energy response is not yet tuned.
237  const G4double e25Scale = 1.03551;
238  z1 *= 9.76972e-01 - 3.85026e-01 * std::tanh(1.82790e+00 * std::log(energy) - 3.66237e+00);
239  p1 *= 0.96;
240 
241  G4double stepLengthLeft = 10000;
242  G4int nSpots_sd = 0; // count total number of spots in SD
243  G4double zInX0 = 0.0; // shower depth in X0 unit
244  G4double deltaZInX0 = 0.0; // segment of depth in X0 unit
245  G4double deltaZ = 0.0; // segment of depth in cm
246  G4double stepLengthLeftInX0 = 0.0; // step length left in X0 unit
247 
248  const G4double divisionStepInX0 = 0.1; //step size in X0 unit
249 
250  Genfun::IncompleteGamma gammaDist;
251 
252  G4double energyInGamma = 0.0; // energy in a specific depth(z) according to Gamma distribution
253  G4double preEnergyInGamma = 0.0; // energy calculated in a previous depth
254  G4double sigmaInGamma = 0.; // sigma of energy in a specific depth(z) according to Gamma distribution
255  G4double preSigmaInGamma = 0.0; // sigma of energy in a previous depth
256 
257  //energy segment in Gamma distribution of shower in each step
258  G4double deltaEnergy = 0.0; // energy in deltaZ
259  G4int spotCounter = 0; // keep track of number of spots generated
260 
261  //step increment along the shower direction
262  G4double deltaStep = 0.0;
263 
264  // Uniqueness of G4Step is important otherwise hits won't be created.
265  G4double timeGlobal = preStepPoint->GetGlobalTime();
266 
267  // this needs to be deleted manually at the end of this loop.
268  // theGflashNavigator = new G4Navigator();
269  theGflashNavigator->SetWorldVolume(
270  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
271 
272  // // loop for longitudinal integration
273 
274 #ifdef EDM_ML_DEBUG
275  edm::LogVerbatim("HFShower") << "HFGflash: Energy = " << energy << " Step Length Left = " << stepLengthLeft;
276 #endif
277  while (energy > 0.0 && stepLengthLeft > 0.0) {
278  stepLengthLeftInX0 = stepLengthLeft / Gflash::radLength[jCalorimeter];
279 
280  if (stepLengthLeftInX0 < divisionStepInX0) {
281  deltaZInX0 = stepLengthLeftInX0;
282  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
283  stepLengthLeft = 0.0;
284  } else {
285  deltaZInX0 = divisionStepInX0;
286  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
287  stepLengthLeft -= deltaZ;
288  }
289 
290  zInX0 += deltaZInX0;
291 
292 #ifdef EDM_ML_DEBUG
293  edm::LogVerbatim("HFShower") << "HFGflash: zInX0 = " << zInX0 << " spotBeta*zInX0 = " << spotBeta * zInX0;
294 #endif
295  if ((!zInX0) || (!(spotBeta * zInX0 != 0)) || (zInX0 < 0.01) || (spotBeta * zInX0 < 0.00001) ||
296  (!(zInX0 * beta != 0)) || (zInX0 * beta < 0.00001))
297  return hit;
298 
299  G4int nSpotsInStep = 0;
300 
301 #ifdef EDM_ML_DEBUG
302  edm::LogVerbatim("HFShower") << "HFGflash: Energy - Energy Cut off = " << energy - energyCutoff;
303 #endif
304 
305  if (energy > energyCutoff) {
306  preEnergyInGamma = energyInGamma;
307  gammaDist.a().setValue(alpha); //alpha
308 
309  energyInGamma = gammaDist(beta * zInX0); //beta
310  G4double energyInDeltaZ = energyInGamma - preEnergyInGamma;
311  deltaEnergy = std::min(energy, energy * energyInDeltaZ);
312 
313  preSigmaInGamma = sigmaInGamma;
314  gammaDist.a().setValue(spotAlpha); //alpha spot
315  sigmaInGamma = gammaDist(spotBeta * zInX0); //beta spot
316  nSpotsInStep = std::max(1, int(nSpots * (sigmaInGamma - preSigmaInGamma)));
317  } else {
318  deltaEnergy = energy;
319  preSigmaInGamma = sigmaInGamma;
320  nSpotsInStep = std::max(1, int(nSpots * (1.0 - preSigmaInGamma)));
321  }
322 
323  if (deltaEnergy > energy || (energy - deltaEnergy) < energyCutoff)
324  deltaEnergy = energy;
325 
326  energy -= deltaEnergy;
327 
328  if (spotCounter + nSpotsInStep > maxNumberOfSpots) {
329  nSpotsInStep = maxNumberOfSpots - spotCounter;
330  if (nSpotsInStep < 1)
331  nSpotsInStep = 1;
332  }
333 
334  // // It begins with 0.5 of deltaZ and then icreases by 1 deltaZ
335  deltaStep += 0.5 * deltaZ;
336  pathLength += deltaStep;
337  deltaStep = 0.5 * deltaZ;
338 
339  //lateral shape and fluctuations for homogenous calo.
340  G4double tScale = tmax * alpha / (alpha - 1.0) * (1.0 - std::exp(-fluctuatedAlpha));
341  G4double tau = std::min(10.0, (zInX0 - 0.5 * deltaZInX0) / tScale);
342  G4double rCore = z1 + z2 * tau;
343  G4double rTail = k1 * (std::exp(k3 * (tau - k2)) + std::exp(k4 * (tau - k2))); // @@ check RT3 sign
344  G4double p23 = (p2 - tau) / p3;
345  G4double probabilityWeight = p1 * std::exp(p23 - std::exp(p23));
346 
347  // Deposition of spots according to lateral distr.
348  // Apply absolute energy scale
349  // Convert into MeV unit
350  G4double emSpotEnergy = deltaEnergy / (nSpotsInStep * e25Scale * GeV);
351 
352 #ifdef EDM_ML_DEBUG
353  edm::LogVerbatim("HFShower") << "HFGflash: nSpotsInStep = " << nSpotsInStep;
354 #endif
355 
356  for (G4int ispot = 0; ispot < nSpotsInStep; ispot++) {
357  spotCounter++;
358  G4double u1 = G4UniformRand();
359  G4double u2 = G4UniformRand();
360  G4double rInRM = 0.0;
361 
362  if (u1 < probabilityWeight) {
363  rInRM = rCore * std::sqrt(u2 / (1.0 - u2));
364  } else {
365  rInRM = rTail * std::sqrt(u2 / (1.0 - u2));
366  }
367 
368  G4double rShower = rInRM * Gflash::rMoliere[jCalorimeter];
369 
370  //Uniform & random rotation of spot along the azimuthal angle
371  G4double azimuthalAngle = twopi * G4UniformRand();
372 
373  //Compute global position of generated spots with taking into account magnetic field
374  //Divide deltaZ into nSpotsInStep and give a spot a global position
375  G4double incrementPath = (deltaZ / nSpotsInStep) * (ispot + 0.5 - 0.5 * nSpotsInStep);
376 
377  // trajectoryPoint give a spot an imaginary point along the shower development
378  GflashTrajectoryPoint trajectoryPoint;
379  theHelix->getGflashTrajectoryPoint(trajectoryPoint, pathLength + incrementPath);
380 
381  G4ThreeVector SpotPosition0 = trajectoryPoint.getPosition() +
382  rShower * std::cos(azimuthalAngle) * trajectoryPoint.getOrthogonalUnitVector() +
383  rShower * std::sin(azimuthalAngle) * trajectoryPoint.getCrossUnitVector();
384 
386  // Convert into mm unit
387  SpotPosition0 *= cm;
388 
389  //---------------------------------------------------
390  // fill a fake step to send it to hit maker
391  //---------------------------------------------------
392 
393  // to make a different time for each fake step. (0.03 nsec is corresponding to 1cm step size)
394  timeGlobal += 0.0001 * nanosecond;
395 
396  //fill equivalent changes to a (fake) step associated with a spot
397 
398  G4double zInX0_spot = std::abs(pathLength + incrementPath - pathLength0) / Gflash::radLength[jCalorimeter];
399 
400 #ifdef EDM_ML_DEBUG
401  edm::LogVerbatim("HFShower") << "HFGflash: zInX0_spot,emSpotEnergy/GeV =" << zInX0_spot << " , "
402  << emSpotEnergy / GeV << "emSpotEnergy/GeV =" << emSpotEnergy / GeV;
403 #endif
404 
405  if (zInX0_spot <= 0)
406  continue;
407  if (emSpotEnergy <= 0)
408  continue;
409  if (rShower / Gflash::rMoliere[jCalorimeter] <= 0)
410  continue;
411 
412  oneHit.depth = 1;
413 
414 #ifdef EDM_ML_DEBUG
415  if (theFillHisto) {
416  em_long->Fill(SpotPosition0.z() / cm, emSpotEnergy * invgev);
417  em_lateral->Fill(rShower / Gflash::rMoliere[jCalorimeter], emSpotEnergy * invgev);
418  em_2d->Fill(SpotPosition0.z() / cm, rShower / Gflash::rMoliere[jCalorimeter], emSpotEnergy * invgev);
419  }
420 #endif
421 
423  //if(SpotPosition0 == 0) continue;
424 
425  double energyratio = emSpotEnergy / (preStepPoint->GetTotalEnergy() / (nSpots * e25Scale));
426 
427  if (emSpotEnergy / GeV < 0.0001)
428  continue;
429  if (energyratio > 80)
430  continue;
431 
432  double zshift = 0;
433  if (SpotPosition0.z() > 0)
434  zshift = 18;
435  if (SpotPosition0.z() < 0)
436  zshift = -18;
437 
438  G4ThreeVector gfshift(0, 0, zshift * (pow(100, 0.1) / pow(energy, 0.1)));
439 
440  G4ThreeVector SpotPosition = gfshift + SpotPosition0;
441 
442  double LengthWeight = std::fabs(std::pow(SpotPosition0.z() / 11370, 1));
443  if (G4UniformRand() > 0.0021 * energyratio * LengthWeight)
444  continue;
445 
446  oneHit.position = SpotPosition;
447  oneHit.time = timeGlobal;
448  oneHit.edep = emSpotEnergy * invgev;
449  hit.push_back(oneHit);
450  nSpots_sd++;
451 
452  } // end of for spot iteration
453 
454  } // end of while for longitudinal integration
455 #ifdef EDM_ML_DEBUG
456  if (theFillHisto) {
457  em_nSpots_sd->Fill(nSpots_sd);
458  }
459 #endif
460  // delete theGflashNavigator;
461  return hit;
462 }
HFGflash::em_ssp_rho
TH1F * em_ssp_rho
Definition: HFGflash.h:64
GflashEMShowerProfile.h
DDAxes::y
Gflash::radLength
const double radLength[kNumberCalorimeter]
Definition: GflashNameSpace.h:57
HFGflash::em_ratio_selected
TH2F * em_ratio_selected
Definition: HFGflash.h:66
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11724
HFGflash::theGflashStep
G4Step * theGflashStep
Definition: HFGflash.h:47
SteppingAction.h
HFGflash.h
HFGflash::Hit::edep
double edep
Definition: HFGflash.h:39
invgev
static const float invgev
Definition: TimingSD.cc:32
GflashTrajectoryPoint.h
metsig::tau
Definition: SignAlgoResolutions.h:49
min
T min(T a, T b)
Definition: MathUtil.h:58
HLT_FULL_cff.beta
beta
Definition: HLT_FULL_cff.py:8651
HFGflash::Hit::time
double time
Definition: HFGflash.h:38
HFGflash::em_ze_ratio
TH1F * em_ze_ratio
Definition: HFGflash.h:65
HFGflash::~HFGflash
virtual ~HFGflash()
Definition: HFGflash.cc:122
alpha
float alpha
Definition: AMPTWrapper.h:105
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
TFileDirectory
Definition: TFileDirectory.h:24
HFGflash::HFGflash
HFGflash(edm::ParameterSet const &p)
Definition: HFGflash.cc:37
GflashTrajectory::initializeTrajectory
void initializeTrajectory(const HepGeom::Vector3D< double > &, const HepGeom::Point3D< double > &, double q, double Field)
Definition: GflashTrajectory.cc:27
GflashTrajectoryPoint::getOrthogonalUnitVector
Gflash3Vector getOrthogonalUnitVector()
Definition: GflashTrajectoryPoint.h:22
tmax
static const double tmax[3]
Definition: CastorTimeSlew.cc:7
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
testProducerWithPsetDescEmpty_cfi.z2
z2
Definition: testProducerWithPsetDescEmpty_cfi.py:41
GflashTrajectory
Definition: GflashTrajectory.h:6
HFGflash::Hit::position
G4ThreeVector position
Definition: HFGflash.h:36
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Service.h
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
tfile
Definition: tfile.py:1
HFGflash::theGflashTouchableHandle
G4TouchableHandle theGflashTouchableHandle
Definition: HFGflash.h:49
p2
double p2[4]
Definition: TauolaWrapper.h:90
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
Gflash::kHF
Definition: GflashNameSpace.h:17
HFGflash::em_long_sd
TH1F * em_long_sd
Definition: HFGflash.h:65
HFGflash::jCalorimeter
Gflash::CalorimeterNumber jCalorimeter
Definition: HFGflash.h:51
HFGflash::gfParameterization
std::vector< Hit > gfParameterization(const G4Step *aStep, bool onlyLong=false)
Definition: HFGflash.cc:128
DDAxes::rho
TFileService.h
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
HFGflash::em_lateral_sd
TH1F * em_lateral_sd
Definition: HFGflash.h:65
HFGflash::theBField
G4double theBField
Definition: HFGflash.h:55
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
listHistos.logY
logY
Definition: listHistos.py:42
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
GflashTrajectoryPoint::getPosition
Gflash3Vector & getPosition()
Definition: GflashTrajectoryPoint.h:20
GeV
const double GeV
Definition: MathUtil.h:16
HFGflash::theWatcherOn
bool theWatcherOn
Definition: HFGflash.h:53
HFGflash::em_incE
TH1F * em_incE
Definition: HFGflash.h:64
HFGflash::Hit
Definition: HFGflash.h:34
HFGflash::em_long
TH1F * em_long
Definition: HFGflash.h:64
edm::Service< TFileService >
testProducerWithPsetDescEmpty_cfi.u1
u1
Definition: testProducerWithPsetDescEmpty_cfi.py:49
HFGflash::em_ssp_z
TH1F * em_ssp_z
Definition: HFGflash.h:64
GflashTrajectoryPoint::getCrossUnitVector
Gflash3Vector getCrossUnitVector()
Definition: GflashTrajectoryPoint.h:23
HFGflash::em_2d_sd
TH2F * em_2d_sd
Definition: HFGflash.h:66
Gflash::rMoliere
const double rMoliere[kNumberCalorimeter]
Definition: GflashNameSpace.h:56
p1
double p1[4]
Definition: TauolaWrapper.h:89
HFGflash::theHelix
GflashTrajectory * theHelix
Definition: HFGflash.h:46
compare.tfile
tfile
Definition: compare.py:324
HFGflash::theGflashNavigator
G4Navigator * theGflashNavigator
Definition: HFGflash.h:48
L1TkHTMissProducer_cfi.deltaZ
deltaZ
Definition: L1TkHTMissProducer_cfi.py:13
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HFGflash::em_ratio
TH2F * em_ratio
Definition: HFGflash.h:66
TFileDirectory::make
T * make(const Args &...args) const
make new ROOT object
Definition: TFileDirectory.h:53
HFGflash::theFillHisto
bool theFillHisto
Definition: HFGflash.h:54
p3
double p3[4]
Definition: TauolaWrapper.h:91
HFGflash::em_lateral
TH1F * em_lateral
Definition: HFGflash.h:64
GflashTrajectory::getPathLengthAtZ
double getPathLengthAtZ(double z) const
Definition: GflashTrajectory.cc:179
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
HFGflash::em_nSpots_sd
TH1F * em_nSpots_sd
Definition: HFGflash.h:65
MetAnalyzer.u2
u2
Definition: MetAnalyzer.py:61
HFGflash::em_2d
TH2F * em_2d
Definition: HFGflash.h:66
HFGflash::Hit::depth
int depth
Definition: HFGflash.h:37
GflashTrajectoryPoint
Definition: GflashTrajectoryPoint.h:8
hit
Definition: SiStripHitEffFromCalibTree.cc:88
HFGflash::em_ze
TH2F * em_ze
Definition: HFGflash.h:66
GflashTrajectory::getGflashTrajectoryPoint
void getGflashTrajectoryPoint(GflashTrajectoryPoint &point, double s) const
Definition: GflashTrajectory.cc:155
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37