17 const DetId coreDet = vdetIds[trkIndex].detIdECAL;
21 std::cout <<
"DetId " << (
EBDetId)(coreDet) <<
" Flag " << vdetIds[trkIndex].okECAL << std::endl;
23 std::cout <<
"DetId " << (
EEDetId)(coreDet) <<
" Flag " << vdetIds[trkIndex].okECAL << std::endl;
26 double maxNearP = -1.0;
27 if (vdetIds[trkIndex].okECAL) {
28 std::vector<DetId> vdets =
spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
30 if (debug)
std::cout <<
"chargeIsolationEcal:: eta/phi/dets " << ieta <<
" " << iphi <<
" " << vdets.size() << std::endl;
33 for (
unsigned int indx=0; indx<vdetIds.size(); ++indx) {
34 if (indx != trkIndex && vdetIds[indx].
ok && vdetIds[indx].okECAL) {
35 const DetId anyCell = vdetIds[indx].detIdECAL;
37 const reco::Track* pTrack = &(*(vdetIds[indx].trkItr));
38 if (maxNearP < pTrack->
p()) maxNearP = pTrack->
p();
51 std::vector<DetId> vdets =
spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
53 if (debug)
std::cout <<
"chargeIsolation:: eta/phi/dets " << ieta <<
" " << iphi <<
" " << vdets.size() << std::endl;
55 double maxNearP = -1.0;
59 reco::TrackCollection::const_iterator trkItr2;
60 for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
65 (pTrack2->
quality(trackQuality_)) :
true;
66 if ( (trkItr2 != trkItr) && trkQuality ) {
69 const GlobalPoint point2(info.first.x(),info.first.y(),info.first.z());
76 if (debug)
std::cout <<
"chargeIsolationEcal Cell " << (
EBDetId)(anyCell) <<
" pt " << pTrack2->
p() << std::endl;
78 if (maxNearP<pTrack2->
p()) maxNearP=pTrack2->
p();
85 if (debug)
std::cout <<
"chargeIsolationEcal Cell " << (
EEDetId)(anyCell) <<
" pt " << pTrack2->
p() << std::endl;
87 if (maxNearP<pTrack2->
p()) maxNearP=pTrack2->
p();
99 std::vector<DetId> dets(1,vdetIds[trkIndex].detIdHCAL);
102 std::cout <<
"DetId " << (
HcalDetId)(dets[0]) <<
" Flag " << vdetIds[trkIndex].okHCAL << std::endl;
105 double maxNearP = -1.0;
106 if (vdetIds[trkIndex].okHCAL) {
107 std::vector<DetId> vdets =
spr::matrixHCALIds(dets, topology, ieta, iphi,
false, debug);
109 if (debug)
std::cout <<
"chargeIsolationHcal:: eta/phi/dets " << ieta <<
" " << iphi <<
" " << vdets.size() << std::endl;
111 for (
unsigned indx = 0; indx<vdetIds.size(); ++indx) {
112 if (indx != trkIndex && vdetIds[indx].
ok && vdetIds[indx].okHCAL) {
113 const DetId anyCell = vdetIds[indx].detIdHCAL;
115 const reco::Track* pTrack = &(*(vdetIds[indx].trkItr));
117 if (debug)
std::cout <<
"chargeIsolationHcal Cell " << (
HcalDetId)(anyCell) <<
" pt " << pTrack->
p() << std::endl;
119 if (maxNearP < pTrack->
p()) maxNearP = pTrack->
p();
131 std::vector<DetId> dets(1,ClosestCell);
135 std::vector<DetId> vdets =
spr::matrixHCALIds(dets, topology, ieta, iphi,
false, debug);
136 std::vector<DetId>::iterator it;
140 for (
unsigned int i=0;
i<vdets.size();
i++) {
141 std::cout <<
"HcalDetId in " <<2*ieta+1 <<
"x" << 2*iphi+1 <<
" " << (
HcalDetId) vdets[
i] << std::endl;
145 double maxNearP = -1.0;
148 reco::TrackCollection::const_iterator trkItr2;
149 for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
154 (pTrack2->
quality(trackQuality_)) :
true;
155 if ( (trkItr2 != trkItr) && trkQuality ) {
157 const GlobalPoint point2(info.first.x(),info.first.y(),info.first.z());
161 std::cout <<
"Track2 (p,eta,phi) " << pTrack2->
p() <<
" " << pTrack2->
eta() <<
" " << pTrack2->
phi() << std::endl;
167 if(maxNearP<pTrack2->
p()) maxNearP=pTrack2->
p();
171 std::cout <<
"maxNearP " << maxNearP <<
" thisCell " 173 << info.first.x() <<
"," << info.first.y() <<
"," 174 << info.first.z() <<
")" << std::endl;
184 bool isIsolated =
true;
185 for (
unsigned int i=0;
i<vdets.size();
i++){
186 if (anyCell == vdets[
i] ) {
194 double coneChargeIsolation(
const edm::Event&
iEvent,
const edm::EventSetup& iSetup, reco::TrackCollection::const_iterator trkItr,
edm::Handle<reco::TrackCollection> trkCollection,
TrackDetectorAssociator&
associator,
TrackAssociatorParameters& parameters_,
const std::string & theTrackQuality,
int &nNearTRKs,
int &nLayers_maxNearP,
int &trkQual_maxNearP,
double &maxNearP_goodTrk,
const GlobalPoint& hpoint1,
const GlobalVector& trackMom,
double dR) {
199 maxNearP_goodTrk = -999.0;
200 double maxNearP = -999.0;
204 reco::TrackCollection::const_iterator trkItr2;
205 for( trkItr2 = trkCollection->begin();
206 trkItr2 != trkCollection->end(); ++trkItr2){
213 (pTrack2->
quality(trackQuality_)) :
true;
214 if (trkQuality) trkQual_maxNearP = 1;
220 if (trkItr2 != trkItr) {
235 if (isConeChargedIso==0) {
237 if(maxNearP<pTrack2->
p()) {
238 maxNearP=pTrack2->
p();
239 if (trkQual_maxNearP>0 && nLayers_maxNearP>7 && maxNearP_goodTrk<pTrack2->
p()) {
240 maxNearP_goodTrk=pTrack2->
p();
251 double chargeIsolationCone(
unsigned int trkIndex, std::vector<spr::propagatedTrackDirection> & trkDirs,
double dR,
int & nNearTRKs,
bool debug) {
253 double maxNearP = -1.0;
255 if (trkDirs[trkIndex].okHCAL) {
257 if (debug)
std::cout <<
"chargeIsolationCone with " << trkDirs.size() <<
" tracks " << std::endl;
259 for (
unsigned int indx=0; indx<trkDirs.size(); ++indx) {
260 if (indx != trkIndex && trkDirs[indx].
ok && trkDirs[indx].okHCAL) {
261 int isConeChargedIso =
spr::coneChargeIsolation(trkDirs[trkIndex].pointHCAL, trkDirs[indx].pointHCAL, trkDirs[trkIndex].directionHCAL, dR);
262 if (isConeChargedIso==0) {
264 const reco::Track* pTrack = &(*(trkDirs[indx].trkItr));
265 if (maxNearP < pTrack->
p()) maxNearP = pTrack->
p();
271 if (debug)
std::cout <<
"chargeIsolationCone Track " << trkDirs[trkIndex].okHCAL <<
" maxNearP " << maxNearP << std::endl;
276 std::pair<double,double>
chargeIsolationCone(
unsigned int trkIndex, std::vector<spr::propagatedTrackDirection> & trkDirs,
double dR,
bool debug) {
278 double maxNearP = -1.0;
280 if (trkDirs[trkIndex].okHCAL) {
282 if (debug)
std::cout <<
"chargeIsolationCone with " << trkDirs.size() <<
" tracks " << std::endl;
284 for (
unsigned int indx=0; indx<trkDirs.size(); ++indx) {
285 if (indx != trkIndex && trkDirs[indx].
ok && trkDirs[indx].okHCAL) {
286 int isConeChargedIso =
spr::coneChargeIsolation(trkDirs[trkIndex].pointHCAL, trkDirs[indx].pointHCAL, trkDirs[trkIndex].directionHCAL, dR);
287 if (isConeChargedIso==0) {
288 const reco::Track* pTrack = &(*(trkDirs[indx].trkItr));
289 sumP += (pTrack->
p());
290 if (maxNearP < pTrack->
p()) maxNearP = pTrack->
p();
296 if (debug)
std::cout <<
"chargeIsolationCone Track " << trkDirs[trkIndex].okHCAL <<
" maxNearP " << maxNearP <<
":" << sumP <<std::endl;
298 return std::pair<double,double>(maxNearP,sumP);
double p() const
momentum vector magnitude
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
double coneChargeIsolation(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::TrackCollection::const_iterator trkItr, edm::Handle< reco::TrackCollection > trkCollection, TrackDetectorAssociator &associator, TrackAssociatorParameters ¶meters_, const std::string &theTrackQuality, int &nNearTRKs, int &nLayers_maxNearP, int &trkQual_maxNearP, double &maxNearP_goodTrk, const GlobalPoint &hpoint1, const GlobalVector &trackMom, double dR)
static const double etaBEEcal
double getDistInPlaneTrackDir(const GlobalPoint &caloPoint, const GlobalVector &caloVector, const GlobalPoint &rechitPoint, bool debug=false)
CaloTopology const * topology(0)
void matrixECALIds(const DetId &det, int ieta, int iphi, const CaloGeometry *geo, const CaloTopology *caloTopology, std::vector< DetId > &vdets, bool debug=false, bool igNoreTransition=true)
bool chargeIsolation(const DetId anyCell, std::vector< DetId > &vdets)
TrackQuality
track quality
std::pair< math::XYZPoint, bool > propagateHCAL(const reco::Track *, const MagneticField *, bool debug=false)
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
double phi() const
azimuthal angle of momentum vector
double chargeIsolationCone(unsigned int trkIndex, std::vector< spr::propagatedTrackDirection > &trkDirs, double dR, int &nNearTRKs, bool debug=false)
int trackerLayersWithMeasurement() const
math::XYZPoint trkGlobPosAtHcal
double chargeIsolationEcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, bool debug=false)
double eta() const
pseudorapidity of momentum vector
Abs< T >::type abs(const T &t)
std::pair< math::XYZPoint, bool > propagateECAL(const reco::Track *, const MagneticField *, bool debug=false)
int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum) ...
virtual DetId getClosestCell(const GlobalPoint &r) const
static TrackQuality qualityByName(const std::string &name)
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
bool quality(const TrackQuality) const
Track quality.
std::vector< DetId > matrixHCALIds(std::vector< DetId > &dets, const HcalTopology *topology, int ieta, int iphi, bool includeHO=false, bool debug=false)
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)