CMS 3D CMS Logo

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