CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BetaCalculatorECAL.cc
Go to the documentation of this file.
3 
11 
14 
16 
17 using namespace susybsm;
18 
20  EBRecHitCollection_ (iConfig.getParameter<edm::InputTag>("EBRecHitCollection")),
21  EERecHitCollection_ (iConfig.getParameter<edm::InputTag>("EERecHitCollection"))
22 {
23  edm::ParameterSet trkParameters = iConfig.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
24  parameters_.loadParameters( trkParameters );
26 
27 }
28 
29 
31  bool setCalo = false;
33 
34  // EcalDetIdAssociator
35  iSetup.get<DetIdAssociatorRecord>().get("EcalDetIdAssociator", ecalDetIdAssociator_);
36  // Get the Bfield
37  iSetup.get<IdealMagneticFieldRecord>().get(bField_);
38  // Geometry
40  const CaloGeometry* theGeometry = theCaloGeometry_.product();
41  // Topology
42  edm::ESHandle<CaloTopology> pCaloTopology;
43  iSetup.get<CaloTopologyRecord>().get(pCaloTopology);
44  const CaloTopology* theCaloTopology = pCaloTopology.product();
45  // EcalRecHits
47  iEvent.getByLabel(EBRecHitCollection_,ebRecHits);
49  iEvent.getByLabel(EERecHitCollection_,eeRecHits);
50 
51  // select the track
52  reco::Track track;
53  if(candidate.hasTrackRef())
54  track = *(candidate.trackRef());
55  else
56  return; // in case there is no track ref, we can't do much
57 
58  // compute the track isolation
59  result.trkIsoDr=100;
60  for(reco::TrackCollection::const_iterator ndTrack = tracks->begin(); ndTrack != tracks->end(); ++ndTrack) {
61  double dr=sqrt(pow((track.outerEta()-ndTrack->outerEta()),2)+pow((track.outerPhi()-ndTrack->outerPhi()),2));
62  if(dr>0.00001 && dr<result.trkIsoDr) result.trkIsoDr=dr;
63  }
64 
65  // use the track associator to propagate to the calo
68  parameters_ );
69 
70  // do a custom propagation through Ecal
71  std::map<int,GlobalPoint> trackExitPositionMap; // rawId to exit position (subtracting cry center)
72  std::map<int,float> trackCrossedXtalCurvedMap; // rawId to trackLength
73 
75  // Build set of points in Ecal (necklace) using the propagator
76  std::vector<SteppingHelixStateInfo> neckLace;
77  neckLace = calcEcalDeposit(&tkInnerState,*ecalDetIdAssociator_);
78  // Initialize variables to be filled by the track-length function
79  double totalLengthCurved = 0.;
80  GlobalPoint internalPointCurved(0., 0., 0.);
81  GlobalPoint externalPointCurved(0., 0., 0.);
82  if(neckLace.size() > 1)
83  {
84  getDetailedTrackLengthInXtals(trackExitPositionMap,
85  trackCrossedXtalCurvedMap,
86  totalLengthCurved,
87  internalPointCurved,
88  externalPointCurved,
89  & (*theGeometry),
90  & (*theCaloTopology),
91  neckLace);
92  }
93 
94  // Make weighted sum of times
95  float sumWeightedTime = 0;
96  float sumTimeErrorSqr = 0;
97  float sumEnergy = 0;
98  float sumTrackLength = 0;
99  std::vector<EcalRecHit> crossedRecHits;
101 
102  std::map<int,GlobalPoint>::const_iterator trackExitMapIt = trackExitPositionMap.begin();
103  for(std::map<int,float>::const_iterator mapIt = trackCrossedXtalCurvedMap.begin();
104  mapIt != trackCrossedXtalCurvedMap.end(); ++mapIt)
105  {
106  if(DetId(mapIt->first).subdetId()==EcalBarrel)
107  {
108  EBDetId ebDetId(mapIt->first);
109  thisHit = ebRecHits->find(ebDetId);
110  if(thisHit == ebRecHits->end())
111  {
112  //std::cout << "\t Could not find crossedEcal detId: " << ebDetId << " in EBRecHitCollection!" << std::endl;
113  continue;
114  }
115  const EcalRecHit hit = *thisHit;
116  // Cut out badly-reconstructed hits
117  if(!hit.isTimeValid())
118  continue;
119  uint32_t rhFlag = hit.recoFlag();
120  if((rhFlag != EcalRecHit::kGood) && (rhFlag != EcalRecHit::kOutOfTime) && (rhFlag != EcalRecHit::kPoorCalib))
121  continue;
122 
123  float errorOnThis = hit.timeError();
124  sumTrackLength+=mapIt->second;
125  sumEnergy+=hit.energy();
126  crossedRecHits.push_back(hit);
127 // result.ecalSwissCrossKs.push_back(EcalSeverityLevelAlgo::spikeFromNeighbours(ebDetId,(*ebRecHits),0.2,EcalSeverityLevelAlgo::kSwissCross));
128 // result.ecalE1OverE9s.push_back(EcalSeverityLevelAlgo::spikeFromNeighbours(ebDetId,(*ebRecHits),0.2,EcalSeverityLevelAlgo::kE1OverE9));
129  result.ecalTrackLengths.push_back(mapIt->second);
130  result.ecalTrackExitPositions.push_back(trackExitMapIt->second);
131  result.ecalEnergies.push_back(hit.energy());
132  result.ecalTimes.push_back(hit.time());
133  result.ecalTimeErrors.push_back(hit.timeError());
134  result.ecalOutOfTimeEnergies.push_back(hit.outOfTimeEnergy());
135  result.ecalOutOfTimeChi2s.push_back(hit.outOfTimeChi2());
136  result.ecalChi2s.push_back(hit.chi2());
137  result.ecalDetIds.push_back(ebDetId);
138  // SIC DEBUG
139  //std::cout << " SIC DEBUG: time error on this crossed RecHit: " << errorOnThis << " energy of hit: "
140  // << hit.energy() << " time of hit: " << hit.time() << " trackLength: " << mapIt->second << std::endl;
141 
142  if(hit.isTimeErrorValid()) // use hit time for weighted time average
143  {
144  sumWeightedTime+=hit.time()/(errorOnThis*errorOnThis);
145  sumTimeErrorSqr+=1/(errorOnThis*errorOnThis);
146  }
147  }
148  trackExitMapIt++;
149  }
150 
151  if(crossedRecHits.size() > 0)
152  {
153  setCalo = true;
154  sort(crossedRecHits.begin(),crossedRecHits.end(),EcalRecHitLess());
155  result.ecalCrossedEnergy = sumEnergy;
156  result.ecalCrysCrossed = crossedRecHits.size();
157  result.ecalDeDx = sumEnergy/sumTrackLength;
158  // replace the below w/o trackassociator quantities?
161 
162  if(sumTimeErrorSqr > 0)
163  {
164  result.ecalTime = sumWeightedTime/sumTimeErrorSqr;
165  result.ecalTimeError = sqrt(1/sumTimeErrorSqr);
166  DetId maxEnergyId = crossedRecHits.begin()->id();
167 
168  if(maxEnergyId != DetId()) // double check
169  {
170  // To get beta, we assume photon propagation time is about the same for muons and e/gamma
171  // Since the typical path length is >> crystal length, this shouldn't be too bad
172  GlobalPoint position = info.getPosition(maxEnergyId); // position of crystal center on front face
173  double frontFaceR = sqrt(pow(position.x(),2)+pow(position.y(),2)+pow(position.z(),2));
174  double muonShowerMax = frontFaceR+11.5; // assume muon "showerMax" is halfway into the crystal
175  double gammaShowerMax = frontFaceR+6.23; // 7 X0 for e/gamma showerMax
176  double speedOfLight = 29.979; // cm/ns
177  result.ecalBeta = (muonShowerMax)/(result.ecalTime*speedOfLight+gammaShowerMax);
178  result.ecalBetaError = (speedOfLight*muonShowerMax*result.ecalTimeError)/pow(speedOfLight*result.ecalTime+gammaShowerMax,2);
179  result.ecalInvBetaError = speedOfLight*result.ecalTimeError/muonShowerMax;
180  }
181  // SIC debug
182  //std::cout << "BetaCalcEcal: CrossedRecHits: " << crossedRecHits.size()
183  // << " ecalTime: " << result.ecalTime << " timeError: " << result.ecalTimeError
184  // << " ecalCrossedEnergy: " << result.ecalCrossedEnergy << " ecalBeta: " << result.ecalBeta
185  // << " ecalBetaError: " << result.ecalBetaError << " ecalDeDx (MeV/cm): " << 1000*result.ecalDeDx << std::endl;
186  }
187  }
188 
189  if(info.crossedHcalRecHits.size() > 0)
190  {
191  // HCAL (not ECAL) info
194  //maxEnergyId = info.findMaxDeposition(TrackDetMatchInfo::HcalRecHits);
197  }
198 
199  if(setCalo)
200  caloInfo = result;
201 }
202 
203 std::vector<SteppingHelixStateInfo> BetaCalculatorECAL::calcEcalDeposit(const FreeTrajectoryState* tkInnerState,
204  const DetIdAssociator& associator)
205 {
206  // Set some parameters
207  double minR = associator.volume().minR () ;
208  double minZ = associator.volume().minZ () ;
209  double maxR = associator.volume().maxR () ;
210  double maxZ = associator.volume().maxZ () ;
211 
212  // Define the TrackOrigin (where the propagation starts)
213  SteppingHelixStateInfo trackOrigin(*tkInnerState);
214 
215  // Define Propagator
217  prop -> setMaterialMode(false);
218  prop -> applyRadX0Correction(true);
219 
220  return propagateThoughFromIP(trackOrigin,prop,associator.volume(), 500,0.1,minR,minZ,maxR,maxZ);
221 }
222 
223 int BetaCalculatorECAL::getDetailedTrackLengthInXtals(std::map<int,GlobalPoint>& trackExitPositionMap,
224  std::map<int,float>& trackCrossedXtalMap,
225  double& totalLengthCurved,
226  GlobalPoint& internalPointCurved,
227  GlobalPoint& externalPointCurved,
228  const CaloGeometry* theGeometry,
229  const CaloTopology* theTopology,
230  const std::vector<SteppingHelixStateInfo>& neckLace)
231 {
232  GlobalPoint origin (0., 0., 0.);
233  internalPointCurved = origin ;
234  externalPointCurved = origin ;
235 
236  bool firstPoint = false;
237  trackCrossedXtalMap.clear();
238 
239  const CaloSubdetectorGeometry* theBarrelSubdetGeometry = theGeometry->getSubdetectorGeometry(DetId::Ecal,1);
240  const CaloSubdetectorGeometry* theEndcapSubdetGeometry = theGeometry->getSubdetectorGeometry(DetId::Ecal,2);
241 
242  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = (neckLace.begin() + 1); itr != neckLace.end(); ++itr)
243  {
244  GlobalPoint probe_gp = (*itr).position();
245  std::vector<DetId> surroundingMatrix;
246 
247  EBDetId closestBarrelDetIdToProbe = ((theBarrelSubdetGeometry -> getClosestCell(probe_gp)).rawId());
248  EEDetId closestEndcapDetIdToProbe = ((theEndcapSubdetGeometry -> getClosestCell(probe_gp)).rawId());
249 
250  // check if the probe is inside the xtal
251  if( (closestEndcapDetIdToProbe) && (theGeometry->getSubdetectorGeometry(closestEndcapDetIdToProbe)->
252  getGeometry(closestEndcapDetIdToProbe)->inside(probe_gp)) )
253  {
254  double step = ((*itr).position() - (*(itr-1)).position()).mag();
255  GlobalPoint point = itr->position();
256  addStepToXtal(trackExitPositionMap, trackCrossedXtalMap, closestEndcapDetIdToProbe, step, point, theEndcapSubdetGeometry);
257  totalLengthCurved += step;
258 
259  if (firstPoint == false)
260  {
261  internalPointCurved = probe_gp ;
262  firstPoint = true ;
263  }
264 
265  externalPointCurved = probe_gp ;
266  }
267 
268  if( (closestBarrelDetIdToProbe) && (theGeometry->getSubdetectorGeometry(closestBarrelDetIdToProbe)->
269  getGeometry(closestBarrelDetIdToProbe)->inside(probe_gp)) )
270  {
271  double step = ((*itr).position() - (*(itr-1)).position()).mag();
272  GlobalPoint point = itr->position();
273  addStepToXtal(trackExitPositionMap, trackCrossedXtalMap, closestBarrelDetIdToProbe, step, point, theBarrelSubdetGeometry);
274  totalLengthCurved += step;
275 
276  if (firstPoint == false)
277  {
278  internalPointCurved = probe_gp ;
279  firstPoint = true ;
280  }
281 
282  externalPointCurved = probe_gp ;
283  }
284  else
285  {
286  // 3x3 matrix surrounding the probe
287  surroundingMatrix = theTopology->getSubdetectorTopology(closestBarrelDetIdToProbe)->getWindow(closestBarrelDetIdToProbe,3,3);
288 
289  for( unsigned int k=0; k<surroundingMatrix.size(); ++k ) {
290  if(theGeometry->getSubdetectorGeometry(surroundingMatrix.at(k))->getGeometry(surroundingMatrix.at(k))->inside(probe_gp))
291  {
292  double step = ((*itr).position() - (*(itr-1)).position()).mag();
293  GlobalPoint point = itr->position();
294  addStepToXtal(trackExitPositionMap, trackCrossedXtalMap, surroundingMatrix[k], step,
295  point, theGeometry->getSubdetectorGeometry(surroundingMatrix.at(k)));
296  totalLengthCurved += step;
297 
298  if (firstPoint == false)
299  {
300  internalPointCurved = probe_gp ;
301  firstPoint = true ;
302  }
303 
304  externalPointCurved = probe_gp ;
305  }
306  }
307 
308  // clear neighborhood matrix
309  surroundingMatrix.clear();
310  }
311  }
312 
313  return 0;
314 }
315 
316 void BetaCalculatorECAL::addStepToXtal(std::map<int,GlobalPoint>& trackExitPositionMap,
317  std::map<int,float>& trackCrossedXtalMap,
318  DetId aDetId,
319  float step,
321  const CaloSubdetectorGeometry* theSubdetGeometry)
322 {
323 
324  const CaloCellGeometry *cell_p = theSubdetGeometry->getGeometry(aDetId);
325  GlobalPoint p = (dynamic_cast <const TruncatedPyramid *> (cell_p))->getPosition(23);
326  GlobalPoint diff(point.x()-p.x(),point.y()-p.y(),point.z()-p.z());
327 
328  std::map<int,GlobalPoint>::iterator xtal = trackExitPositionMap.find(aDetId.rawId());
329  if (xtal!=trackExitPositionMap.end())
330  ((*xtal).second)=diff;
331  else
332  trackExitPositionMap.insert(std::pair<int,GlobalPoint>(aDetId.rawId(),diff));
333 
334  std::map<int,float>::iterator xtal2 = trackCrossedXtalMap.find(aDetId.rawId());
335  if (xtal2!= trackCrossedXtalMap.end())
336  ((*xtal2).second)+=step;
337  else
338  trackCrossedXtalMap.insert(std::pair<int,float>(aDetId.rawId(),step));
339 }
340 
341 
342 
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
GlobalPoint getPosition(const DetId &)
int getDetailedTrackLengthInXtals(std::map< int, GlobalPoint > &trackExitPositionMap, std::map< int, float > &trackCrossedXtalMap, double &totalLengthCurved, GlobalPoint &internalPointCurved, GlobalPoint &externalPointCurved, const CaloGeometry *theGeometry, const CaloTopology *theTopology, const std::vector< SteppingHelixStateInfo > &neckLace)
std::vector< float > ecalTrackLengths
Definition: HSCPCaloInfo.h:28
std::vector< float > ecalEnergies
Definition: HSCPCaloInfo.h:29
bool isTimeErrorValid() const
Definition: EcalRecHit.cc:121
std::vector< const HBHERecHit * > crossedHcalRecHits
list step
Definition: launcher.py:15
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
FreeTrajectoryState innerFreeState(const reco::Track &tk, const MagneticField *field)
double maxZ(bool withTolerance=true) const
std::vector< GlobalPoint > ecalTrackExitPositions
Definition: HSCPCaloInfo.h:36
void useDefaultPropagator()
use the default propagator
std::vector< EcalRecHit >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:62
BetaCalculatorECAL(const edm::ParameterSet &iConfig)
Flags recoFlag() const
DEPRECATED provided for temporary backward compatibility.
Definition: EcalRecHit.cc:133
void addInfoToCandidate(susybsm::HSCParticle &candidate, edm::Handle< reco::TrackCollection > &tracks, edm::Event &iEvent, const edm::EventSetup &iSetup, susybsm::HSCPCaloInfo &caloInfo)
float time() const
Definition: CaloRecHit.h:21
double minZ(bool withTolerance=true) const
bool hasTrackRef() const
Definition: HSCParticle.h:56
double nXnEnergy(const DetId &, EnergyType, int gridSize=1)
get energy of the NxN shape (N = 2*gridSize + 1) around given detector element
bool isTimeValid() const
Definition: EcalRecHit.cc:112
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
std::vector< DetId > ecalDetIds
Definition: HSCPCaloInfo.h:35
std::vector< float > ecalTimeErrors
Definition: HSCPCaloInfo.h:31
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
std::vector< float > ecalTimes
Definition: HSCPCaloInfo.h:30
int iEvent
Definition: GenABIO.cc:243
float outOfTimeEnergy() const
Definition: EcalRecHit.cc:40
float outOfTimeChi2() const
Definition: EcalRecHit.cc:34
float energy() const
Definition: CaloRecHit.h:19
T sqrt(T t)
Definition: SSEVec.h:46
float chi2() const
Definition: EcalRecHit.cc:26
T z() const
Definition: PV3DBase.h:63
tuple result
Definition: query.py:137
edm::ESHandle< DetIdAssociator > ecalDetIdAssociator_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double crossedEnergy(EnergyType)
energy in detector elements crossed by the track by types
float timeError() const
Definition: EcalRecHit.cc:95
std::vector< float > ecalChi2s
Definition: HSCPCaloInfo.h:32
int k[5][pyjets_maxn]
double minR(bool withTolerance=true) const
TrackAssociatorParameters parameters_
Definition: DetId.h:20
double outerEta() const
pseudorapidity of the momentum vector at the outermost hit position
Definition: Track.h:89
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
edm::InputTag EERecHitCollection_
T const * product() const
Definition: ESHandle.h:62
void addStepToXtal(std::map< int, GlobalPoint > &trackExitPositionMap, std::map< int, float > &trackCrossedXtalMap, DetId aDetId, float step, GlobalPoint point, const CaloSubdetectorGeometry *theSubdetGeometry)
const CaloSubdetectorTopology * getSubdetectorTopology(const DetId &id) const
access the subdetector Topology for the given subdetector directly
Definition: CaloTopology.cc:26
A base class to handle the particular shape of Ecal Xtals. Taken from ORCA Calorimetry Code...
reco::TrackRef trackRef() const
Definition: HSCParticle.h:70
TrackDetectorAssociator trackAssociator_
static int position[264][3]
Definition: ReadPGInfo.cc:509
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
std::vector< float > ecalOutOfTimeChi2s
Definition: HSCPCaloInfo.h:33
const FiducialVolume & volume() const
get active detector volume
std::vector< SteppingHelixStateInfo > propagateThoughFromIP(const SteppingHelixStateInfo &state, const Propagator *prop, const FiducialVolume &volume, int nsteps, float step, float minR, float minZ, float maxR, float maxZ)
edm::ESHandle< MagneticField > bField_
std::vector< float > ecalOutOfTimeEnergies
Definition: HSCPCaloInfo.h:34
T x() const
Definition: PV3DBase.h:61
double outerPhi() const
azimuthal angle of the momentum vector at the outermost hit position
Definition: Track.h:87
edm::ESHandle< CaloGeometry > theCaloGeometry_
double maxR(bool withTolerance=true) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
void loadParameters(const edm::ParameterSet &)
std::vector< SteppingHelixStateInfo > calcEcalDeposit(const FreeTrajectoryState *tkInnerState, const DetIdAssociator &associator)
const double speedOfLight
edm::InputTag EBRecHitCollection_