CMS 3D CMS Logo

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

#include <ECalSD.h>

Inheritance diagram for ECalSD:
CaloSD SensitiveCaloDetector Observer< const BeginOfRun * > Observer< const BeginOfEvent * > Observer< const BeginOfTrack * > Observer< const EndOfTrack * > Observer< const EndOfEvent * > SensitiveDetector

Public Member Functions

 ECalSD (G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
virtual uint16_t getDepth (G4Step *)
 
virtual double getEnergyDeposit (G4Step *)
 
virtual uint16_t getRadiationLength (G4Step *)
 
virtual int getTrackID (G4Track *)
 
virtual uint32_t setDetUnitId (G4Step *)
 
void setNumberingScheme (EcalNumberingScheme *)
 
virtual ~ECalSD ()
 
- Public Member Functions inherited from CaloSD
 CaloSD (G4String aSDname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, int tSlice=1, bool ignoreTkID=false)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *eventHC)
 
void fillHits (edm::PCaloHitContainer &, std::string n)
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
virtual void PrintAll ()
 
virtual bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory)
 
virtual bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *)
 
virtual ~CaloSD ()
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (G4ThreeVector point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 
- Public Member Functions inherited from Observer< const BeginOfRun * >
 Observer ()
 
void slotForUpdate (const BeginOfRun *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack * >
 Observer ()
 
void slotForUpdate (const BeginOfTrack *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack * >
 Observer ()
 
void slotForUpdate (const EndOfTrack *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent * >
 Observer ()
 
void slotForUpdate (const EndOfEvent *iT)
 
virtual ~Observer ()
 

Private Member Functions

double crystalLength (G4LogicalVolume *)
 
double curve_LY (G4Step *)
 
void getBaseNumber (const G4Step *)
 
double getBirkL3 (G4Step *)
 
std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)
 
std::vector< std::string > getStringArray (const std::string &, const DDsvalues_type &)
 
void initMap (G4String, const DDCompactView &)
 

Private Attributes

double birk1
 
double birk2
 
double birk3
 
double birkCut
 
double birkSlope
 
std::string crystalMat
 
std::string depth1Name
 
std::string depth2Name
 
std::vector< G4LogicalVolume * > noWeight
 
EcalNumberingSchemenumberingScheme
 
double slopeLY
 
bool storeRL
 
bool storeTrack
 
EcalBaseNumber theBaseNumber
 
bool useBirk
 
bool useBirkL3
 
std::vector< G4LogicalVolume * > useDepth1
 
std::vector< G4LogicalVolume * > useDepth2
 
bool useWeight
 
std::map< G4LogicalVolume
*, double > 
xtalLMap
 

Additional Inherited Members

- Public Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Member Functions inherited from CaloSD
G4bool checkHit ()
 
virtual void clearHits ()
 
CaloG4HitcreateNewHit ()
 
virtual bool filterHit (CaloG4Hit *, double)
 
double getAttenuation (G4Step *aStep, double birk1, double birk2, double birk3)
 
double getResponseWt (G4Track *)
 
virtual G4bool getStepInfo (G4Step *aStep)
 
G4bool hitExists ()
 
virtual void initRun ()
 
void resetForNewPrimary (G4ThreeVector, double)
 
G4ThreeVector setToLocal (G4ThreeVector, const G4VTouchable *)
 
virtual void update (const BeginOfRun *)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const BeginOfEvent *)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const BeginOfTrack *trk)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const EndOfTrack *trk)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const ::EndOfEvent *)
 
void updateHit (CaloG4Hit *)
 
- Protected Member Functions inherited from Observer< const EndOfEvent * >
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 
- Protected Attributes inherited from CaloSD
int checkHits
 
double correctT
 
bool corrTOFBeam
 
CaloG4HitcurrentHit
 
CaloHitID currentID
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double eminHitD
 
G4int emPDG
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
G4int epPDG
 
bool forceSave
 
G4int gammaPDG
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
const SimTrackManagerm_trackManager
 
G4ThreeVector posGlobal
 
G4StepPoint * preStepPoint
 
