CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
HDetIdAssociator Class Referenceabstract

#include <DetIdAssociator.h>

Inheritance diagram for HDetIdAssociator:
HCaloDetIdAssociator HEcalDetIdAssociator HHcalDetIdAssociator

Public Member Functions

virtual std::set< DetIdgetCrossedDetIds (const std::set< DetId > &, const std::vector< GlobalPoint > &trajectory)
 
virtual std::set< DetIdgetDetIdsCloseToAPoint (const GlobalPoint &, const int idR=0)
 
virtual std::set< DetIdgetDetIdsCloseToAPoint (const GlobalPoint &point, const double dR=0)
 
virtual std::set< DetIdgetDetIdsInACone (const std::set< DetId > &, const std::vector< GlobalPoint > &trajectory, const double)
 
virtual std::set< DetIdgetMaxEDetId (const std::set< DetId > &, edm::Handle< CaloTowerCollection > caloTowers)
 
virtual std::set< DetIdgetMaxEDetId (const std::set< DetId > &, edm::Handle< HBHERecHitCollection > recHits)
 
virtual std::vector< GlobalPointgetTrajectory (const FreeTrajectoryState &, const std::vector< GlobalPoint > &)
 
 HDetIdAssociator ()
 
 HDetIdAssociator (const int nPhi, const int nEta, const double etaBinSize)
 
virtual int iEta (const GlobalPoint &)
 
virtual int iPhi (const GlobalPoint &)
 
virtual void setPropagator (Propagator *ptr)
 
virtual ~HDetIdAssociator ()
 

Protected Member Functions

virtual void buildMap ()
 
virtual void check_setup ()
 
virtual std::set< DetIdgetASetOfValidDetIds ()=0
 
virtual std::vector< GlobalPointgetDetIdPoints (const DetId &)=0
 
virtual GlobalPoint getPosition (const DetId &)=0
 
virtual bool insideElement (const GlobalPoint &, const DetId &)=0
 
virtual bool nearElement (const GlobalPoint &point, const DetId &id, const double distance)
 

Protected Attributes

const double etaBinSize_
 
PropagatorivProp_
 
const int nEta_
 
const int nPhi_
 
std::vector< std::vector
< std::set< DetId > > > * 
theMap_
 

Detailed Description

Definition at line 30 of file DetIdAssociator.h.

Constructor & Destructor Documentation

HDetIdAssociator::HDetIdAssociator ( )
inline

Definition at line 32 of file DetIdAssociator.h.

32 :theMap_(0),nPhi_(0),nEta_(0),etaBinSize_(0),ivProp_(0){};
const double etaBinSize_
std::vector< std::vector< std::set< DetId > > > * theMap_
Propagator * ivProp_
HDetIdAssociator::HDetIdAssociator ( const int  nPhi,
const int  nEta,
const double  etaBinSize 
)
inline

Definition at line 33 of file DetIdAssociator.h.

const double etaBinSize_
std::vector< std::vector< std::set< DetId > > > * theMap_
Propagator * ivProp_
virtual HDetIdAssociator::~HDetIdAssociator ( )
inlinevirtual

Definition at line 36 of file DetIdAssociator.h.

36 {};

Member Function Documentation

void HDetIdAssociator::buildMap ( )
protectedvirtual

Definition at line 249 of file DetIdAssociator.cc.

References assert(), check_setup(), PV3DBase< T, PVType, FrameType >::eta(), etaBinSize_, HLT_25ns14e33_v1_cff::etaMax, HLT_25ns14e33_v1_cff::etaMin, getASetOfValidDetIds(), getPosition(), iEta(), iPhi(), LogTrace, nEta_, nPhi_, point, std::swap(), and theMap_.

Referenced by getDetIdsCloseToAPoint().

