CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFShowerParam.cc
Go to the documentation of this file.
1 // File: HFShowerParam.cc
3 // Description: Parametrized version of HF hits
5 
10 
14 
15 #include "G4VPhysicalVolume.hh"
16 #include "G4Step.hh"
17 #include "G4Track.hh"
18 #include "G4NavigationHistory.hh"
19 #include "Randomize.hh"
20 
21 #include "CLHEP/Units/GlobalPhysicalConstants.h"
22 #include "CLHEP/Units/GlobalSystemOfUnits.h"
23 
24 #include<iostream>
25 
26 //#define DebugLog
27 //#define plotDebug
28 //#define mkdebug
29 
31  edm::ParameterSet const & p) : showerLibrary(0),
32  fibre(0), gflash(0),
33  fillHisto(false) {
34  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
35  pePerGeV = m_HF.getParameter<double>("PEPerGeV");
36  trackEM = m_HF.getParameter<bool>("TrackEM");
37  bool useShowerLibrary = m_HF.getParameter<bool>("UseShowerLibrary");
38  bool useGflash = m_HF.getParameter<bool>("UseHFGflash");
39  edMin = m_HF.getParameter<double>("EminLibrary");
40  onlyLong = m_HF.getParameter<bool>("OnlyLong");
41  ref_index = m_HF.getParameter<double>("RefIndex");
42  double lambdaMean = m_HF.getParameter<double>("LambdaMean");
43  aperture = cos(asin(m_HF.getParameter<double>("Aperture")));
44  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
45  parametrizeLast = m_HF.getUntrackedParameter<bool>("ParametrizeLast",false);
46  edm::LogInfo("HFShower") << "HFShowerParam::Use of shower library is set to "
47  << useShowerLibrary << " Use of Gflash is set to "
48  << useGflash << " P.E. per GeV " << pePerGeV
49  << ", ref. index of fibre " << ref_index
50  << ", Track EM Flag " << trackEM << ", edMin "
51  << edMin << " GeV, use of Short fibre info in"
52  << " shower library set to " << !(onlyLong)
53  << ", use of parametrization for last part set to "
54  << parametrizeLast << ", Mean lambda " << lambdaMean
55  << ", aperture (cutoff) " << aperture
56  << ", Application of Fiducial Cut " << applyFidCut;
57 
58 #ifdef plotDebug
60  if (tfile.isAvailable()) {
61  fillHisto = true;
62  LogDebug("HFShower") << "HFShowerParam::Save histos in directory "
63  << "ProfileFromParam";
64  TFileDirectory showerDir = tfile->mkdir("ProfileFromParam");
65  hzvem = showerDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part);Number of PE",330,0.0,1650.0);
66  hzvhad = showerDir.make<TH1F>("hzvhad","Longitudinal Profile (Had Part);Number of PE",330,0.0,1650.0);
67  em_2d_1 = showerDir.make<TH2F>("em_2d_1","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
68  em_long_1 = showerDir.make<TH1F>("em_long_1","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
69  em_long_1_tuned = showerDir.make<TH1F>("em_long_1_tuned","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
70  em_lateral_1 = showerDir.make<TH1F>("em_lateral_1","Lateral Profile;cm;Events",100,50.0,150.0);
71  em_2d_2 = showerDir.make<TH2F>("em_2d_2","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
72  em_long_2 = showerDir.make<TH1F>("em_long_2","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
73  em_lateral_2 = showerDir.make<TH1F>("em_lateral_2","Lateral Profile;cm;Events",100,50.0,150.0);
74  em_long_gflash = showerDir.make<TH1F>("em_long_gflash","Longitudinal Profile From GFlash;cm;Number of Spots",800,800.0,1600.0);
75  em_long_sl = showerDir.make<TH1F>("em_long_sl","Longitudinal Profile From Shower Library;cm;Number of Spots",800,800.0,1600.0);
76  } else {
77  fillHisto = false;
78  edm::LogInfo("HFShower") << "HFShowerParam::No file is available for "
79  << "saving histos so the flag is set to false";
80  }
81 #endif
82 
83  G4String attribute = "ReadOutName";
84  G4String value = name;
86  DDValue ddv(attribute,value,0);
88  DDFilteredView fv(cpv);
89  fv.addFilter(filter);
90  bool dodet = fv.firstChild();
91  if (dodet) {
93  //Special Geometry parameters
94  gpar = getDDDArray("gparHF",sv);
95  edm::LogInfo("HFShower") << "HFShowerParam: " <<gpar.size() <<" gpar (cm)";
96  for (unsigned int ig=0; ig<gpar.size(); ig++)
97  edm::LogInfo("HFShower") << "HFShowerParam: gpar[" << ig << "] = "
98  << gpar[ig]/cm << " cm";
99  } else {
100  edm::LogError("HFShower") << "HFShowerParam: cannot get filtered "
101  << " view for " << attribute << " matching " << name;
102  throw cms::Exception("Unknown", "HFShowerParam") << "cannot match " << attribute
103  << " to " << name <<"\n";
104  }
105 
106  if (useShowerLibrary) showerLibrary = new HFShowerLibrary(name, cpv, p);
107  if (useGflash) gflash = new HFGflash(p);
108  fibre = new HFFibre(name, cpv, p);
109  attLMeanInv = fibre->attLength(lambdaMean);
110  edm::LogInfo("HFShower") << "att. length used for (lambda=" << lambdaMean
111  << ") = " << 1/(attLMeanInv*cm) << " cm";
112 }
113 
115  if (fibre) delete fibre;
116  if (gflash) delete gflash;
117  if (showerLibrary) delete showerLibrary;
118 }
119 
120 void HFShowerParam::initRun(G4ParticleTable * theParticleTable) {
121  emPDG = theParticleTable->FindParticle("e-")->GetPDGEncoding();
122  epPDG = theParticleTable->FindParticle("e+")->GetPDGEncoding();
123  gammaPDG = theParticleTable->FindParticle("gamma")->GetPDGEncoding();
124 #ifdef DebugLog
125  edm::LogInfo("HFShower") << "HFShowerParam: Particle code for e- = " << emPDG
126  << " for e+ = " << epPDG << " for gamma = " << gammaPDG;
127 #endif
128  if (showerLibrary) showerLibrary->initRun(theParticleTable);
129 }
130 
131 std::vector<HFShowerParam::Hit> HFShowerParam::getHits(G4Step * aStep,
132  double weight) {
133  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
134  G4Track * track = aStep->GetTrack();
135  G4ThreeVector hitPoint = preStepPoint->GetPosition();
136  G4int particleCode = track->GetDefinition()->GetPDGEncoding();
137  double zv = std::abs(hitPoint.z()) - gpar[4] - 0.5*gpar[1];
138  G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(),hitPoint.y(),zv);
139 
140  double pin = (preStepPoint->GetTotalEnergy())/GeV;
141  double zint = hitPoint.z();
142  double zz = std::abs(zint) - gpar[4];
143 
144 #ifdef DebugLog
145  edm::LogInfo("HFShower") << "HFShowerParam: getHits "
146  << track->GetDefinition()->GetParticleName()
147  << " of energy " << pin << " GeV"
148  << " Pos x,y,z = " << hitPoint.x() << ","
149  << hitPoint.y() << "," << zint << " (" << zz << ","
150  << localPoint.z() << ", "
151  << (localPoint.z()+0.5*gpar[1]) << ") Local "
152  << localPoint;
153 #endif
154  std::vector<HFShowerParam::Hit> hits;
156  hit.position = hitPoint;
157 
158  // look for other charged particles
159  bool other = false;
160  double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
161  double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).z();
162  if (hitPoint.z() < 0) dirz *= -1.;
163 #ifdef DebugLog
164  edm::LogInfo("HFShower") << "HFShowerParam: getHits Momentum "
165  <<track->GetDynamicParticle()->GetMomentumDirection()
166  << " HitPoint " << hitPoint << " dirz " << dirz;
167 #endif
168  if (particleCode != emPDG && particleCode != epPDG && particleCode != gammaPDG ) {
169  if (track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1/ref_index) &&
170  aStep->GetTotalEnergyDeposit() > 0) other = true;
171  }
172 
173  // take only e+-/gamma/or special particles
174  if (particleCode == emPDG || particleCode == epPDG ||
175  particleCode == gammaPDG || other) {
176  // Leave out the last part
177  double edep = 0.;
178  bool kill = false;
179  if ((!trackEM) && ((zz<(gpar[1]-gpar[2])) || parametrizeLast) && (!other)){
180  edep = pin;
181  kill = true;
182  } else if ((track->GetDefinition()->GetPDGCharge() != 0) &&
183  (pBeta > (1/ref_index)) && (dirz > aperture)) {
184  edep = (aStep->GetTotalEnergyDeposit())/GeV;
185  }
186  std::string path = "ShowerLibrary";
187 #ifdef DebugLog
188  edm::LogInfo("HFShower") << "HFShowerParam: getHits edep = " << edep
189  << " weight " << weight << " final " <<edep*weight
190  << ", Kill = " << kill << ", pin = " << pin
191  << ", edMin = " << edMin << " Other " << other;
192 #endif
193  edep *= weight;
194  if (edep > 0) {
195  if ((showerLibrary || gflash) && kill && pin > edMin && (!other)) {
196  if (showerLibrary) {
197  std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary->getHits(aStep,kill,weight,onlyLong);
198  for (unsigned int i=0; i<hitSL.size(); i++) {
199  bool ok = true;
200 #ifdef DebugLog
201  edm::LogInfo("HFShower") << "HFShowerParam: getHits applyFidCut = " << applyFidCut;
202 #endif
203  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
204  int npmt = HFFibreFiducial:: PMTNumber(hitSL[i].position);
205  if (npmt <= 0) ok = false;
206  }
207  if (ok) {
208  hit.position = hitSL[i].position;
209  hit.depth = hitSL[i].depth;
210  hit.time = hitSL[i].time;
211  hit.edep = 1;
212  hits.push_back(hit);
213 #ifdef plotDebug
214  if (fillHisto) {
215  double zv = std::abs(hit.position.z()) - gpar[4];
216  hzvem->Fill(zv);
217  em_long_sl->Fill(hit.position.z()/cm);
218  double sq = sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2));
219  double zp = hit.position.z()/cm;
220  if (hit.depth == 1) {
221  em_2d_1->Fill(zp, sq);
222  em_lateral_1->Fill(sq);
223  em_long_1->Fill(zp);
224  } else if (hit.depth == 2) {
225  em_2d_2->Fill(zp, sq);
226  em_lateral_2->Fill(sq);
227  em_long_2->Fill(zp);
228  }
229  }
230 #endif
231 #ifdef DebugLog
232  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth "
233  << hit.depth << " with edep " << hit.edep
234  << " Time " << hit.time;
235 #endif
236  }
237  }
238  } else { // GFlash clusters with known z
239  std::vector<HFGflash::Hit>hitSL=gflash->gfParameterization(aStep,kill, onlyLong);
240  for (unsigned int i=0; i<hitSL.size(); ++i) {
241  bool ok = true;
242  G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(),
243  hitSL[i].position.z());
244  double zv = std::abs(pe_effect.z()) - gpar[4];
245  //depth
246  int depth = 1;
247  int npmt = 0;
248  if (zv < 0. || zv > gpar[1]) {
249 #ifdef mkdebug
250  std::cout<<"-#Zcut-HFShowerParam::getHits:z="<<zv<<",m="<<gpar[1]<<std::endl;
251 #endif
252  ok = false;
253  }
254  if (ok && applyFidCut) {
255  npmt = HFFibreFiducial:: PMTNumber(pe_effect);
256 #ifdef DebugLog
257  edm::LogInfo("HFShower") << "HFShowerParam::getHits:#PMT= "
258  << npmt << ",z = " << zv;
259 #endif
260  if (npmt <= 0) {
261 #ifdef DebugLog
262  edm::LogInfo("HFShower") << "-#PMT=0 cut-HFShowerParam::"
263  << "getHits: npmt = " << npmt;
264 #endif
265  ok = false;
266  } else if (npmt > 24) { // a short fibre
267  if (zv > gpar[0]) {
268  depth = 2;
269  } else {
270 #ifdef DebugLog
271  edm::LogInfo("HFShower") << "-SHORT cut-HFShowerParam::"
272  << "getHits:zMin=" << gpar[0];
273 #endif
274  ok = false;
275  }
276  }
277 #ifdef DebugLog
278  edm::LogInfo("HFShower") << "HFShowerParam: npmt " << npmt
279  << " zv " << std::abs(pe_effect.z())
280  << ":" << gpar[4] << ":" << zv << ":"
281  << gpar[0] << " ok " << ok << " depth "
282  << depth;
283 #endif
284  } else {
285  if (G4UniformRand() > 0.5) depth = 2;
286  if (depth == 2 && zv < gpar[0]) ok = false;
287  }
288  //attenuation
289  double dist = fibre->zShift(localPoint,depth,0); // distance to PMT
290  double r1 = G4UniformRand();
291 #ifdef DebugLog
292  edm::LogInfo("HFShower") << "HFShowerParam:Distance to PMT (" <<npmt
293  << ") " << dist << ", exclusion flag "
294  << (r1 > exp(-attLMeanInv*zv));
295 #endif
296  if (r1 > exp(-attLMeanInv*dist)) ok = false;
297  if (ok) {
298  double r2 = G4UniformRand();
299 #ifdef DebugLog
300  edm::LogInfo("HFShower") << "HFShowerParam:Extra exclusion "
301  << r2 << ">" << weight << " "
302  << (r2 > weight);
303 #endif
304  if (r2 < weight) {
305  double time = fibre->tShift(localPoint,depth,0);
306 
307  hit.position = hitSL[i].position;
308  hit.depth = depth;
309  hit.time = time + hitSL[i].time;
310  hit.edep = 1;
311  hits.push_back(hit);
312 #ifdef plotDebug
313  if (fillHisto) {
314  em_long_gflash->Fill(pe_effect.z()/cm, hitSL[i].edep);
315  hzvem->Fill(zv);
316  double sq = sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2));
317  double zp = hit.position.z()/cm;
318  if (hit.depth == 1) {
319  em_2d_1->Fill(zp, sq);
320  em_lateral_1->Fill(s);
321  em_long_1->Fill(zp);
322  } else if (hit.depth == 2) {
323  em_2d_2->Fill(zp, sq);
324  em_lateral_2->Fill(sq);
325  em_long_2->Fill(zp);
326  }
327  }
328 #endif
329 #ifdef DebugLog
330  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth "
331  << hit.depth << " with edep "
332  << hit.edep << " Time " << hit.time;
333 #endif
334  }
335  }
336  }
337  }
338  } else {
339  path = "Rest";
340  edep *= pePerGeV;
341  double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
342  double time = fibre->tShift(localPoint,1,0); // remaining part
343  bool ok = true;
344  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
345  int npmt = HFFibreFiducial:: PMTNumber(hitPoint);
346  if (npmt <= 0) ok = false;
347  }
348 #ifdef DebugLog
349  edm::LogInfo("HFShower") << "HFShowerParam: getHits hitPoint " << hitPoint << " flag " << ok;
350 #endif
351  if (ok) {
352  hit.depth = 1;
353  hit.time = tSlice+time;
354  hit.edep = edep;
355  hits.push_back(hit);
356 #ifdef DebugLog
357  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth 1 with edep "
358  << edep << " Time " << tSlice << ":" << time
359  << ":" << hit.time;
360 #endif
361 #ifdef plotDebug
362  double zv = std::abs(hitPoint.z()) - gpar[4];
363  if (fillHisto) {
364  hzvhad->Fill(zv);
365  }
366 #endif
367  if (zz >= gpar[0]) {
368  time = fibre->tShift(localPoint,2,0);
369  hit.depth = 2;
370  hit.time = tSlice+time;
371  hits.push_back(hit);
372 #ifdef DebugLog
373  edm::LogInfo("HFShower") <<"HFShowerParam: Hit at depth 2 with edep "
374  << edep << " Time " << tSlice << ":" << time
375  << hit.time;
376 #endif
377 #ifdef plotDebug
378  if (fillHisto) {
379  hzvhad->Fill(zv);
380  }
381 #endif
382  }
383  }
384  }
385 #ifdef DebugLog
386  for (unsigned int ii=0; ii<hits.size(); ++ii) {
387  double zv = std::abs(hits[ii].position.z());
388  if (zv > 12790) edm::LogInfo("HFShower")<< "HFShowerParam: Abnormal hit along "
389  << path << " in "
390  << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
391  << " at " << hits[ii].position << " zz "
392  << zv << " Edep " << edep << " due to "
393  <<track->GetDefinition()->GetParticleName()
394  << " time " << hit.time;
395  }
396 #endif
397  if (kill) {
398  track->SetTrackStatus(fStopAndKill);
399  G4TrackVector tv = *(aStep->GetSecondary());
400  for (unsigned int kk=0; kk<tv.size(); ++kk) {
401  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
402  tv[kk]->SetTrackStatus(fStopAndKill);
403  }
404  }
405 #ifdef DebugLog
406  edm::LogInfo("HFShower") << "HFShowerParam: getHits kill (" << kill
407  << ") track " << track->GetTrackID()
408  << " at " << hitPoint
409  << " and deposit " << edep << " " << hits.size()
410  << " times" << " ZZ " << zz << " " << gpar[0];
411 #endif
412  }
413  }
414  return hits;
415 }
416 
417 std::vector<double> HFShowerParam::getDDDArray(const std::string & str,
418  const DDsvalues_type & sv)
419 {
420 #ifdef DebugLog
421  LogDebug("HFShower") << "HFShowerParam:getDDDArray called for " << str;
422 #endif
423  DDValue value(str);
424  if (DDfetch(&sv,value))
425  {
426 #ifdef DebugLog
427  LogDebug("HFShower") << value;
428 #endif
429  const std::vector<double> & fvec = value.doubles();
430  int nval = fvec.size();
431  if (nval < 2) {
432  edm::LogError("HFShower") << "HFShowerParam : # of " << str
433  << " bins " << nval << " < 2 ==> illegal";
434  throw cms::Exception("Unknown", "HFShowerParam") << "nval < 2 for array "
435  << str << "\n";
436  }
437  return fvec;
438  } else {
439  edm::LogError("HFShower") << "HFShowerParam : cannot get array " << str;
440  throw cms::Exception("Unknown", "HFShowerParam") << "cannot get array "
441  << str << "\n";
442  }
443 }
#define LogDebug(id)
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:132
std::vector< Hit > gfParameterization(G4Step *aStep, bool &ok, bool onlyLong=false)
Definition: HFGflash.cc:91
T getParameter(std::string const &) const
double attLMeanInv
Definition: HFShowerParam.h:54
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:134
bool parametrizeLast
Definition: HFShowerParam.h:55
void addFilter(const DDFilter &, log_op op=AND)
virtual ~HFShowerParam()
std::vector< Hit > getHits(G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
double attLength(double lambda)
Definition: HFFibre.cc:114
std::vector< double > gpar
Definition: HFShowerParam.h:57
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:144
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:51
int ii
Definition: cuy.py:588
HFShowerParam(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
type of data representation of DDCompactView
Definition: DDCompactView.h:77
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
TH1F * em_lateral_2
Definition: HFShowerParam.h:59
HFFibre * fibre
Definition: HFShowerParam.h:52
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
static int PMTNumber(const G4ThreeVector &pe_effect)
tuple path
else: Piece not in the list, fine.
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:19
bool isAvailable() const
Definition: Service.h:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T * make(const Args &...args) const
make new ROOT object
HFGflash * gflash
Definition: HFShowerParam.h:53
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
G4ThreeVector position
Definition: HFShowerParam.h:38
static const G4LogicalVolume * GetVolume(const std::string &name)
double ref_index
Definition: HFShowerParam.h:54
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
std::vector< Hit > getHits(G4Step *aStep, double weight)
DDsvalues_type mergedSpecifics() const
void initRun(G4ParticleTable *theParticleTable)
TH1F * em_long_gflash
Definition: HFShowerParam.h:60
TH1F * em_long_sl
Definition: HFShowerParam.h:61
TH1F * em_long_1
Definition: HFShowerParam.h:59
void initRun(G4ParticleTable *)
TH1F * em_long_1_tuned
Definition: HFShowerParam.h:60
bool firstChild()
set the current node to the first child ...
TH1F * em_lateral_1
Definition: HFShowerParam.h:59
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
tuple cout
Definition: gather_cfg.py:121
TH1F * em_long_2
Definition: HFShowerParam.h:59
volatile std::atomic< bool > shutdown_flag false
int weight
Definition: histoStyle.py:50
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37