CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HGCNumberingScheme.cc
Go to the documentation of this file.
1 // File: HGCNumberingScheme.cc
3 // Description: Numbering scheme for High Granularity Calorimeter
8 
9 #include "CLHEP/Units/GlobalSystemOfUnits.h"
10 #include <iostream>
11 
13  std::string & name, bool check,
14  int verbose) :
15  CaloNumberingScheme(0), check_(check), verbosity(verbose),
16  hgcons(new HGCalDDDConstants(cpv,name)) {
17  edm::LogInfo("HGCSim") << "Creating HGCNumberingScheme for " << name;
18 }
19 
21  edm::LogInfo("HGCSim") << "Deleting HGCNumberingScheme";
22 }
23 
24 //
25 uint32_t HGCNumberingScheme::getUnitID(ForwardSubdetector subdet, int layer, int sector, int iz, const G4ThreeVector &pos) {
26 
27  std::pair<int,int> phicell = hgcons->assignCell(pos.x(),pos.y(),layer,0,false);
28  int phiSector = phicell.first;
29  int icell = phicell.second;
30 
31  //build the index
32  uint32_t index = HGCalDetId(subdet,iz,layer,sector,phiSector,icell).rawId();
33 
34  //check if it fits
35  if ((!HGCalDetId::isValid(subdet,iz,layer,sector,phiSector,icell)) ||
36  (!hgcons->isValid(layer,sector,icell,false))) {
37  index = 0;
38  if (check_ && icell != -1) {
39  edm::LogError("HGCSim") << "[HGCNumberingScheme] ID out of bounds :"
40  << " Subdet= " << subdet << " Zside= " << iz
41  << " Layer= " << layer << " Sector= " << sector
42  << " SubSector= " << phiSector << " Cell= "
43  << icell << " Local position: (" << pos.x()
44  << "," << pos.y() << "," << pos.z() << ")";
45  }
46  }
47  if (verbosity > 0)
48  std::cout << "HGCNumberingScheme::i/p " << subdet << ":" << layer << ":"
49  << sector << ":" << iz << ":" << pos << " o/p " << phiSector
50  << ":" << icell << ":" << std::hex << index << std::dec << " "
51  << HGCalDetId(index) << std::endl;
52  return index;
53 }
54 
55 //
56 int HGCNumberingScheme::assignCell(float x, float y, int layer) {
57 
58  std::pair<int,int> phicell = hgcons->assignCell(x,y,layer,0,false);
59  return phicell.second;
60 }
61 
62 //
63 std::pair<float,float> HGCNumberingScheme::getLocalCoords(int cell, int layer){
64 
65  return hgcons->locateCell(cell,layer,0,false);
66 }
bool isValid(int lay, int mod, int cell, bool reco) const
int assignCell(float x, float y, int layer)
maps a hit position to a sequential cell in a trapezoid surface defined by h,b,t
std::pair< float, float > getLocalCoords(int cell, int layer)
inverts the cell number in a trapezoid surface to local coordinates
virtual uint32_t getUnitID(ForwardSubdetector subdet, int layer, int module, int iz, const G4ThreeVector &pos)
assigns the det id to a hit
ForwardSubdetector
type of data representation of DDCompactView
Definition: DDCompactView.h:77
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
T x() const
Cartesian x coordinate.
bool check(const std::string &)
HGCalDDDConstants * hgcons
static bool isValid(ForwardSubdetector subdet, int zp, int lay, int mod, int subsec, int cell)
Definition: HGCalDetId.cc:47
std::pair< int, int > assignCell(float x, float y, int lay, int subSec, bool reco) const
tuple cout
Definition: gather_cfg.py:121
std::pair< float, float > locateCell(int cell, int lay, int subSec, bool reco) const