CaloHitID previousID
 
int primIDSaved
 
bool runInit
 
bool suppressHeavy
 
G4Track * theTrack
 
double tmaxHit
 
bool useMap
 

Detailed Description

Definition at line 24 of file ECalSD.h.

Constructor & Destructor Documentation

ECalSD::ECalSD ( G4String  name,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 31 of file ECalSD.cc.

References DDFilteredView::addFilter(), birk1, birk2, birk3, birkCut, birkSlope, crystalMat, depth1Name, depth2Name, DDSpecificsFilter::equals, align_tpl::filter, DDFilteredView::firstChild(), g, getDDDArray(), edm::ParameterSet::getParameter(), getStringArray(), edm::ParameterSet::getUntrackedParameter(), initMap(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, DDFilteredView::mergedSpecifics(), DDSpecificsFilter::setCriteria(), setNumberingScheme(), slopeLY, storeRL, storeTrack, CaloSD::suppressHeavy, useBirk, useBirkL3, and useWeight.

33  :
34  CaloSD(name, cpv, clg, p, manager,
35  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<int>("TimeSliceUnit"),
36  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<bool>("IgnoreTrackID")),
37  numberingScheme(0) {
38 
39  // static SimpleConfigurable<bool> on1(false, "ECalSD:UseBirkLaw");
40  // static SimpleConfigurable<double> bk1(0.00463,"ECalSD:BirkC1");
41  // static SimpleConfigurable<double> bk2(-0.03, "ECalSD:BirkC2");
42  // static SimpleConfigurable<double> bk3(1.0, "ECalSD:BirkC3");
43  // Values from NIM A484 (2002) 239-244: as implemented in Geant3
44  // useBirk = on1.value();
45  // birk1 = bk1.value()*(g/(MeV*cm2));
46  // birk2 = bk2.value()*(g/(MeV*cm2))*(g/(MeV*cm2));
47 
48  edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
49  useBirk = m_EC.getParameter<bool>("UseBirkLaw");
50  useBirkL3 = m_EC.getParameter<bool>("BirkL3Parametrization");
51  birk1 = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
52  birk2 = m_EC.getParameter<double>("BirkC2");
53  birk3 = m_EC.getParameter<double>("BirkC3");
54  birkSlope = m_EC.getParameter<double>("BirkSlope");
55  birkCut = m_EC.getParameter<double>("BirkCut");
56  slopeLY = m_EC.getParameter<double>("SlopeLightYield");
57  storeTrack = m_EC.getParameter<bool>("StoreSecondary");
58  crystalMat = m_EC.getUntrackedParameter<std::string>("XtalMat","E_PbWO4");
59  bool isItTB = m_EC.getUntrackedParameter<bool>("TestBeam", false);
60  bool nullNS = m_EC.getUntrackedParameter<bool>("NullNumbering", false);
61  storeRL = m_EC.getUntrackedParameter<bool>("StoreRadLength", false);
62 
63  //Material list for HB/HE/HO sensitive detectors
64  std::string attribute = "ReadOutName";
66  DDValue ddv(attribute,name,0);
68  DDFilteredView fv(cpv);
69  fv.addFilter(filter);
70  fv.firstChild();
71  DDsvalues_type sv(fv.mergedSpecifics());
72  // Use of Weight
73  useWeight= true;
74  std::vector<double> tempD = getDDDArray("EnergyWeight",sv);
75  if (tempD.size() > 0) { if (tempD[0] < 0.1) useWeight = false; }
76  std::vector<std::string> tempS = getStringArray("Depth1Name",sv);
77  if (tempS.size() > 0) depth1Name = tempS[0];
78  else depth1Name = " ";
79  tempS = getStringArray("Depth2Name",sv);
80  if (tempS.size() > 0) depth2Name = tempS[0];
81  else depth2Name = " ";
82 
83  EcalNumberingScheme* scheme=0;
84  if (nullNS) scheme = 0;
85  else if (name == "EcalHitsEB") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
86  else if (name == "EcalHitsEE") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
87  else if (name == "EcalHitsES") {
88  if (isItTB) scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
89  else scheme = dynamic_cast<EcalNumberingScheme*>(new EcalPreshowerNumberingScheme());
90  useWeight = false;
91  } else {edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported\n";}
92 
93  if (scheme) setNumberingScheme(scheme);
94 #ifdef DebugLog
95  LogDebug("EcalSim") << "Constructing a ECalSD with name " << GetName();
96 #endif
97  if (useWeight) {
98  edm::LogInfo("EcalSim") << "ECalSD:: Use of Birks law is set to "
99  << useBirk << " with three constants kB = "
100  << birk1 << ", C1 = " << birk2 << ", C2 = "
101  << birk3 <<"\n Use of L3 parametrization "
102  << useBirkL3 << " with slope " << birkSlope
103  << " and cut off " << birkCut << "\n"
104  << " Slope for Light yield is set to "
105  << slopeLY;
106  } else {
107  edm::LogInfo("EcalSim") << "ECalSD:: energy deposit is not corrected "
108  << " by Birk or light yield curve";
109  }
110 
111  edm::LogInfo("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy
112  << " protons below " << kmaxProton << " MeV,"
113  << " neutrons below " << kmaxNeutron << " MeV and"
114  << " ions below " << kmaxIon << " MeV\n"
115  << " Depth1 Name = " << depth1Name
116  << " and Depth2 Name = " << depth2Name;
117 
118  if (useWeight) initMap(name,cpv);
119 
120 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool useBirkL3
Definition: ECalSD.h:52
double kmaxNeutron
Definition: CaloSD.h:132
std::vector< std::string > getStringArray(const std::string &, const DDsvalues_type &)
Definition: ECalSD.cc:493
double birkSlope
Definition: ECalSD.h:53
double slopeLY
Definition: ECalSD.h:54
double birk1
Definition: ECalSD.h:53
bool useWeight
Definition: ECalSD.h:51
double kmaxProton
Definition: CaloSD.h:132
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
bool storeRL
Definition: ECalSD.h:51
CaloSD(G4String aSDname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, int tSlice=1, bool ignoreTkID=false)
Definition: CaloSD.cc:21
std::string depth1Name
Definition: ECalSD.h:55
bool storeTrack
Definition: ECalSD.h:51
double kmaxIon
Definition: CaloSD.h:132
bool suppressHeavy
Definition: CaloSD.h:130
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
void initMap(G4String, const DDCompactView &)
Definition: ECalSD.cc:245
std::string depth2Name
Definition: ECalSD.h:55
bool useBirk
Definition: ECalSD.h:52
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: ECalSD.cc:474
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
double birk2
Definition: ECalSD.h:53
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:50
void setNumberingScheme(EcalNumberingScheme *)
Definition: ECalSD.cc:236
double birkCut
Definition: ECalSD.h:53
std::string crystalMat
Definition: ECalSD.h:55
double birk3
Definition: ECalSD.h:53
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
ECalSD::~ECalSD ( )
virtual

Definition at line 122 of file ECalSD.cc.

References numberingScheme.

122  {
123  if (numberingScheme) delete numberingScheme;
124 }
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:50

Member Function Documentation

double ECalSD::crystalLength ( G4LogicalVolume *  lv)
private

Definition at line 420 of file ECalSD.cc.

References xtalLMap.

Referenced by curve_LY(), and getRadiationLength().

420  {
421 
422  double length= 230.;
423  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.find(lv);
424  if (ite != xtalLMap.end()) length = ite->second;
425  return length;
426 }
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:56
double ECalSD::curve_LY ( G4Step *  aStep)
private

Definition at line 390 of file ECalSD.cc.

References crystalLength(), LogDebug, CaloSD::setToLocal(), slopeLY, and CommonMethods::weight().

Referenced by getEnergyDeposit().

390  {
391 
392  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
393  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
394 
395  double weight = 1.;
396  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
397  stepPoint->GetTouchable());
398  double crlength = crystalLength(lv);
399  double dapd = 0.5 * crlength - localPoint.z();
400  if (dapd >= -0.1 || dapd <= crlength+0.1) {
401  if (dapd <= 100.)
402  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
403  } else {
404  edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
405  << "to APD " << dapd << " crlength = "
406  << crlength <<" crystal name = " <<lv->GetName()
407  << " z of localPoint = " << localPoint.z()
408  << " take weight = " << weight;
409  }
410 #ifdef DebugLog
411  LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd
412  << " crlength = " << crlength
413  << " crystal name = " << lv->GetName()
414  << " z of localPoint = " << localPoint.z()
415  << " take weight = " << weight;
416 #endif
417  return weight;
418 }
G4ThreeVector setToLocal(G4ThreeVector, const G4VTouchable *)
Definition: CaloSD.cc:310
#define LogDebug(id)
double slopeLY
Definition: ECalSD.h:54
double crystalLength(G4LogicalVolume *)
Definition: ECalSD.cc:420
void ECalSD::getBaseNumber ( const G4Step *  aStep)
private

Definition at line 428 of file ECalSD.cc.

References EcalBaseNumber::addLevel(), EcalBaseNumber::getCapacity(), LogDebug, EcalBaseNumber::reset(), EcalBaseNumber::setSize(), and theBaseNumber.

Referenced by setDetUnitId().

428  {
429 
431  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
432  int theSize = touch->GetHistoryDepth()+1;
433  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
434  //Get name and copy numbers
435  if ( theSize > 1 ) {
436  for (int ii = 0; ii < theSize ; ii++) {
437  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
438 #ifdef DebugLog
439  LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
440  << ": " << touch->GetVolume(ii)->GetName() << "["
441  << touch->GetReplicaNumber(ii) << "]";
442 #endif
443  }
444  }
445 
446 }
#define LogDebug(id)
void addLevel(const std::string &name, const int &copyNumber)
EcalBaseNumber theBaseNumber
Definition: ECalSD.h:58
void setSize(const int &size)
double ECalSD::getBirkL3 ( G4Step *  aStep)
private

Definition at line 448 of file ECalSD.cc.

References birk1, birkCut, birkSlope, DeDxDiscriminatorTools::charge(), funct::log(), LogDebug, and CommonMethods::weight().

Referenced by getEnergyDeposit().

448  {
449 
450  double weight = 1.;
451  double charge = aStep->GetPreStepPoint()->GetCharge();
452 
453  if (charge != 0. && aStep->GetStepLength() > 0) {
454  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
455  double density = mat->GetDensity();
456  double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
457  double rkb = birk1/density;
458  if (dedx > 0) {
459  weight = 1. - birkSlope*log(rkb*dedx);
460  if (weight < birkCut) weight = birkCut;
461  else if (weight > 1.) weight = 1.;
462  }
463 #ifdef DebugLog
464  LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
465  << " Charge " << charge << " dE/dx " << dedx
466  << " Birk Const " << rkb << " Weight = " << weight
467  << " dE " << aStep->GetTotalEnergyDeposit();
468 #endif
469  }
470  return weight;
471 
472 }
#define LogDebug(id)
double birkSlope
Definition: ECalSD.h:53
double birk1
Definition: ECalSD.h:53
double charge(const std::vector< uint8_t > &Ampls)
Log< T >::type log(const T &t)
Definition: Log.h:22
double birkCut
Definition: ECalSD.h:53
std::vector< double > ECalSD::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 474 of file ECalSD.cc.

References DDfetch(), DDValue::doubles(), LogDebug, and relativeConstraints::value.

Referenced by ECalSD().

475  {
476 
477 #ifdef DebugLog
478  LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
479 #endif
480  DDValue value(str);
481  if (DDfetch(&sv,value)) {
482 #ifdef DebugLog
483  LogDebug("EcalSim") << value;
484 #endif
485  const std::vector<double> & fvec = value.doubles();
486  return fvec;
487  } else {
488  std::vector<double> fvec;
489  return fvec;
490  }
491 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
uint16_t ECalSD::getDepth ( G4Step *  aStep)
virtual

Reimplemented from CaloSD.

Definition at line 196 of file ECalSD.cc.

References prof2calltree::count, getRadiationLength(), LogDebug, runTheMatrix::ret, storeRL, useDepth1, and useDepth2.

196  {
197  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
198  uint16_t ret = 0;
199  if (std::count(useDepth1.begin(),useDepth1.end(),lv) != 0) ret = 1;
200  else if (std::count(useDepth2.begin(),useDepth2.end(),lv) != 0) ret = 2;
201  else if (storeRL) ret = getRadiationLength(aStep);
202 #ifdef DebugLog
203  LogDebug("EcalSim") << "Volume " << lv->GetName() << " Depth " << ret;
204 #endif
205  return ret;
206 }
#define LogDebug(id)
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:57
bool storeRL
Definition: ECalSD.h:51
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:57
virtual uint16_t getRadiationLength(G4Step *)
Definition: ECalSD.cc:208
double ECalSD::getEnergyDeposit ( G4Step *  aStep)
virtual

Reimplemented from CaloSD.

Definition at line 126 of file ECalSD.cc.

References birk1, birk2, birk3, prof2calltree::count, curve_LY(), CaloSD::getAttenuation(), getBirkL3(), CaloSD::getResponseWt(), TrackInformation::isPrimary(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, noWeight, NULL, CaloSD::preStepPoint, CaloSD::suppressHeavy, CaloSD::theTrack, useBirk, useBirkL3, useWeight, and CommonMethods::weight().

126  {
127 
128  if (aStep == NULL) {
129  return 0;
130  } else {
131  preStepPoint = aStep->GetPreStepPoint();
132  G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
133 
134  // take into account light collection curve for crystals
135  double weight = 1.;
136  if (suppressHeavy) {
137  G4Track* theTrack = aStep->GetTrack();
138  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
139  if (trkInfo) {
140  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
141  if (!(trkInfo->isPrimary())) { // Only secondary particles
142  double ke = theTrack->GetKineticEnergy()/MeV;
143  if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
144  ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
145  if ((pdg == 2212) && (ke < kmaxProton)) weight = 0;
146  if ((pdg == 2112) && (ke < kmaxNeutron)) weight = 0;
147 #ifdef DebugLog
148  if (weight == 0)
149  LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
150  << " Type " << theTrack->GetDefinition()->GetParticleName()
151  << " Kinetic Energy " << ke << " MeV";
152 #endif
153  }
154  }
155  }
156  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
157  if (useWeight && std::count(noWeight.begin(),noWeight.end(),lv) == 0) {
158  weight *= curve_LY(aStep);
159  if (useBirk) {
160  if (useBirkL3) weight *= getBirkL3(aStep);
161  else weight *= getAttenuation(aStep, birk1, birk2, birk3);
162  }
163  }
164  double wt1 = getResponseWt(theTrack);
165  double edep = aStep->GetTotalEnergyDeposit() * weight * wt1;
166 #ifdef DebugLog
167  LogDebug("EcalSim") << "ECalSD:: " << nameVolume
168  <<" Light Collection Efficiency " <<weight << ":" <<wt1
169  << " Weighted Energy Deposit " << edep/MeV << " MeV";
170 #endif
171  return edep;
172  }
173 }
#define LogDebug(id)
bool useBirkL3
Definition: ECalSD.h:52
double kmaxNeutron
Definition: CaloSD.h:132
double birk1
Definition: ECalSD.h:53
bool useWeight
Definition: ECalSD.h:51
#define NULL
Definition: scimark2.h:8
double kmaxProton
Definition: CaloSD.h:132
std::vector< G4LogicalVolume * > noWeight
Definition: ECalSD.h:57
double kmaxIon
Definition: CaloSD.h:132
bool suppressHeavy
Definition: CaloSD.h:130
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:480
G4Track * theTrack
Definition: CaloSD.h:116
double curve_LY(G4Step *)
Definition: ECalSD.cc:390
G4StepPoint * preStepPoint
Definition: CaloSD.h:118
bool useBirk
Definition: ECalSD.h:52
double getBirkL3(G4Step *)
Definition: ECalSD.cc:448
bool isPrimary() const
double getResponseWt(G4Track *)
Definition: CaloSD.cc:629
double birk2
Definition: ECalSD.h:53
double birk3
Definition: ECalSD.h:53
uint16_t ECalSD::getRadiationLength ( G4Step *  aStep)
virtual

Definition at line 208 of file ECalSD.cc.

References crystalLength(), NULL, CaloSD::setToLocal(), and useWeight.

Referenced by getDepth().

208  {
209 
210  uint16_t thisX0 = 0;
211  if (aStep != NULL) {
212  G4StepPoint* hitPoint = aStep->GetPreStepPoint();
213  G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
214 
215  if (useWeight) {
216  G4ThreeVector localPoint = setToLocal(hitPoint->GetPosition(),
217  hitPoint->GetTouchable());
218  double crlength = crystalLength(lv);
219  double radl = hitPoint->GetMaterial()->GetRadlen();
220  double detz = (float)(0.5*crlength + localPoint.z());
221  thisX0 = (uint16_t)floor(detz/radl);
222  }
223  }
224  return thisX0;
225 }
G4ThreeVector setToLocal(G4ThreeVector, const G4VTouchable *)
Definition: CaloSD.cc:310
bool useWeight
Definition: ECalSD.h:51
#define NULL
Definition: scimark2.h:8
double crystalLength(G4LogicalVolume *)
Definition: ECalSD.cc:420
std::vector< std::string > ECalSD::getStringArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 493 of file ECalSD.cc.

References DDfetch(), LogDebug, DDValue::strings(), and relativeConstraints::value.

Referenced by ECalSD().

494  {
495 
496 #ifdef DebugLog
497  LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
498 #endif
499  DDValue value(str);
500  if (DDfetch(&sv,value)) {
501 #ifdef DebugLog
502  LogDebug("EcalSim") << value;
503 #endif
504  const std::vector<std::string> & fvec = value.strings();
505  return fvec;
506  } else {
507  std::vector<std::string> fvec;
508  return fvec;
509  }
510 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
int ECalSD::getTrackID ( G4Track *  aTrack)
virtual

Reimplemented from CaloSD.

Definition at line 175 of file ECalSD.cc.

References prof2calltree::count, CaloSD::forceSave, CaloSD::getTrackID(), CaloSD::preStepPoint, storeTrack, useDepth1, and useDepth2.

175  {
176 
177  int primaryID(0);
178  bool flag(false);
179  if (storeTrack) {
180  G4LogicalVolume* lv = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
181  if (std::count(useDepth1.begin(),useDepth1.end(),lv) != 0) {
182  flag = true;
183  } else if (std::count(useDepth2.begin(),useDepth2.end(),lv) != 0) {
184  flag = true;
185  }
186  }
187  if (flag) {
188  forceSave = true;
189  primaryID = aTrack->GetTrackID();
190  } else {
191  primaryID = CaloSD::getTrackID(aTrack);
192  }
193  return primaryID;
194 }
long int flag
Definition: mlp_lapack.h:47
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:57
virtual int getTrackID(G4Track *)
Definition: CaloSD.cc:596
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:57
bool forceSave
Definition: CaloSD.h:135
bool storeTrack
Definition: ECalSD.h:51
G4StepPoint * preStepPoint
Definition: CaloSD.h:118
void ECalSD::initMap ( G4String  sd,
const DDCompactView cpv 
)
private

Definition at line 245 of file ECalSD.cc.

References DDFilteredView::addFilter(), prof2calltree::count, crystalMat, ddtrap, depth1Name, depth2Name, DDSpecificsFilter::equals, align_tpl::filter, DDFilteredView::firstChild(), i, LogDebug, DDFilteredView::logicalPart(), DDLogicalPart::material(), DDName::name(), SensitiveDetector::name, DDBase< N, C >::name(), DDFilteredView::next(), noWeight, DDSolid::parameters(), DDSpecificsFilter::setCriteria(), DDSolid::shape(), DDLogicalPart::solid(), useDepth1, useDepth2, and xtalLMap.

Referenced by ECalSD().

245  {
246 
247  G4String attribute = "ReadOutName";
249  DDValue ddv(attribute,sd,0);
251  DDFilteredView fv(cpv);
252  fv.addFilter(filter);
253  fv.firstChild();
254 
255  std::vector<G4LogicalVolume*> lvused;
256  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
257  std::vector<G4LogicalVolume *>::const_iterator lvcite;
258  bool dodet=true;
259  std::string lvnamx, lvnamy, lvname;
260  while (dodet) {
261  const std::string &matname = fv.logicalPart().material().name().name();
262  lvname = fv.logicalPart().name().name();
263  G4LogicalVolume* lv=0;
264  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
265  if (!strcmp((*lvcite)->GetName().c_str(), lvname.c_str())) {
266  lv = (*lvcite);
267  break;
268  }
269  }
270  if (depth1Name != " ") {
271  if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
272  if (std::count(useDepth1.begin(),useDepth1.end(),lv) == 0) {
273  useDepth1.push_back(lv);
274 #ifdef DebugLog
275  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
276  <<" in Depth 1 volume list";
277 #endif
278  }
279  lvnamx = lvname + "_refl";
280  G4LogicalVolume* lvr = 0;
281  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
282  if (!strcmp((*lvcite)->GetName().c_str(), lvnamx.c_str())) {
283  lvr = (*lvcite);
284  break;
285  }
286  }
287  if (lvr != 0 && std::count(useDepth1.begin(),useDepth1.end(),lvr)==0) {
288  useDepth1.push_back(lvr);
289 #ifdef DebugLog
290  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvnamx
291  <<" in Depth 1 volume list";
292 #endif
293  }
294  }
295  }
296  if (depth2Name != " ") {
297  if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
298  if (std::count(useDepth2.begin(),useDepth2.end(),lv) == 0) {
299  useDepth2.push_back(lv);
300 #ifdef DebugLog
301  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
302  <<" in Depth 2 volume list";
303 #endif
304  }
305  lvnamy = lvname + "_refl";
306  G4LogicalVolume* lvr = 0;
307  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
308  if (!strcmp((*lvcite)->GetName().c_str(), lvnamy.c_str())) {
309  lvr = (*lvcite);
310  break;
311  }
312  }
313  if (lvr != 0 && std::count(useDepth2.begin(),useDepth2.end(),lvr)==0) {
314  useDepth2.push_back(lvr);
315 #ifdef DebugLog
316  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvnamy
317  <<" in Depth 2 volume list";
318 #endif
319  }
320  }
321  }
322  if (lv != 0) {
323  if (strcmp(crystalMat.c_str(), matname.c_str()) == 0) {
324  if (std::count(lvused.begin(),lvused.end(),lv) == 0) {
325  lvused.push_back(lv);
326  const DDSolid & sol = fv.logicalPart().solid();
327  const std::vector<double> & paras = sol.parameters();
328 #ifdef DebugLog
329  LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid "
330  << lvname << " Shape " << sol.shape()
331  << " Parameter 0 = "<< paras[0]
332  << " Logical Volume " << lv;
333 #endif
334  if (sol.shape() == ddtrap) {
335  double dz = 2*paras[0];
336  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
337  lvname += "_refl";
338  lv = 0;
339  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
340  if (!strcmp((*lvcite)->GetName().c_str(), lvname.c_str())) {
341  lv = (*lvcite);
342  break;
343  }
344  }
345  if (lv != 0)
346  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
347  }
348  }
349  } else {
350  if (std::count(noWeight.begin(),noWeight.end(),lv) == 0) {
351  noWeight.push_back(lv);
352 #ifdef DebugLog
353  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
354  << " Material " << matname <<" in noWeight list";
355 #endif
356  }
357  lvname += "_refl";
358  lv = 0;
359  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
360  if (!strcmp((*lvcite)->GetName().c_str(), lvname.c_str())) {
361  lv = (*lvcite);
362  break;
363  }
364  }
365  if (lv != 0 && std::count(noWeight.begin(),noWeight.end(),lv) == 0) {
366  noWeight.push_back(lv);
367 #ifdef DebugLog
368  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
369  << " Material " << matname <<" in noWeight list";
370 #endif
371  }
372  }
373  }
374  dodet = fv.next();
375  }
376 #ifdef DebugLog
377  LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = "
378  << sd << ":";
379  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
380  int i=0;
381  for (; ite != xtalLMap.end(); ite++, i++) {
382  G4String name = "Unknown";
383  if (ite->first != 0) name = (ite->first)->GetName();
384  LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
385  << " L = " << ite->second;
386  }
387 #endif
388 }
#define LogDebug(id)
const std::vector< double > & parameters() const
Don&#39;t use (only meant to be used by DDbox(), DDtub(), ...)
Definition: DDSolid.cc:153
int i
Definition: DBlmapReader.cc:9
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:57
DDSolidShape shape() const
The type of the solid.
Definition: DDSolid.cc:147
std::vector< G4LogicalVolume * > noWeight
Definition: ECalSD.h:57
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:57
A DDSolid represents the shape of a part.
Definition: DDSolid.h:42
std::string depth1Name
Definition: ECalSD.h:55
std::string depth2Name
Definition: ECalSD.h:55
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:56
std::string crystalMat
Definition: ECalSD.h:55
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
uint32_t ECalSD::setDetUnitId ( G4Step *  aStep)
virtual

