CMS 3D CMS Logo

FastTimeGeometry.cc
Go to the documentation of this file.
7 
8 #include <cmath>
9 
10 #include <Math/Transform3D.h>
11 #include <Math/EulerAngles.h>
12 
14 typedef std::vector<float> ParmVec;
15 
16 //#define EDM_ML_DEBUG
17 
19  : m_topology(topology_),
20  m_cellVec(topology_.totalGeomModules()),
21  m_validGeomIds(topology_.totalGeomModules()),
22  m_Type(topology_.detectorType()),
23  m_subdet(topology_.subDetector()) {
24  m_validIds.reserve(topology().totalModules());
25 #ifdef EDM_ML_DEBUG
26  edm::LogVerbatim("FastTimeGeom") << "Expected total # of Geometry Modules " << topology().totalGeomModules();
27 #endif
28 }
29 
31 
33 
35 
36 void FastTimeGeometry::localCorners(Pt3DVec& lc, const CCGFloat* pv, unsigned int i, Pt3D& ref) {
37  FlatTrd::localCorners(lc, pv, ref);
38 }
39 
41  const GlobalPoint& f1, const GlobalPoint& f2, const GlobalPoint& f3, const CCGFloat* parm, const DetId& detId) {
43  DetId geomId = (DetId)(FastTimeDetId(detId).geometryCell());
44  int nEtaZ = topology().dddConstants().numberEtaZ(m_Type);
46 
47  const uint32_t cellIndex(topology().detId2denseGeomId(detId));
48 
49  m_cellVec.at(cellIndex) = FlatTrd(cornersMgr(), f1, f2, f3, parm);
50  m_validGeomIds.at(cellIndex) = geomId;
51 
52 #ifdef EDM_ML_DEBUG
53  unsigned int nOld = m_validIds.size();
54 #endif
55  for (int etaZ = 1; etaZ <= nEtaZ; ++etaZ) {
56  id.iEtaZ = etaZ;
57  for (int phi = 1; phi <= nPhi; ++phi) {
58  id.iPhi = phi;
59  DetId idc = topology().encode(id);
60  if (topology().valid(idc)) {
61  m_validIds.emplace_back(idc);
62  }
63  }
64  }
65 
66 #ifdef EDM_ML_DEBUG
67  edm::LogVerbatim("FastTimeGeom") << "FastTimeGeometry::newCell-> [" << cellIndex << "] front:" << f1.x() << '/'
68  << f1.y() << '/' << f1.z() << " back:" << f2.x() << '/' << f2.y() << '/' << f2.z()
69  << " eta|phi " << m_cellVec[cellIndex].etaPos() << ":"
70  << m_cellVec[cellIndex].phiPos() << " id:" << FastTimeDetId(detId)
71  << " with valid DetId from " << nOld << " to " << m_validIds.size();
72  edm::LogVerbatim("FastTimeGeom") << "Cell[" << cellIndex << "] " << std::hex << geomId.rawId() << ":"
73  << m_validGeomIds[cellIndex].rawId() << std::dec;
74 #endif
75 }
76 
77 std::shared_ptr<const CaloCellGeometry> FastTimeGeometry::getGeometry(const DetId& id) const {
78  if (id == DetId())
79  return nullptr; // nothing to get
80  DetId geoId = (DetId)(FastTimeDetId(id).geometryCell());
81  const uint32_t cellIndex(topology().detId2denseGeomId(geoId));
82  const GlobalPoint pos = (id != geoId) ? getPosition(id) : GlobalPoint();
83  return cellGeomPtr(cellIndex, pos);
84 }
85 
86 bool FastTimeGeometry::present(const DetId& id) const {
87  if (id == DetId())
88  return false;
89  DetId geoId = (DetId)(FastTimeDetId(id).geometryCell());
90  const uint32_t index(topology().detId2denseGeomId(geoId));
91  return (nullptr != getGeometryRawPtr(index));
92 }
93 
95  FastTimeDetId id_ = FastTimeDetId(id);
96  auto pos = topology().dddConstants().getPosition(m_Type, id_.ieta(), id_.iphi(), id_.zside());
97  return GlobalPoint(0.1 * pos.x(), 0.1 * pos.y(), 0.1 * pos.z());
98 }
99 
101  FastTimeDetId id_ = FastTimeDetId(id);
102  auto corners = topology().dddConstants().getCorners(m_Type, id_.ieta(), id_.iphi(), id_.zside());
104  for (const auto& corner : corners) {
105  out.emplace_back(0.1 * corner.x(), 0.1 * corner.y(), 0.1 * corner.z());
106  }
107  return out;
108 }
109 
111  int zside = (r.z() > 0) ? 1 : -1;
112  std::pair<int, int> etaZPhi;
113  if (m_Type == 1) {
114  double zz = (zside > 0) ? r.z() : -r.z();
115  etaZPhi = topology().dddConstants().getZPhi(zz, r.phi());
116  } else {
117  double phi = (zside > 0) ? static_cast<double>(r.phi()) : atan2(r.y(), -r.x());
118  // Cast needed to resolve compile-time ambiguity of ? operator between
119  // convertible Phi class and atan2 template function.
120 
121  etaZPhi = topology().dddConstants().getEtaPhi(r.perp(), phi);
122  }
123  FastTimeDetId id = FastTimeDetId(m_Type, etaZPhi.first, etaZPhi.second, zside);
124 #ifdef EDM_ML_DEBUG
125  edm::LogVerbatim("FastTimeGeom") << "getClosestCell: for (" << r.x() << ", " << r.y() << ", " << r.z() << ") Id "
126  << id.type() << ":" << id.zside() << ":" << id.ieta() << ":" << id.iphi();
127 #endif
128 
129  return (topology().valid(id) ? DetId(id) : DetId());
130 }
131 
134  return dss;
135 }
136 
138  if (m_Type == 1)
139  return "FastTimeBarrel";
140  else if (m_Type == 2)
141  return "FastTimeEndcap";
142  else
143  return "Unknown";
144 }
145 
146 unsigned int FastTimeGeometry::indexFor(const DetId& id) const {
147  unsigned int cellIndex = m_cellVec.size();
148  if (id != DetId()) {
149  DetId geoId = (DetId)(FastTimeDetId(id).geometryCell());
150  cellIndex = topology().detId2denseGeomId(geoId);
151 #ifdef EDM_ML_DEBUG
152  edm::LogVerbatim("FastTimeGeom") << "indexFor " << std::hex << id.rawId() << ":" << geoId.rawId() << std::dec
153  << " index " << cellIndex;
154 #endif
155  }
156  return cellIndex;
157 }
158 
160 
162  // Modify the RawPtr class
163  const CaloCellGeometry* cell(&m_cellVec[index]);
164  return (m_cellVec.size() < index || nullptr == cell->param() ? nullptr : cell);
165 }
166 
167 std::shared_ptr<const CaloCellGeometry> FastTimeGeometry::cellGeomPtr(uint32_t index) const {
168  if ((index >= m_cellVec.size()) || (m_validGeomIds[index].rawId() == 0))
169  return nullptr;
170  static const auto do_not_delete = [](const void*) {};
171  auto cell = std::shared_ptr<const CaloCellGeometry>(&m_cellVec[index], do_not_delete);
172  if (nullptr == cell->param())
173  return nullptr;
174  return cell;
175 }
176 
177 std::shared_ptr<const CaloCellGeometry> FastTimeGeometry::cellGeomPtr(uint32_t index, const GlobalPoint& pos) const {
178  if ((index >= m_cellVec.size()) || (m_validGeomIds[index].rawId() == 0))
179  return nullptr;
180  if (pos == GlobalPoint())
181  return cellGeomPtr(index);
182  auto cell = std::make_shared<FlatTrd>(m_cellVec[index]);
183  cell->setPosition(pos);
184 #ifdef EDM_ML_DEBUG
185  edm::LogVerbatim("FastTimeGeomX") << "cellGeomPtr " << pos << ":" << cell;
186 #endif
187  if (nullptr == cell->param())
188  return nullptr;
189  return cell;
190 }
191 
193  edm::LogError("FastTimeGeom") << "FastTimeGeometry::addValidID is not implemented";
194 }
195 
196 // FIXME: Change sorting algorithm if needed
197 namespace {
198  struct rawIdSort {
199  bool operator()(const DetId& a, const DetId& b) { return (a.rawId() < b.rawId()); }
200  };
201 } // namespace
202 
204  m_validIds.shrink_to_fit();
205  std::sort(m_validIds.begin(), m_validIds.end(), rawIdSort());
206 }
207 
211  CaloSubdetectorGeometry::IVec& dinsVector) const {
212  unsigned int numberOfCells = topology().totalGeomModules(); // total Geom Modules both sides
215 
216  trVector.reserve(numberOfCells * numberOfTransformParms());
217  iVector.reserve(numberOfCells);
218  dimVector.reserve(numberOfShapes * numberOfParametersPerShape);
219  dinsVector.reserve(numberOfCells);
220 
221  for (unsigned int k = 0; k < topology().totalGeomModules(); ++k) {
224  params[1] = params[2] = 0;
228  params[6] = params[10] = 0;
229  params[11] = (k == 0) ? 1.0 : -1.0;
230  dimVector.insert(dimVector.end(), params.begin(), params.end());
231  }
232 
233  for (unsigned int i(0); i < numberOfCells; ++i) {
234  DetId detId = m_validGeomIds[i];
235  dinsVector.emplace_back(topology().detId2denseGeomId(detId));
236  iVector.emplace_back(1);
237 
238  Tr3D tr;
239  auto ptr(cellGeomPtr(i));
240  if (nullptr != ptr) {
241  ptr->getTransform(tr, (Pt3DVec*)nullptr);
242 
243  if (Tr3D() == tr) { // there is no rotation
244  const GlobalPoint& gp(ptr->getPosition());
245  tr = HepGeom::Translate3D(gp.x(), gp.y(), gp.z());
246  }
247 
248  const CLHEP::Hep3Vector tt(tr.getTranslation());
249  trVector.emplace_back(tt.x());
250  trVector.emplace_back(tt.y());
251  trVector.emplace_back(tt.z());
252  if (6 == numberOfTransformParms()) {
253  const CLHEP::HepRotation rr(tr.getRotation());
254  const ROOT::Math::Transform3D rtr(
255  rr.xx(), rr.xy(), rr.xz(), tt.x(), rr.yx(), rr.yy(), rr.yz(), tt.y(), rr.zx(), rr.zy(), rr.zz(), tt.z());
257  rtr.GetRotation(ea);
258  trVector.emplace_back(ea.Phi());
259  trVector.emplace_back(ea.Theta());
260  trVector.emplace_back(ea.Psi());
261  }
262  }
263  }
264 }
265 
267 
Log< level::Info, true > LogVerbatim
double getRin(int type) const
void localCorners(Pt3DVec &lc, const CCGFloat *pv, unsigned int i, Pt3D &ref)
std::set< DetId > DetIdSet
bool present(const DetId &id) const override
is this detid present in the geometry?
double getZHalf(int type) const
A base class to handle the particular shape of HGCal volumes.
Definition: FlatTrd.h:19
unsigned int indexFor(const DetId &id) const override
std::vector< CCGFloat > DimVec
int numberEtaZ(int type) const
std::vector< GlobalPoint > CornersVec
unsigned int totalGeomModules() const
CornersVec getCorners(const DetId &id) const
Returns the corner points of this cell&#39;s volume.
CaloCellGeometry::Pt3DVec Pt3DVec
virtual unsigned int numberOfTransformParms() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< unsigned int > IVec
virtual void fillNamedParams(DDFilteredView fv)
void addValidID(const DetId &id)
std::vector< float > ParmVec
std::vector< CCGFloat > TrVec
HepGeom::Transform3D Tr3D
virtual unsigned int numberOfShapes() const
int zside(DetId const &)
Log< level::Error, false > LogError
CaloCellGeometry::CCGFloat CCGFloat
CaloCellGeometry::Pt3D Pt3D
FastTimeGeometry(const FastTimeTopology &topology)
double getRout(int type) const
CaloCellGeometry::Tr3D Tr3D
int iphi() const
get the absolute value of the cell #&#39;s along y-axis (EC) | phi (Barrel)
Definition: FastTimeDetId.h:45
~FastTimeGeometry() override
std::vector< DetId > m_validIds
def pv(vc)
Definition: MetAnalyzer.py:7
int ieta() const
get the absolute value of the cell #&#39;s along x-axis (EC) | z-axis (Barel)
Definition: FastTimeDetId.h:41
std::pair< int, int > getZPhi(double z, double phi) const
std::string cellElement() const
const FastTimeDDDConstants & dddConstants() const
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
Definition: FlatTrd.cc:153
DetIdSet getCells(const GlobalPoint &r, double dR) const override
Get a list of all cells within a dR of the given cell.
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
std::shared_ptr< const CaloCellGeometry > cellGeomPtr(uint32_t index) const override
Definition: DetId.h:17
int zside() const
get the z-side of the cell (1/-1)
Definition: FastTimeDetId.h:48
AlgebraicVector EulerAngles
Definition: Definitions.h:34
std::pair< int, int > getEtaPhi(double r, double phi) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void getSummary(CaloSubdetectorGeometry::TrVec &trVector, CaloSubdetectorGeometry::IVec &iVector, CaloSubdetectorGeometry::DimVec &dimVector, CaloSubdetectorGeometry::IVec &dinsVector) const override
void initializeParms() override
CaloCellGeometry::CornersMgr * cornersMgr()
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const override
Get the cell geometry of a given detector id. Should return false if not found.
const FastTimeTopology & topology() const
double b
Definition: hdecay.h:118
GlobalPoint getPosition(int type, int izeta, int iphi, int zside) const
DecodedDetId decode(const DetId &id) const
virtual uint32_t detId2denseGeomId(const DetId &id) const
unsigned int sizeForDenseIndex() const
void newCell(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId) override
DetId getClosestCell(const GlobalPoint &r) const override
CaloCellGeometry::Tr3D Tr3D
std::vector< GlobalPoint > getCorners(int type, int izeta, int iphi, int zside) const
double a
Definition: hdecay.h:119
GlobalPoint getPosition(const DetId &id) const
FastTimeDetId geometryCell() const
Definition: FastTimeDetId.h:32
const CCGFloat * param() const
std::vector< DetId > m_validGeomIds
DetId encode(const DecodedDetId &id_) const
const CaloCellGeometry * getGeometryRawPtr(uint32_t index) const override
virtual unsigned int numberOfParametersPerShape() const
int numberPhi(int type) const