250 {
251 // modified version: take only detector central position
252  check_setup();
253  LogTrace("HDetIdAssociator")<<"building map" << "\n";
254  if(theMap_) delete theMap_;
255  theMap_ = new std::vector<std::vector<std::set<DetId> > >(nEta_,std::vector<std::set<DetId> >(nPhi_));
256  int numberOfDetIdsOutsideEtaRange = 0;
257  int numberOfDetIdsActive = 0;
258  std::set<DetId> validIds = getASetOfValidDetIds();
259  for (std::set<DetId>::const_iterator id_itr = validIds.begin(); id_itr!=validIds.end(); id_itr++) {
260 // std::vector<GlobalPoint> points = getDetIdPoints(*id_itr);
261  GlobalPoint point = getPosition(*id_itr);
262 // reject fake DetIds (eta=0 - what are they anyway???)
263  if(point.eta()==0)continue;
264 
265  int ieta = iEta(point);
266  int iphi = iPhi(point);
267  int etaMax(-1);
268  int etaMin(nEta_);
269  int phiMax(-1);
270  int phiMin(nPhi_);
271  if ( iphi >= nPhi_ ) iphi = iphi % nPhi_;
272  assert (iphi>=0);
273  if ( etaMin > ieta) etaMin = ieta;
274  if ( etaMax < ieta) etaMax = ieta;
275  if ( phiMin > iphi) phiMin = iphi;
276  if ( phiMax < iphi) phiMax = iphi;
277 // for abs(eta)>1.8 one tower covers two phi segments
278  if ((ieta>54||ieta<15) && iphi%2==0) phiMax++;
279  if ((ieta>54||ieta<15) && iphi%2==1) phiMin--;
280 
281  if (etaMax<0||phiMax<0||etaMin>=nEta_||phiMin>=nPhi_) {
282  LogTrace("HDetIdAssociator")<<"Out of range: DetId:" << id_itr->rawId() <<
283  "\n\teta (min,max): " << etaMin << "," << etaMax <<
284  "\n\tphi (min,max): " << phiMin << "," << phiMax <<
285  "\nTower id: " << id_itr->rawId() << "\n";
286  numberOfDetIdsOutsideEtaRange++;
287  continue;
288  }
289 
290  if (phiMax-phiMin > phiMin+nPhi_-phiMax){
291  phiMin += nPhi_;
292  std::swap(phiMin,phiMax);
293  }
294  for(int ieta = etaMin; ieta <= etaMax; ieta++)
295  for(int iphi = phiMin; iphi <= phiMax; iphi++)
296  (*theMap_)[ieta][iphi%nPhi_].insert(*id_itr);
297  numberOfDetIdsActive++;
298  }
299  LogTrace("HDetIdAssociator") << "Number of elements outside the allowed range ( |eta|>"<<
300  nEta_/2*etaBinSize_ << "): " << numberOfDetIdsOutsideEtaRange << "\n";
301  LogTrace("HDetIdAssociator") << "Number of active DetId's mapped: " <<
302  numberOfDetIdsActive << "\n";
303 }
virtual void check_setup()
const double etaBinSize_
virtual GlobalPoint getPosition(const DetId &)=0
std::vector< std::vector< std::set< DetId > > > * theMap_
assert(m_qm.get())
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
virtual int iPhi(const GlobalPoint &)
virtual std::set< DetId > getASetOfValidDetIds()=0
#define LogTrace(id)
virtual int iEta(const GlobalPoint &)
T eta() const
Definition: PV3DBase.h:76
*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
virtual void HDetIdAssociator::check_setup ( )
inlineprotectedvirtual

Reimplemented in HCaloDetIdAssociator.

Definition at line 70 of file DetIdAssociator.h.

References etaBinSize_, Exception, ivProp_, nEta_, and nPhi_.

Referenced by buildMap(), HCaloDetIdAssociator::check_setup(), getCrossedDetIds(), getDetIdsCloseToAPoint(), getDetIdsInACone(), getMaxEDetId(), and getTrajectory().

71  {
72  if (nEta_==0) throw cms::Exception("FatalError") << "Number of eta bins is not set.\n";
73  if (nPhi_==0) throw cms::Exception("FatalError") << "Number of phi bins is not set.\n";
74  if (ivProp_==0) throw cms::Exception("FatalError") << "Track propagator is not defined\n";
75  if (etaBinSize_==0) throw cms::Exception("FatalError") << "Eta bin size is not set.\n";
76  }
const double etaBinSize_
Propagator * ivProp_
virtual std::set<DetId> HDetIdAssociator::getASetOfValidDetIds ( )
protectedpure virtual
std::set< DetId > HDetIdAssociator::getCrossedDetIds ( const std::set< DetId > &  inset,
const std::vector< GlobalPoint > &  trajectory 
)
virtual

Definition at line 351 of file DetIdAssociator.cc.

