CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Attributes
HFGflash Class Reference

#include <HFGflash.h>

Classes

struct  Hit
 

Public Member Functions

std::vector< HitgfParameterization (G4Step *aStep, bool &ok, bool onlyLong=false)
 
 HFGflash (edm::ParameterSet const &p)
 
virtual ~HFGflash ()
 

Private Attributes

TH2F * em_2d
 
TH2F * em_2d_sd
 
TH1F * em_incE
 
TH1F * em_lateral
 
TH1F * em_lateral_sd
 
TH1F * em_long
 
TH1F * em_long_sd
 
TH1F * em_nSpots_sd
 
TH2F * em_ratio
 
TH2F * em_ratio_selected
 
TH1F * em_ssp_rho
 
TH1F * em_ssp_z
 
TH2F * em_ze
 
TH1F * em_ze_ratio
 
G4double energyScale [Gflash::kNumberCalorimeter]
 
G4double energyToDeposit
 
Gflash::CalorimeterNumber jCalorimeter
 
G4double lateralPar [4]
 
G4double longEcal [Gflash::NPar]
 
G4double longHcal [Gflash::NPar]
 
G4int showerType
 
G4double theBField
 
bool theFillHisto
 
G4Navigator * theGflashNavigator
 
G4Step * theGflashStep
 
G4TouchableHandle theGflashTouchableHandle
 
GflashTrajectorytheHelix
 
bool theWatcherOn
 

Detailed Description

Definition at line 29 of file HFGflash.h.

Constructor & Destructor Documentation

HFGflash::HFGflash ( edm::ParameterSet const &  p)

Definition at line 37 of file HFGflash.cc.

References em_2d, em_2d_sd, em_incE, em_lateral, em_lateral_sd, em_long, em_long_sd, em_nSpots_sd, em_ratio, em_ratio_selected, em_ssp_rho, em_ssp_z, em_ze, em_ze_ratio, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), edm::Service< T >::isAvailable(), jCalorimeter, Gflash::kHF, TFileDirectory::make(), TFileService::mkdir(), theBField, theFillHisto, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theHelix, and theWatcherOn.

