CMS 3D CMS Logo

MTDCPEBase.cc
Go to the documentation of this file.
4 
6 
7 // MessageLogger
9 
10 // Magnetic field
12 
13 #include <iostream>
14 
15 using namespace std;
16 
17 //-----------------------------------------------------------------------------
18 // A constructor run for generic and templates
19 //
20 //-----------------------------------------------------------------------------
21 MTDCPEBase::MTDCPEBase(edm::ParameterSet const& conf, const MTDGeometry& geom) : geom_(geom) { fillDetParams(); }
22 
23 //-----------------------------------------------------------------------------
24 // Fill all variables which are constant for an event (geometry)
25 //-----------------------------------------------------------------------------
27  auto const& dus = geom_.detUnits();
28  unsigned detectors = dus.size();
29  m_DetParams.resize(detectors);
30  LogDebug("MTDCPEBase::fillDetParams():") << "caching " << detectors << "MTD detectors" << endl;
31  for (unsigned i = 0; i != detectors; ++i) {
32  auto& p = m_DetParams[i];
33  p.theDet = dynamic_cast<const MTDGeomDetUnit*>(dus[i]);
34  assert(p.theDet);
35 
36  p.theOrigin = p.theDet->surface().toLocal(GlobalPoint(0, 0, 0));
37 
38  //--- p.theDet->type() returns a GeomDetType, which implements subDetector()
39  p.thePart = p.theDet->type().subDetector();
40 
41  //--- bounds() is implemented in BoundSurface itself.
42  p.theThickness = p.theDet->surface().bounds().thickness();
43 
44  // Cache the det id for templates and generic erros
45  p.theTopol = &(static_cast<const ProxyMTDTopology&>(p.theDet->topology()));
46  assert(p.theTopol);
47  p.theRecTopol = &(static_cast<const RectangularMTDTopology&>(p.theTopol->specificTopology()));
48  assert(p.theRecTopol);
49 
50  //--- The geometrical description of one module/plaquette
51  std::pair<float, float> pitchxy = p.theRecTopol->pitch();
52  p.thePitchX = pitchxy.first; // pitch along x
53  p.thePitchY = pitchxy.second; // pitch along y
54 
55  LogDebug("MTDCPEBase::fillDetParams()") << "***** MTD LAYOUT *****"
56  << " thePart = " << p.thePart << " theThickness = " << p.theThickness
57  << " thePitchX = " << p.thePitchX << " thePitchY = " << p.thePitchY;
58  }
59 }
60 
61 //-----------------------------------------------------------------------------
62 // One function to cache the variables common for one DetUnit.
63 //-----------------------------------------------------------------------------
65 
66 //------------------------------------------------------------------------
67 MTDCPEBase::DetParam const& MTDCPEBase::detParam(const GeomDetUnit& det) const { return m_DetParams.at(det.index()); }
68 
70  //remember measurement point is row(col)+0.5f
72  return dp.theTopol->localPosition(pos);
73 }
74 
76  constexpr double one_over_twelve = 1. / 12.;
78  MeasurementError simpleRect(one_over_twelve, 0, one_over_twelve);
79  return dp.theTopol->localError(pos, simpleRect);
80 }
81 
83  return cp.theCluster->time();
84 }
85 
87  return cp.theCluster->timeError();
88 }
#define LogDebug(id)
virtual LocalError localError(DetParam const &dp, ClusterParam &cp) const
Definition: MTDCPEBase.cc:75
LocalPoint localPosition(const MeasurementPoint &) const override
float y() const
Definition: FTLCluster.h:126
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
Definition: MTDGeometry.h:27
LocalError localError(const MeasurementPoint &, const MeasurementError &) const override
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const FTLCluster * theCluster
Definition: MTDCPEBase.h:51
float time() const
Definition: FTLCluster.h:131
void fillDetParams()
Definition: MTDCPEBase.cc:26
virtual TimeValue clusterTime(DetParam const &dp, ClusterParam &cp) const
Definition: MTDCPEBase.cc:82
float timeError() const
Definition: FTLCluster.h:136
float x() const
Definition: FTLCluster.h:121
def detectors(dt=True, csc=True, me42=False, chambers=True, superlayers=False, layers=False)
DetParams m_DetParams
Definition: MTDCPEBase.h:104
virtual TimeValueError clusterTimeError(DetParam const &dp, ClusterParam &cp) const
Definition: MTDCPEBase.cc:86
int index() const
Definition: GeomDet.h:83
DetParam const & detParam(const GeomDetUnit &det) const
Definition: MTDCPEBase.cc:67
const ProxyMTDTopology * theTopol
Definition: MTDCPEBase.h:36
virtual LocalPoint localPosition(DetParam const &dp, ClusterParam &cp) const
Definition: MTDCPEBase.cc:69
void setTheClu(DetParam const &dp, ClusterParam &cp) const
Definition: MTDCPEBase.cc:64
const MTDGeometry & geom_
Definition: MTDCPEBase.h:92
#define constexpr
MTDCPEBase(edm::ParameterSet const &conf, const MTDGeometry &geom)
Definition: MTDCPEBase.cc:21