References check_setup(), and insideElement().

353 {
354  check_setup();
355  std::set<DetId> outset;
356  for(std::set<DetId>::const_iterator id_iter = inset.begin(); id_iter != inset.end(); id_iter++)
357  for(std::vector<GlobalPoint>::const_iterator point_iter = trajectory.begin(); point_iter != trajectory.end(); point_iter++)
358  if (insideElement(*point_iter, *id_iter)) outset.insert(*id_iter);
359  return outset;
360 }
virtual void check_setup()
virtual bool insideElement(const GlobalPoint &, const DetId &)=0
virtual std::vector<GlobalPoint> HDetIdAssociator::getDetIdPoints ( const DetId )
protectedpure virtual

Implemented in HCaloDetIdAssociator.

std::set< DetId > HDetIdAssociator::getDetIdsCloseToAPoint ( const GlobalPoint direction,
const int  idR = 0 
)
virtual

Definition at line 143 of file DetIdAssociator.cc.

References begin, buildMap(), check_setup(), end, PV3DBase< T, PVType, FrameType >::eta(), i, iEta(), iPhi(), j, LogTrace, nEta_, nPhi_, PV3DBase< T, PVType, FrameType >::phi(), and theMap_.

Referenced by getDetIdsCloseToAPoint().

145 {
146  std::set<DetId> set;
147  check_setup();
148  int nDets=0;
149  if (! theMap_) buildMap();
150  LogTrace("MatchPoint") << "point (eta,phi): " << direction.eta() << "," << direction.phi() << "\n";
151  int ieta = iEta(direction);
152  int iphi = iPhi(direction);
153 
154  LogTrace("MatchPoint") << "(ieta,iphi): " << ieta << "," << iphi << "\n";
155 
156  if (ieta>=0 && ieta<nEta_ && iphi>=0 && iphi<nPhi_){
157  set = (*theMap_)[ieta][iphi];
158  nDets++;
159  if (idR>0){
160  LogTrace("MatchPoint") << "Add neighbors (ieta,iphi): " << ieta << "," << iphi << "\n";
161  //add neighbors
162  int maxIEta = ieta+idR;
163  int minIEta = ieta-idR;
164  if(maxIEta>=nEta_) maxIEta = nEta_-1;
165  if(minIEta<0) minIEta = 0;
166  int maxIPhi = iphi+idR;
167  int minIPhi = iphi-idR;
168  if(minIPhi<0) {
169  minIPhi+=nPhi_;
170  maxIPhi+=nPhi_;
171  }
172  LogTrace("MatchPoint") << "\tieta (min,max): " << minIEta << "," << maxIEta<< "\n";
173  LogTrace("MatchPoint") << "\tiphi (min,max): " << minIPhi << "," << maxIPhi<< "\n";
174  for (int i=minIEta;i<=maxIEta;i++)
175  for (int j=minIPhi;j<=maxIPhi;j++) {
176  if( i==ieta && j==iphi) continue; // already in the set
177  set.insert((*theMap_)[i][j%nPhi_].begin(),(*theMap_)[i][j%nPhi_].end());
178  nDets++;
179  }
180  }
181  }
182 // if(set.size() > 0) {
183 // if (ieta+idR<55 && ieta-idR>14 && set.size() != (2*idR+1)*(2*idR+1)){
184 // std::cout<<" RRRA: "<<set.size()<<" DetIds in region "<<ieta<<" "<<iphi<<std::endl;
185 // for( std::set<DetId>::const_iterator itr=set.begin(); itr!=set.end(); itr++) {
186 // GlobalPoint point = getPosition(*itr);
187 // std::cout << "DetId: " <<itr->rawId() <<" (eta,phi): " << point.eta() << "," << point.phi()<<" "<<iEta(point)<<" "<<iPhi(point)<<std::endl;
188 // }
189 // }
190 // else {
191 // std::cout <<" HDetIdAssociator::getDetIdsCloseToAPoint::There are strange days "<<std::endl;
192 // }
193  return set;
194 }
virtual void check_setup()
virtual void buildMap()
int i
Definition: DBlmapReader.cc:9
std::vector< std::vector< std::set< DetId > > > * theMap_
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
virtual int iPhi(const GlobalPoint &)
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
#define LogTrace(id)
virtual int iEta(const GlobalPoint &)
T eta() const
Definition: PV3DBase.h:76
#define begin
Definition: vmac.h:30
virtual std::set<DetId> HDetIdAssociator::getDetIdsCloseToAPoint ( const GlobalPoint point,
const double  dR = 0 
)
inlinevirtual