37  {
38 
39  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFGflash");
40  theBField = m_HF.getUntrackedParameter<double>("BField", 3.8);
41  theWatcherOn = m_HF.getUntrackedParameter<bool>("WatcherOn",true);
42  theFillHisto = m_HF.getUntrackedParameter<bool>("FillHisto",true);
43  edm::LogInfo("HFShower") << "HFGFlash:: Set B-Field to " << theBField
44  << ", WatcherOn to " << theWatcherOn
45  << " and FillHisto to " << theFillHisto;
46 
48  theGflashStep = new G4Step();
49  theGflashNavigator = new G4Navigator();
50  // theGflashNavigator = 0;
51  theGflashTouchableHandle = new G4TouchableHistory();
52 
53 #ifdef DebugLog
54  if (theFillHisto) {
56  if ( tfile.isAvailable() ) {
57  TFileDirectory showerDir = tfile->mkdir("GflashEMShowerProfile");
58 
59  em_incE = showerDir.make<TH1F>("em_incE","Incoming energy (GeV)",500,0,500.);
60  em_ssp_rho = showerDir.make<TH1F>("em_ssp_rho","Shower starting position;#rho (cm);Number of Events",100,100.0,200.0);
61  em_ssp_z = showerDir.make<TH1F>("em_ssp_z","Shower starting position;z (cm);Number of Events",2000,0.0,2000.0);
62  em_long = showerDir.make<TH1F>("em_long","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
63  em_lateral = showerDir.make<TH1F>("em_lateral","Lateral Profile;Radiation Length;Moliere Radius",100,0.0,5.0);
64  em_2d = showerDir.make<TH2F>("em_2d","Lateral Profile vs. Shower Depth;Radiation Length;Moliere Radius",800,800.0,1600.0,100,0.0,5.0);
65  em_long_sd = showerDir.make<TH1F>("em_long_sd","Longitudinal Profile in Sensitive Detector;Radiation Length;Number of Spots",800,800.0,1600.0);
66  em_lateral_sd = showerDir.make<TH1F>("em_lateral_sd","Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",100,0.0,5.0);
67  em_2d_sd = showerDir.make<TH2F>("em_2d_sd","Lateral Profile vs. Shower Depth in Sensitive Detector;Radiation Length;Moliere Radius",800,800.0,1600.0,100,0.0,5.0);
68  em_ze = showerDir.make<TH2F>("em_ze","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,1.0);
69  em_ratio = showerDir.make<TH2F>("em_ratio","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,100.0);
70  em_ratio_selected = showerDir.make<TH2F>("em_ratio_selected","Profile vs. Energy;Radiation Length;Moliere Radius",800,800.0,1600.0,1000,0.0,100.0);
71  em_nSpots_sd = showerDir.make<TH1F>("em_nSpots_sd","Number of Gflash Spots in Sensitive Detector;Number of Spots;Number of Events",100,0.0,100);
72  em_ze_ratio = showerDir.make<TH1F>("em_ze_ratio","Ratio of Energy and Z Position",1000,0.0,0.001);
73  } else {
74  theFillHisto = false;
75  edm::LogInfo("HFShower") << "HFGFlash::No file is available for saving"
76  << " histos so the flag is set to false";
77  }
78  }
79 #endif
81 
82 }
G4Step * theGflashStep
Definition: HFGflash.h:50
T getUntrackedParameter(std::string const &, T const &) const
TH2F * em_ze
Definition: HFGflash.h:69
G4Navigator * theGflashNavigator
Definition: HFGflash.h:51
TH1F * em_ssp_z
Definition: HFGflash.h:67
Gflash::CalorimeterNumber jCalorimeter
Definition: HFGflash.h:54
G4double theBField
Definition: HFGflash.h:58
bool theFillHisto
Definition: HFGflash.h:57
TH2F * em_ratio_selected
Definition: HFGflash.h:69
bool isAvailable() const
Definition: Service.h:46
G4TouchableHandle theGflashTouchableHandle
Definition: HFGflash.h:52
T * make(const Args &...args) const
make new ROOT object
bool theWatcherOn
Definition: HFGflash.h:56
TH2F * em_ratio
Definition: HFGflash.h:69
TH1F * em_ze_ratio
Definition: HFGflash.h:68
TH1F * em_ssp_rho
Definition: HFGflash.h:67
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
TH1F * em_lateral_sd
Definition: HFGflash.h:68
TH1F * em_nSpots_sd
Definition: HFGflash.h:68
TH1F * em_lateral
Definition: HFGflash.h:67
TH1F * em_long
Definition: HFGflash.h:67
TH1F * em_incE
Definition: HFGflash.h:67
GflashTrajectory * theHelix
Definition: HFGflash.h:49
TH1F * em_long_sd
Definition: HFGflash.h:68
TH2F * em_2d
Definition: HFGflash.h:69
TH2F * em_2d_sd
Definition: HFGflash.h:69
HFGflash::~HFGflash ( )
virtual

Definition at line 84 of file HFGflash.cc.

References theGflashNavigator, theGflashStep, and theHelix.

84  {
85  if (theHelix) delete theHelix;
86  if (theGflashStep) delete theGflashStep;
88 }
G4Step * theGflashStep
Definition: HFGflash.h:50
G4Navigator * theGflashNavigator
Definition: HFGflash.h:51
GflashTrajectory * theHelix
Definition: HFGflash.h:49

Member Function Documentation

std::vector< HFGflash::Hit > HFGflash::gfParameterization ( G4Step *  aStep,
bool &  ok,
bool  onlyLong = false 
)

V.Ivanchenko - not clear if it is correct

V.Ivanchenko what this means??

Definition at line 91 of file HFGflash.cc.

References funct::abs(), alpha, beta, RecoTauCleanerPlugins::charge, funct::cos(), HFGflash::Hit::depth, alignCSCRings::e, HFGflash::Hit::edep, em_2d, em_incE, em_lateral, em_long, em_nSpots_sd, em_ssp_rho, em_ssp_z, relval_parameters_module::energy, create_public_lumi_plots::exp, GflashTrajectoryPoint::getCrossUnitVector(), GflashTrajectory::getGflashTrajectoryPoint(), GflashTrajectoryPoint::getOrthogonalUnitVector(), GflashTrajectory::getPathLengthAtZ(), GflashTrajectoryPoint::getPosition(), GeV, GflashTrajectory::initializeTrajectory(), jCalorimeter, relval_steps::k2, Gflash::kHF, cmsBatch::log, LogDebug, bookConverter::max, min(), p1, p2, p3, HFGflash::Hit::position, funct::pow(), Gflash::radLength, rho, Gflash::rMoliere, funct::sin(), mathSSE::sqrt(), metsig::tau, theBField, theFillHisto, theGflashNavigator, theHelix, HFGflash::Hit::time, tmax, MetAnalyzer::u1, MetAnalyzer::u2, and y.

Referenced by HFShowerParam::getHits().

91  {
92  double tempZCalo = 26; // Gflash::Z[jCalorimeter]
93  double hfcriticalEnergy = 0.021; // Gflash::criticalEnergy
94 
95  std::vector<HFGflash::Hit> hit;
96  HFGflash::Hit oneHit;
97 
98  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
99  //G4StepPoint * postStepPoint = aStep->GetPostStepPoint();
100  G4Track * track = aStep->GetTrack();
101  // Get Z-direction
102  const G4DynamicParticle *aParticle = track->GetDynamicParticle();
103  G4ThreeVector momDir = aParticle->GetMomentumDirection();
104 
105  G4ThreeVector hitPoint = preStepPoint->GetPosition();
106  G4String partType = track->GetDefinition()->GetParticleName();
107  // int parCode = track->GetDefinition()->GetPDGEncoding();
108 
109  // This part of code is copied from the original GFlash Fortran code.
110  // reference : hep-ex/0001020v1
111 
112  const G4double energyCutoff = 1;
113  const G4int maxNumberOfSpots = 10000000;
114 
115  G4ThreeVector showerStartingPosition = track->GetPosition()/cm;
116  G4ThreeVector showerMomentum = preStepPoint->GetMomentum()/GeV;
118 
119  G4double logEinc = std::log((preStepPoint->GetTotalEnergy())/GeV);
120 
121  G4double y = ((preStepPoint->GetTotalEnergy())/GeV) / hfcriticalEnergy; // y = E/Ec, criticalEnergy is in GeV
122  G4double logY = std::log(y);
123 
124 
125  G4double nSpots = 93.0 * std::log(tempZCalo) * ((preStepPoint->GetTotalEnergy())/GeV); // total number of spot due linearization
126  if(preStepPoint->GetTotalEnergy()/GeV < 1.6) nSpots = 140.4 * std::log(tempZCalo) * std::pow(((preStepPoint->GetTotalEnergy())/GeV),0.876);
127 
128 
129  // // implementing magnetic field effects
130  double charge = track->GetStep()->GetPreStepPoint()->GetCharge();
131  theHelix->initializeTrajectory(showerMomentum,showerStartingPosition,charge,theBField);
132 
133  G4double pathLength0 = theHelix->getPathLengthAtZ(showerStartingPosition.getZ());
134  G4double pathLength = pathLength0; // this will grow along the shower development
135 
136  //--- 2.2 Fix intrinsic properties of em. showers.
137 
138  G4double fluctuatedTmax = std::log(logY - 0.7157);
139  G4double fluctuatedAlpha= std::log(0.7996 +(0.4581 + 1.8628/tempZCalo)*logY);
140 
141  G4double sigmaTmax = 1.0/( -1.4 + 1.26 * logY);
142  G4double sigmaAlpha = 1.0/( -0.58 + 0.86 * logY);
143  G4double rho = 0.705 - 0.023 * logY;
144  G4double sqrtPL = std::sqrt((1.0+rho)/2.0);
145  G4double sqrtLE = std::sqrt((1.0-rho)/2.0);
146 
147  G4double norm1 = G4RandGauss::shoot();
148  G4double norm2 = G4RandGauss::shoot();
149  G4double tempTmax = fluctuatedTmax + sigmaTmax*(sqrtPL*norm1 + sqrtLE*norm2);
150  G4double tempAlpha = fluctuatedAlpha + sigmaAlpha*(sqrtPL*norm1 - sqrtLE*norm2);
151 
152  // tmax, alpha, beta : parameters of gamma distribution
153  G4double tmax = std::exp(tempTmax);
154  G4double alpha = std::exp(tempAlpha);
155  G4double beta = (alpha - 1.0)/tmax;
156 
157  if (!alpha) return hit;
158  if (!beta) return hit;
159  if (alpha < 0.00001) return hit;
160  if (beta < 0.00001) return hit;
161 
162  // spot fluctuations are added to tmax, alpha, beta
163  G4double averageTmax = logY-0.858;
164  G4double averageAlpha = 0.21+(0.492+2.38/tempZCalo)*logY;
165  G4double spotTmax = averageTmax * (0.698 + .00212*tempZCalo);
166  G4double spotAlpha= averageAlpha * (0.639 + .00334*tempZCalo);
167  G4double spotBeta = (spotAlpha-1.0)/spotTmax;
168 
169  if (!spotAlpha) return hit;
170  if (!spotBeta) return hit;
171  if (spotAlpha < 0.00001) return hit;
172  if (spotBeta < 0.00001) return hit;
173 
174 #ifdef DebugLog
175  LogDebug("HFShower") << "Incoming energy = " << ((preStepPoint->GetTotalEnergy())/GeV) << " Position (rho,z) = (" << showerStartingPosition.rho() << ", " << showerStartingPosition.z() << ")";
176 
177  if(theFillHisto) {
178  em_incE->Fill(((preStepPoint->GetTotalEnergy())/GeV));
179  em_ssp_rho->Fill(showerStartingPosition.rho());
180  em_ssp_z->Fill(std::abs(showerStartingPosition.z()));
181  }
182 #endif
183  // parameters for lateral distribution and fluctuation
184  G4double z1=0.0251+0.00319*logEinc;
185  G4double z2=0.1162-0.000381*tempZCalo;
186 
187  G4double k1=0.659 - 0.00309 * tempZCalo;
188  G4double k2=0.645;
189  G4double k3=-2.59;
190  G4double k4=0.3585+ 0.0421*logEinc;
191 
192  G4double p1=2.623 -0.00094*tempZCalo;
193  G4double p2=0.401 +0.00187*tempZCalo;
194  G4double p3=1.313 -0.0686*logEinc;
195 
196  // // @@@ dwjang, intial tuning by comparing 20-150GeV TB data
197  // // the width of energy response is not yet tuned.
198  G4double e25Scale = 1.03551;
199  z1 *= 9.76972e-01 - 3.85026e-01 * std::tanh(1.82790e+00*std::log(((preStepPoint->GetTotalEnergy())/GeV)) - 3.66237e+00);
200  p1 *= 0.96;
201 
202  G4double stepLengthLeft = 10000;
203  G4int nSpots_sd = 0; // count total number of spots in SD
204  G4double zInX0 = 0.0; // shower depth in X0 unit
205  G4double deltaZInX0 = 0.0; // segment of depth in X0 unit
206  G4double deltaZ = 0.0; // segment of depth in cm
207  G4double stepLengthLeftInX0 = 0.0; // step length left in X0 unit
208 
209  const G4double divisionStepInX0 = 0.1; //step size in X0 unit
210  G4double energy = ((preStepPoint->GetTotalEnergy())/GeV); // energy left in GeV
211 
212  Genfun::IncompleteGamma gammaDist;
213 
214  G4double energyInGamma = 0.0; // energy in a specific depth(z) according to Gamma distribution
215  G4double preEnergyInGamma = 0.0; // energy calculated in a previous depth
216  G4double sigmaInGamma = 0.; // sigma of energy in a specific depth(z) according to Gamma distribution
217  G4double preSigmaInGamma = 0.0; // sigma of energy in a previous depth
218 
219  //energy segment in Gamma distribution of shower in each step
220  G4double deltaEnergy =0.0 ; // energy in deltaZ
221  G4int spotCounter = 0; // keep track of number of spots generated
222 
223  //step increment along the shower direction
224  G4double deltaStep = 0.0;
225 
226  // Uniqueness of G4Step is important otherwise hits won't be created.
227  G4double timeGlobal = track->GetStep()->GetPreStepPoint()->GetGlobalTime();
228 
229  // this needs to be deleted manually at the end of this loop.
230  // theGflashNavigator = new G4Navigator();
231  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
232 
233  // // loop for longitudinal integration
234 
235 #ifdef DebugLog
236  LogDebug("HFShower") << " Energy = " << energy << " Step Length Left = " << stepLengthLeft;
237 #endif
238  while(energy > 0.0 && stepLengthLeft > 0.0) {
239  stepLengthLeftInX0 = stepLengthLeft / Gflash::radLength[jCalorimeter];
240 
241  if ( stepLengthLeftInX0 < divisionStepInX0 ) {
242  deltaZInX0 = stepLengthLeftInX0;
243  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
244  stepLengthLeft = 0.0;
245  } else {
246  deltaZInX0 = divisionStepInX0;
247  deltaZ = deltaZInX0 * Gflash::radLength[jCalorimeter];
248  stepLengthLeft -= deltaZ;
249  }
250 
251  zInX0 += deltaZInX0;
252 
253 #ifdef DebugLog
254  LogDebug("HFShower") << " zInX0 = " << zInX0 << " spotBeta*zInX0 = " << spotBeta*zInX0;
255 #endif
256  if ((!zInX0) || (!spotBeta*zInX0) || (zInX0 < 0.01) ||
257  (spotBeta*zInX0 < 0.00001) || (!zInX0*beta) || (zInX0*beta < 0.00001))
258  return hit;
259 
260  G4int nSpotsInStep = 0;
261 
262 #ifdef DebugLog
263  LogDebug("HFShower") << " Energy - Energy Cut off = " << energy - energyCutoff;
264 #endif
265 
266  if ( energy > energyCutoff ) {
267  preEnergyInGamma = energyInGamma;
268  gammaDist.a().setValue(alpha); //alpha
269 
270  energyInGamma = gammaDist(beta*zInX0); //beta
271  G4double energyInDeltaZ = energyInGamma - preEnergyInGamma;
272  deltaEnergy = std::min(energy,((preStepPoint->GetTotalEnergy())/GeV)*energyInDeltaZ);
273 
274  preSigmaInGamma = sigmaInGamma;
275  gammaDist.a().setValue(spotAlpha); //alpha spot
276  sigmaInGamma = gammaDist(spotBeta*zInX0); //beta spot
277  nSpotsInStep = std::max(1,int(nSpots * (sigmaInGamma - preSigmaInGamma)));
278  } else {
279  deltaEnergy = energy;
280  preSigmaInGamma = sigmaInGamma;
281  nSpotsInStep = std::max(1,int(nSpots * (1.0 - preSigmaInGamma)));
282  }
283 
284  if ( deltaEnergy > energy || (energy-deltaEnergy) < energyCutoff ) deltaEnergy = energy;
285 
286  energy -= deltaEnergy;
287 
288  if ( spotCounter+nSpotsInStep > maxNumberOfSpots ) {
289  nSpotsInStep = maxNumberOfSpots - spotCounter;
290  if (nSpotsInStep < 1) nSpotsInStep = 1;
291  }
292 
293 
294  // // It begins with 0.5 of deltaZ and then icreases by 1 deltaZ
295  deltaStep += 0.5*deltaZ;
296  pathLength += deltaStep;
297  deltaStep = 0.5*deltaZ;
298 
299 
300  //lateral shape and fluctuations for homogenous calo.
301  G4double tScale = tmax *alpha/(alpha-1.0) * (std::exp(fluctuatedAlpha)-1.0)/std::exp(fluctuatedAlpha);
302  G4double tau = std::min(10.0,(zInX0 - 0.5*deltaZInX0)/tScale);
303  G4double rCore = z1 + z2 * tau;
304  G4double rTail = k1 *( std::exp(k3*(tau-k2)) + std::exp(k4*(tau-k2))); // @@ check RT3 sign
305  G4double p23 = (p2 - tau)/p3;
306  G4double probabilityWeight = p1 * std::exp( p23 - std::exp(p23) );
307 
308 
309  // Deposition of spots according to lateral distr.
310  // Apply absolute energy scale
311  // Convert into MeV unit
312  G4double emSpotEnergy = deltaEnergy / nSpotsInStep * e25Scale * GeV;
313 
314 
315 #ifdef DebugLog
316  LogDebug("HFShower") << " nSpotsInStep = " << nSpotsInStep;
317 #endif
318 
319 
320 
321  for (G4int ispot = 0 ; ispot < nSpotsInStep ; ispot++) {
322  spotCounter++;
323  G4double u1 = G4UniformRand();
324  G4double u2 = G4UniformRand();
325  G4double rInRM = 0.0;
326 
327  if (u1 < probabilityWeight) {
328  rInRM = rCore * std::sqrt( u2/(1.0-u2) );
329  } else {
330  rInRM = rTail * std::sqrt( u2/(1.0-u2) );
331  }
332 
333  G4double rShower = rInRM * Gflash::rMoliere[jCalorimeter];
334 
335  //Uniform & random rotation of spot along the azimuthal angle
336  G4double azimuthalAngle = twopi*G4UniformRand();
337 
338  //Compute global position of generated spots with taking into account magnetic field
339  //Divide deltaZ into nSpotsInStep and give a spot a global position
340  G4double incrementPath = (deltaZ/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
341 
342  // trajectoryPoint give a spot an imaginary point along the shower development
343  GflashTrajectoryPoint trajectoryPoint;
344  theHelix->getGflashTrajectoryPoint(trajectoryPoint,pathLength+incrementPath);
345 
346 
347  G4ThreeVector SpotPosition0 = trajectoryPoint.getPosition() + rShower*std::cos(azimuthalAngle)*trajectoryPoint.getOrthogonalUnitVector() + rShower*std::sin(azimuthalAngle)*trajectoryPoint.getCrossUnitVector();
348 
350  // Convert into mm unit
351  SpotPosition0 *= cm;
352 
353 
354  //---------------------------------------------------
355  // fill a fake step to send it to hit maker
356  //---------------------------------------------------
357 
358  // to make a different time for each fake step. (0.03 nsec is corresponding to 1cm step size)
359  timeGlobal += 0.0001*nanosecond;
360 
361  //fill equivalent changes to a (fake) step associated with a spot
362 
363  G4double zInX0_spot = std::abs(pathLength+incrementPath - pathLength0)/Gflash::radLength[jCalorimeter];
364 
365 #ifdef DebugLog
366  LogDebug("HFShower") << "zInX0_spot,emSpotEnergy/GeV =" << zInX0_spot << " , " << emSpotEnergy/GeV << "emSpotEnergy/GeV =" << emSpotEnergy/GeV;
367 #endif
368 
369  if ((!zInX0_spot) || (zInX0_spot < 0)) continue;
370  if ((!emSpotEnergy/GeV) || (emSpotEnergy < 0)) continue;
371  if ((!rShower/Gflash::rMoliere[jCalorimeter]) || (rShower/Gflash::rMoliere[jCalorimeter] < 0)) continue;
372 
373  oneHit.depth = 1;
374 
375 #ifdef DebugLog
376  if (theFillHisto) {
377  em_long->Fill(SpotPosition0.z()/cm,emSpotEnergy/GeV);
378  em_lateral->Fill(rShower/Gflash::rMoliere[jCalorimeter],emSpotEnergy/GeV);
379  em_2d->Fill(SpotPosition0.z()/cm,rShower/Gflash::rMoliere[jCalorimeter],emSpotEnergy/GeV);
380  }
381 #endif
382 
384  //if(SpotPosition0 == 0) continue;
385 
386  double energyratio = emSpotEnergy/(preStepPoint->GetTotalEnergy()/(nSpots*e25Scale));
387 
388  if (emSpotEnergy/GeV < 0.0001) continue;
389  if (energyratio > 80) continue;
390 
391  double zshift =0;
392  if(SpotPosition0.z() > 0) zshift = 18;
393  if(SpotPosition0.z() < 0) zshift = -18;
394 
395 
396  G4ThreeVector gfshift(0,0,zshift*(pow(100,0.1)/pow(preStepPoint->GetTotalEnergy()/GeV,0.1)));
397 
398  G4ThreeVector SpotPosition = gfshift + SpotPosition0;
399 
400  double LengthWeight = std::fabs(std::pow(SpotPosition0.z()/11370,1));
401  if (G4UniformRand()> 0.0021 * energyratio * LengthWeight) continue;
402 
403  oneHit.position = SpotPosition;
404  oneHit.time = timeGlobal;
405  oneHit.edep = emSpotEnergy/GeV;
406  hit.push_back(oneHit);
407  nSpots_sd++;
408 
409  } // end of for spot iteration
410 
411  } // end of while for longitudinal integration
412 #ifdef DebugLog
413  if (theFillHisto) {
414  em_nSpots_sd->Fill(nSpots_sd);
415  }
416 #endif
417  // delete theGflashNavigator;
418  return hit;
419 }
#define LogDebug(id)
const double beta
G4Navigator * theGflashNavigator
Definition: HFGflash.h:51
float alpha
Definition: AMPTWrapper.h:95
const double GeV
Definition: MathUtil.h:16
double getPathLengthAtZ(double z) const
Gflash3Vector getCrossUnitVector()
TH1F * em_ssp_z
Definition: HFGflash.h:67
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Gflash::CalorimeterNumber jCalorimeter
Definition: HFGflash.h:54
G4double theBField
Definition: HFGflash.h:58
bool theFillHisto
Definition: HFGflash.h:57
const double rMoliere[kNumberCalorimeter]
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double edep
Definition: HFGflash.h:41
T min(T a, T b)
Definition: MathUtil.h:58
double p2[4]
Definition: TauolaWrapper.h:90
Gflash3Vector getOrthogonalUnitVector()
void getGflashTrajectoryPoint(GflashTrajectoryPoint &point, double s) const
static const double tmax[3]
TH1F * em_ssp_rho
Definition: HFGflash.h:67
const double radLength[kNumberCalorimeter]
double time
Definition: HFGflash.h:40
TH1F * em_nSpots_sd
Definition: HFGflash.h:68
TH1F * em_lateral
Definition: HFGflash.h:67
TH1F * em_long
Definition: HFGflash.h:67
Gflash3Vector & getPosition()
double p1[4]
Definition: TauolaWrapper.h:89
G4ThreeVector position
Definition: HFGflash.h:38
TH1F * em_incE
Definition: HFGflash.h:67
GflashTrajectory * theHelix
Definition: HFGflash.h:49
TH2F * em_2d
Definition: HFGflash.h:69
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double p3[4]
Definition: TauolaWrapper.h:91
tuple log
Definition: cmsBatch.py:341
void initializeTrajectory(const HepGeom::Vector3D< double > &, const HepGeom::Point3D< double > &, double q, double Field)

Member Data Documentation

TH2F* HFGflash::em_2d
private

Definition at line 69 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH2F * HFGflash::em_2d_sd
private

Definition at line 69 of file HFGflash.h.

Referenced by HFGflash().

TH1F* HFGflash::em_incE
private

Definition at line 67 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH1F * HFGflash::em_lateral
private

Definition at line 67 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH1F * HFGflash::em_lateral_sd
private

Definition at line 68 of file HFGflash.h.

Referenced by HFGflash().

TH1F * HFGflash::em_long
private

Definition at line 67 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH1F* HFGflash::em_long_sd
private

Definition at line 68 of file HFGflash.h.

Referenced by HFGflash().

TH1F * HFGflash::em_nSpots_sd
private

Definition at line 68 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH2F * HFGflash::em_ratio
private

Definition at line 69 of file HFGflash.h.

Referenced by HFGflash().

TH2F * HFGflash::em_ratio_selected
private

Definition at line 69 of file HFGflash.h.

Referenced by HFGflash().

TH1F * HFGflash::em_ssp_rho
private

Definition at line 67 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH1F * HFGflash::em_ssp_z
private

Definition at line 67 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

TH2F * HFGflash::em_ze
private

Definition at line 69 of file HFGflash.h.

Referenced by HFGflash().

TH1F * HFGflash::em_ze_ratio
private

Definition at line 68 of file HFGflash.h.

Referenced by HFGflash().

G4double HFGflash::energyScale[Gflash::kNumberCalorimeter]
private

Definition at line 62 of file HFGflash.h.

G4double HFGflash::energyToDeposit
private

Definition at line 61 of file HFGflash.h.

Gflash::CalorimeterNumber HFGflash::jCalorimeter
private

Definition at line 54 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

G4double HFGflash::lateralPar[4]
private

Definition at line 65 of file HFGflash.h.

G4double HFGflash::longEcal[Gflash::NPar]
private

Definition at line 64 of file HFGflash.h.

G4double HFGflash::longHcal[Gflash::NPar]
private

Definition at line 63 of file HFGflash.h.

G4int HFGflash::showerType
private

Definition at line 60 of file HFGflash.h.

G4double HFGflash::theBField
private

Definition at line 58 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

bool HFGflash::theFillHisto
private

Definition at line 57 of file HFGflash.h.

Referenced by gfParameterization(), and HFGflash().

G4Navigator* HFGflash::theGflashNavigator
private

Definition at line 51 of file HFGflash.h.

Referenced by gfParameterization(), HFGflash(), and ~HFGflash().

G4Step* HFGflash::theGflashStep
private

Definition at line 50 of file HFGflash.h.

Referenced by HFGflash(), and ~HFGflash().

G4TouchableHandle HFGflash::theGflashTouchableHandle
private

Definition at line 52 of file HFGflash.h.

Referenced by HFGflash().

GflashTrajectory* HFGflash::theHelix
private

Definition at line 49 of file HFGflash.h.

Referenced by gfParameterization(), HFGflash(), and ~HFGflash().

bool HFGflash::theWatcherOn
private

Definition at line 56 of file HFGflash.h.

Referenced by HFGflash().