Implements CaloSD.

Definition at line 227 of file ECalSD.cc.

References getBaseNumber(), EcalNumberingScheme::getUnitID(), numberingScheme, and theBaseNumber.

227  {
228  if (numberingScheme == 0) {
229  return EBDetId(1,1)();
230  } else {
231  getBaseNumber(aStep);
233  }
234 }
void getBaseNumber(const G4Step *)
Definition: ECalSD.cc:428
EcalBaseNumber theBaseNumber
Definition: ECalSD.h:58
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:50
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
void ECalSD::setNumberingScheme ( EcalNumberingScheme scheme)

Definition at line 236 of file ECalSD.cc.

References numberingScheme.

Referenced by ECalSD(), and HcalTB04Analysis::update().

236  {
237  if (scheme != 0) {
238  edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for "
239  << GetName() << "\n";
240  if (numberingScheme) delete numberingScheme;
241  numberingScheme = scheme;
242  }
243 }
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:50

Member Data Documentation

double ECalSD::birk1
private

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), getBirkL3(), and getEnergyDeposit().

double ECalSD::birk2
private

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birk3
private

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birkCut
private

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

double ECalSD::birkSlope
private

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

std::string ECalSD::crystalMat
private

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth1Name
private

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth2Name
private

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::vector<G4LogicalVolume*> ECalSD::noWeight
private