Definition at line 44 of file DetIdAssociator.h.

References PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, etaBinSize_, getDetIdsCloseToAPoint(), and nPhi_.

46  {
47  int etaIdR = int(dR/etaBinSize_);
48  int phiIdR = int(dR/(2*3.1416)*nPhi_);
49  if (etaIdR>phiIdR)
50  return getDetIdsCloseToAPoint(point, 1+etaIdR);
51  else
52  return getDetIdsCloseToAPoint(point, 1+phiIdR);
53  }
virtual std::set< DetId > getDetIdsCloseToAPoint(const GlobalPoint &, const int idR=0)
const double etaBinSize_
std::set< DetId > HDetIdAssociator::getDetIdsInACone ( const std::set< DetId > &  inset,
const std::vector< GlobalPoint > &  trajectory,
const double  dR 
)
virtual

Definition at line 306 of file DetIdAssociator.cc.

References begin, check_setup(), end, getPosition(), i, iEta(), iPhi(), j, nearElement(), nEta_, nPhi_, point, findQualityFiles::size, and theMap_.

309 {
310 // modified version: if dR<0, returns 3x3 towers around the input one (Michal)
311  check_setup();
312  std::set<DetId> outset;
313 
314  if(dR>=0) {
315  for(std::set<DetId>::const_iterator id_iter = inset.begin(); id_iter != inset.end(); id_iter++)
316  for(std::vector<GlobalPoint>::const_iterator point_iter = trajectory.begin(); point_iter != trajectory.end(); point_iter++)
317  if (nearElement(*point_iter,*id_iter,dR)) outset.insert(*id_iter);
318  }
319  else {
320  if (inset.size()!=1) return outset;
321  std::set<DetId>::const_iterator id_inp = inset.begin();
322  int ieta;
323  int iphi;
324  GlobalPoint point = getPosition(*id_inp);
325  ieta = iEta(point);
326  iphi = iPhi(point);
327  for (int i=ieta-1;i<=ieta+1;i++) {
328  for (int j=iphi-1;j<=iphi+1;j++) {
329 // if( i==ieta && j==iphi) continue;
330  if( i<0 || i>=nEta_) continue;
331  int j2fill = j%nPhi_;
332  if(j2fill<0) j2fill+=nPhi_;
333  if((*theMap_)[i][j2fill].size()==0)continue;
334  outset.insert((*theMap_)[i][j2fill].begin(),(*theMap_)[i][j2fill].end());
335  }
336  }
337  }
338 
339 // if (outset.size() > 0) {
340 // std::cout<<" RRRA: DetIds in cone:"<<std::endl;
341 // for( std::set<DetId>::const_iterator itr=outset.begin(); itr!=outset.end(); itr++) {
342 // GlobalPoint point = getPosition(*itr);
343 // std::cout << "DetId: " <<itr->rawId() <<" (eta,phi): " << point.eta() << "," << point.phi()<<std::endl;
344 // }
345 // }
346 
347  return outset;
348 }
virtual void check_setup()
int i
Definition: DBlmapReader.cc:9
virtual GlobalPoint getPosition(const DetId &)=0
std::vector< std::vector< std::set< DetId > > > * theMap_
virtual bool nearElement(const GlobalPoint &point, const DetId &id, const double distance)
virtual int iPhi(const GlobalPoint &)
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
virtual int iEta(const GlobalPoint &)
#define begin
Definition: vmac.h:30
tuple size
Write out results.
*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
std::set< DetId > HDetIdAssociator::getMaxEDetId ( const std::set< DetId > &  inset,
edm::Handle< CaloTowerCollection caloTowers 
)
virtual

Definition at line 363 of file DetIdAssociator.cc.

References check_setup().