Definition at line 57 of file ECalSD.h.

Referenced by getEnergyDeposit(), and initMap().

EcalNumberingScheme* ECalSD::numberingScheme
private

Definition at line 50 of file ECalSD.h.

Referenced by setDetUnitId(), setNumberingScheme(), and ~ECalSD().

double ECalSD::slopeLY
private

Definition at line 54 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

bool ECalSD::storeRL
private

Definition at line 51 of file ECalSD.h.

Referenced by ECalSD(), and getDepth().

bool ECalSD::storeTrack
private

Definition at line 51 of file ECalSD.h.

Referenced by ECalSD(), and getTrackID().

EcalBaseNumber ECalSD::theBaseNumber
private

Definition at line 58 of file ECalSD.h.

Referenced by getBaseNumber(), and setDetUnitId().

bool ECalSD::useBirk
private

Definition at line 52 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

bool ECalSD::useBirkL3
private

Definition at line 52 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

std::vector<G4LogicalVolume*> ECalSD::useDepth1
private

Definition at line 57 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

std::vector<G4LogicalVolume*> ECalSD::useDepth2
private

Definition at line 57 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

bool ECalSD::useWeight
private

Definition at line 51 of file ECalSD.h.

Referenced by ECalSD(), getEnergyDeposit(), and getRadiationLength().

std::map<G4LogicalVolume*,double> ECalSD::xtalLMap
private

Definition at line 56 of file ECalSD.h.

Referenced by crystalLength(), and initMap().