365 {
366 // returns the most energetic tower in the NxN box (Michal)
367  check_setup();
368  std::set<DetId> outset;
369  std::set<DetId>::const_iterator id_max = inset.begin();
370  double Ehadmax=0;
371 
372  for(std::set<DetId>::const_iterator id_iter = inset.begin(); id_iter != inset.end(); id_iter++) {
373  DetId id(*id_iter);
374 // GlobalPoint point = getPosition(*id_iter);
375 // int ieta = iEta(point);
376 // int iphi = iPhi(point);
377  CaloTowerCollection::const_iterator tower = (*caloTowers).find(id);
378  if(tower != (*caloTowers).end() && tower->hadEnergy()>Ehadmax) {
379  id_max = id_iter;
380  Ehadmax = tower->hadEnergy();
381  }
382  }
383 
384  if (Ehadmax > 0) outset.insert(*id_max);
385 
386 // if (outset.size() > 0) {
387 // std::cout<<" RRRA: Most energetic DetId:"<<std::endl;
388 // for( std::set<DetId>::const_iterator itr=outset.begin(); itr!=outset.end(); itr++) {
389 // GlobalPoint point = getPosition(*itr);
390 // std::cout << "DetId: " <<itr->rawId() <<" (eta,phi): " << point.eta() << "," << point.phi()<<std::endl;
391 // }
392 // }
393 
394  return outset;
395 }
virtual void check_setup()
std::vector< CaloTower >::const_iterator const_iterator
Definition: DetId.h:18
std::set< DetId > HDetIdAssociator::getMaxEDetId ( const std::set< DetId > &  inset,
edm::Handle< HBHERecHitCollection recHits 
)
virtual

Definition at line 398 of file DetIdAssociator.cc.

References check_setup().

400 {
401 // returns the most energetic tower in the NxN box - from RecHits (Michal)
402  check_setup();
403  std::set<DetId> outset;
404  std::set<DetId>::const_iterator id_max = inset.begin();
405  double Ehadmax=0;
406 
407  for(std::set<DetId>::const_iterator id_iter = inset.begin(); id_iter != inset.end(); id_iter++) {
408  DetId id(*id_iter);
409 // GlobalPoint point = getPosition(*id_iter);
410 // int ieta = iEta(point);
411 // int iphi = iPhi(point);
412  HBHERecHitCollection::const_iterator hit = (*recHits).find(id);
413  if(hit != (*recHits).end() && hit->energy()>Ehadmax) {
414  id_max = id_iter;
415  Ehadmax = hit->energy();
416  }
417  }
418 
419  if (Ehadmax > 0) outset.insert(*id_max);
420 
421 // if (outset.size() > 0) {
422 // std::cout<<" RRRA: Most energetic DetId:"<<std::endl;
423 // for( std::set<DetId>::const_iterator itr=outset.begin(); itr!=outset.end(); itr++) {
424 // GlobalPoint point = getPosition(*itr);
425 // std::cout << "DetId: " <<itr->rawId() <<" (eta,phi): " << point.eta() << "," << point.phi()<<std::endl;
426 // }
427 // }
428 
429  return outset;
430 }
virtual void check_setup()
std::vector< HBHERecHit >::const_iterator const_iterator
Definition: DetId.h:18
virtual GlobalPoint HDetIdAssociator::getPosition ( const DetId )
protectedpure virtual

Implemented in HCaloDetIdAssociator.

Referenced by buildMap(), getDetIdsInACone(), and nearElement().

std::vector< GlobalPoint > HDetIdAssociator::getTrajectory ( const FreeTrajectoryState ftsStart,
const std::vector< GlobalPoint > &  surfaces 
)
virtual

Definition at line 25 of file DetIdAssociator.cc.

References funct::abs(), check_setup(), PV3DBase< T, PVType, FrameType >::eta(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), ivProp_, LogTrace, PV3DBase< T, PVType, FrameType >::mag(), FreeTrajectoryState::momentum(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), point, FreeTrajectoryState::position(), Propagator::propagate(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

27 {
28  check_setup();
29  std::vector<GlobalPoint> trajectory;
30  TrajectoryStateOnSurface tSOSDest;
31  FreeTrajectoryState ftsCurrent = ftsStart;
32 
33  for(std::vector<GlobalPoint>::const_iterator surface_iter = surfaces.begin();
34  surface_iter != surfaces.end(); surface_iter++) {
35  // this stuff is some weird pointer, which destroy itself
36  Cylinder *cylinder = new Cylinder(surface_iter->perp(),
37  Surface::PositionType(0,0,0),
39  Plane *forwardEndcap = new Plane(Surface::PositionType(0,0,surface_iter->z()),
41  Plane *backwardEndcap = new Plane(Surface::PositionType(0,0,-surface_iter->z()),
43 
44 
45  LogTrace("StartingPoint")<< "Propagate from "<< "\n"
46  << "\tx: " << ftsStart.position().x()<< "\n"
47  << "\ty: " << ftsStart.position().y()<< "\n"
48  << "\tz: " << ftsStart.position().z()<< "\n"
49  << "\tmomentum eta: " << ftsStart.momentum().eta()<< "\n"
50  << "\tmomentum phi: " << ftsStart.momentum().phi()<< "\n"
51  << "\tmomentum: " << ftsStart.momentum().mag()<< "\n";
52 
53  float tanTheta = ftsCurrent.momentum().perp()/ftsCurrent.momentum().z();
54  float corner = surface_iter->perp()/surface_iter->z();
55 /*
56  std::cout<<"Propagate from "<< "\n"
57  << "\tx: " << ftsCurrent.position().x()<< "\n"
58  << "\ty: " << ftsCurrent.position().y()<< "\n"
59  << "\tz: " << ftsCurrent.position().z()<< "\n"
60  << "\tz: " << ftsCurrent.position().perp()<< "\n"
61  << "\tz: " << tanTheta<<" "<< corner <<"\n"
62  << "\tmomentum eta: " << ftsCurrent.momentum().eta()<< "\n"
63  << "\tmomentum phi: " << ftsCurrent.momentum().phi()<< "\n"
64  << "\tmomentum: " << ftsCurrent.momentum().mag()<<std::endl;
65 */
66  // First propage the track to the cylinder if |eta|<1, othewise to the encap
67  // and correct depending on the result
68  int ibar = 0;
69  if (fabs(tanTheta) > corner)
70  {
71  tSOSDest = ivProp_->propagate(ftsCurrent, *cylinder);
72  // std::cout<<" Propagate to cylinder "<<std::endl;
73  }
74  else if(tanTheta > 0.)
75  {tSOSDest = ivProp_->propagate(ftsCurrent, *forwardEndcap); ibar=1; }
76  else
77  {tSOSDest = ivProp_->propagate(ftsCurrent, *backwardEndcap); ibar=-1; }
78 
79 // std::cout<<" Trajectory valid? "<<tSOSDest.isValid()<<" First propagation in "<<ibar<<std::endl;
80 
81  if(! tSOSDest.isValid() )
82  {
83 // barrel
84  if(ibar == 0){
85  if (tanTheta < 0 ) tSOSDest = ivProp_->propagate( ftsCurrent,*forwardEndcap);
86  if (tanTheta >= 0 ) tSOSDest = ivProp_->propagate( ftsCurrent,*backwardEndcap);
87  }
88  else
89  {
90  tSOSDest = ivProp_->propagate(ftsCurrent, *cylinder);
91  }
92  } else
93  {
94 // missed target
95  if(abs(ibar) > 0)
96  {
97  if(tSOSDest.globalPosition().perp() > surface_iter->perp())
98  {
99  tSOSDest = ivProp_->propagate(ftsCurrent, *cylinder);
100  }
101  }
102  else
103  {
104  if (tanTheta < 0 ) tSOSDest = ivProp_->propagate( ftsCurrent,*forwardEndcap);
105  if (tanTheta >= 0 ) tSOSDest = ivProp_->propagate( ftsCurrent,*backwardEndcap);
106  }
107  }
108 
109 
110  // If missed the target, propagate to again
111 // if ((!tSOSDest.isValid()) && point.perp() > surface_iter->perp())
112 // {tSOSDest = ivProp_->propagate(ftsCurrent, *cylinder);std::cout<<" Propagate again 1 "<<std::endl;}
113 // std::cout<<" Track is ok after repropagation to cylinder or not? "<<tSOSDest.isValid()<<std::endl;
114 // if ((!tSOSDest.isValid()) && ftsStart.momentum().eta()>0. && fabs(ftsStart.momentum().eta())>1.)
115 // {tSOSDest = ivProp_->propagate(ftsStart, *forwardEndcap);std::cout<<" Propagate again 2 "<<std::endl;}
116 // std::cout<<" Track is ok after repropagation forward or not? "<<tSOSDest.isValid()<<std::endl;
117 // if ((!tSOSDest.isValid()) && ftsStart.momentum().eta()<0.&&fabs(ftsStart.momentum().eta())>1.)
118 // {tSOSDest = ivProp_->propagate(ftsStart, *backwardEndcap);std::cout<<" Propagate again 3 "<<std::endl;}
119 // std::cout<<" Track is after repropagation backward ok or not? "<<tSOSDest.isValid()<<std::endl;
120 
121 
122  if (! tSOSDest.isValid()) return trajectory;
123 
124 // std::cout<<" Propagate reach something"<<std::endl;
125  LogTrace("SuccessfullPropagation") << "Great, I reached something." << "\n"
126  << "\tx: " << tSOSDest.freeState()->position().x() << "\n"
127  << "\ty: " << tSOSDest.freeState()->position().y() << "\n"
128  << "\tz: " << tSOSDest.freeState()->position().z() << "\n"
129  << "\teta: " << tSOSDest.freeState()->position().eta() << "\n"
130  << "\tphi: " << tSOSDest.freeState()->position().phi() << "\n";
131 
132 // std::cout<<" The position of trajectory "<<tSOSDest.freeState()->position().perp()<<" "<<tSOSDest.freeState()->position().z()<<std::endl;
133 
134  GlobalPoint point = tSOSDest.freeState()->position();
135  point = tSOSDest.freeState()->position();
136  ftsCurrent = *tSOSDest.freeState();
137  trajectory.push_back(point);
138  }
139  return trajectory;
140 }
virtual void check_setup()
virtual FreeTrajectoryState propagate(const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest) const final
Definition: Propagator.h:119
T perp() const
Definition: PV3DBase.h:72
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
Definition: Plane.h:17
T mag() const
Definition: PV3DBase.h:67
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GlobalVector momentum() const
#define LogTrace(id)
Propagator * ivProp_
GlobalPoint position() const
Point3DBase< float, GlobalTag > PositionType
T eta() const
Definition: PV3DBase.h:76
TkRotation< float > RotationType
T x() const
Definition: PV3DBase.h:62
*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
int HDetIdAssociator::iEta ( const GlobalPoint point)
virtual

Definition at line 197 of file DetIdAssociator.cc.

References PV3DBase< T, PVType, FrameType >::eta(), etaBinSize_, and nEta_.

Referenced by buildMap(), getDetIdsCloseToAPoint(), and getDetIdsInACone().

198 {
199 // unequal bin sizes for endcap, following HCAL geometry
200  int iEta1 = int(point.eta()/etaBinSize_ + nEta_/2);
201  if (point.eta()>1.827 && point.eta()<=1.830) return iEta1-1;
202  else if (point.eta()>1.914 && point.eta()<=1.930) return iEta1-1;
203  else if (point.eta()>2.001 && point.eta()<=2.043) return iEta1-1;
204  else if (point.eta()>2.088 && point.eta()<=2.172) return iEta1-1;
205  else if (point.eta()>2.175 && point.eta()<=2.262) return iEta1-1;
206  else if (point.eta()>2.262 && point.eta()<=2.332) return iEta1-2;
207  else if (point.eta()>2.332 && point.eta()<=2.349) return iEta1-1;
208  else if (point.eta()>2.349 && point.eta()<=2.436) return iEta1-2;
209  else if (point.eta()>2.436 && point.eta()<=2.500) return iEta1-3;
210  else if (point.eta()>2.500 && point.eta()<=2.523) return iEta1-2;
211  else if (point.eta()>2.523 && point.eta()<=2.610) return iEta1-3;
212  else if (point.eta()>2.610 && point.eta()<=2.650) return iEta1-4;
213  else if (point.eta()>2.650 && point.eta()<=2.697) return iEta1-3;
214  else if (point.eta()>2.697 && point.eta()<=2.784) return iEta1-4;
215  else if (point.eta()>2.784 && point.eta()<=2.868) return iEta1-5;
216  else if (point.eta()>2.868 && point.eta()<=2.871) return iEta1-4;
217  else if (point.eta()>2.871 && point.eta()<=2.958) return iEta1-5;
218  else if (point.eta()>2.958) return iEta1-6;
219  else if (point.eta()<-1.827 && point.eta()>=-1.830) return iEta1+1;
220  else if (point.eta()<-1.914 && point.eta()>=-1.930) return iEta1+1;
221  else if (point.eta()<-2.001 && point.eta()>=-2.043) return iEta1+1;
222  else if (point.eta()<-2.088 && point.eta()>=-2.172) return iEta1+1;
223  else if (point.eta()<-2.175 && point.eta()>=-2.262) return iEta1+1;
224  else if (point.eta()<-2.262 && point.eta()>=-2.332) return iEta1+2;
225  else if (point.eta()<-2.332 && point.eta()>=-2.349) return iEta1+1;
226  else if (point.eta()<-2.349 && point.eta()>=-2.436) return iEta1+2;
227  else if (point.eta()<-2.436 && point.eta()>=-2.500) return iEta1+3;
228  else if (point.eta()<-2.500 && point.eta()>=-2.523) return iEta1+2;
229  else if (point.eta()<-2.523 && point.eta()>=-2.610) return iEta1+3;
230  else if (point.eta()<-2.610 && point.eta()>=-2.650) return iEta1+4;
231  else if (point.eta()<-2.650 && point.eta()>=-2.697) return iEta1+3;
232  else if (point.eta()<-2.697 && point.eta()>=-2.784) return iEta1+4;
233  else if (point.eta()<-2.784 && point.eta()>=-2.868) return iEta1+5;
234  else if (point.eta()<-2.868 && point.eta()>=-2.871) return iEta1+4;
235  else if (point.eta()<-2.871 && point.eta()>=-2.958) return iEta1+5;
236  else if (point.eta()<-2.349) return iEta1+6;
237  else return iEta1;
238 }
const double etaBinSize_
T eta() const
Definition: PV3DBase.h:76
virtual bool HDetIdAssociator::insideElement ( const GlobalPoint ,
const DetId  
)
protectedpure virtual

Implemented in HCaloDetIdAssociator.

Referenced by getCrossedDetIds().

int HDetIdAssociator::iPhi ( const GlobalPoint point)
virtual

Definition at line 241 of file DetIdAssociator.cc.

References nPhi_, PV3DBase< T, PVType, FrameType >::phi(), and pi.

Referenced by buildMap(), getDetIdsCloseToAPoint(), and getDetIdsInACone().

242 {
243  double pi=4*atan(1.);
244  int iPhi1 = int((double(point.phi())+pi)/(2*pi)*nPhi_);
245  return iPhi1;
246 }
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
const Double_t pi
virtual bool HDetIdAssociator::nearElement ( const GlobalPoint point,
const DetId id,
const double  distance 
)
inlineprotectedvirtual

Definition at line 84 of file DetIdAssociator.h.

References HLT_25ns14e33_v1_cff::distance, PV3DBase< T, PVType, FrameType >::eta(), getPosition(), PV3DBase< T, PVType, FrameType >::phi(), funct::pow(), and mathSSE::sqrt().

Referenced by getDetIdsInACone().

85  {
86  GlobalPoint center = getPosition(id);
87  return sqrt(pow(point.eta()-center.eta(),2)+pow(point.phi()-center.phi(),2)) < distance;
88  };
virtual GlobalPoint getPosition(const DetId &)=0
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T sqrt(T t)
Definition: SSEVec.h:48
T eta() const
Definition: PV3DBase.h:76
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual void HDetIdAssociator::setPropagator ( Propagator ptr)
inlinevirtual

Definition at line 67 of file DetIdAssociator.h.

References ivProp_.

67 { ivProp_ = ptr; };
Propagator * ivProp_

Member Data Documentation

const double HDetIdAssociator::etaBinSize_
protected

Definition at line 93 of file DetIdAssociator.h.

Referenced by buildMap(), check_setup(), getDetIdsCloseToAPoint(), and iEta().

Propagator* HDetIdAssociator::ivProp_
protected

Definition at line 94 of file DetIdAssociator.h.

Referenced by check_setup(), getTrajectory(), and setPropagator().

const int HDetIdAssociator::nEta_
protected
const int HDetIdAssociator::nPhi_
protected
std::vector<std::vector<std::set<DetId> > >* HDetIdAssociator::theMap_
protected

Definition at line 88 of file DetIdAssociator.h.

Referenced by buildMap(), getDetIdsCloseToAPoint(), and getDetIdsInACone().