CMS 3D CMS Logo

HeterogeneousHGCalHEFCellPositionsFiller.cc
Go to the documentation of this file.
2 
5  edm::ESInputTag{"", "HGCalHESiliconSensitive"});
7 }
8 
10 
11 //the geometry is not required if the layer offset is hardcoded (potential speed-up)
13  //fill the CPU position structure from the geometry
14  posmap_->zLayer.clear();
15  posmap_->nCellsLayer.clear();
16  posmap_->nCellsWaferUChunk.clear();
17  posmap_->nCellsHexagon.clear();
18  posmap_->detid.clear();
19 
20  int nlayers = ddd_->lastLayer(true) - ddd_->firstLayer() + 1;
21  int upper_estimate_wafer_number_1D = 2 * nlayers * (ddd_->waferMax() - ddd_->waferMin());
22  int upper_estimate_wafer_number_2D = upper_estimate_wafer_number_1D * (ddd_->waferMax() - ddd_->waferMin());
23  int upper_estimate_cell_number = upper_estimate_wafer_number_2D * 3 * 12 * 12;
24  posmap_->zLayer.resize(nlayers * 2);
25  posmap_->nCellsLayer.reserve(nlayers * 2);
26  posmap_->nCellsWaferUChunk.reserve(upper_estimate_wafer_number_1D);
27  posmap_->nCellsHexagon.reserve(upper_estimate_wafer_number_2D);
28  posmap_->detid.reserve(upper_estimate_cell_number);
29  //set position-related variables
30  posmap_->waferSize = static_cast<float>(params_->waferSize_);
31  posmap_->sensorSeparation = static_cast<float>(params_->sensorSeparation_);
33  assert(posmap_->firstLayer == 1); //otherwise the loop over the layers has to be changed
34  posmap_->lastLayer = ddd_->lastLayer(true);
37 
38  unsigned sumCellsLayer, sumCellsWaferUChunk;
39 
40  //store detids following a geometry ordering
41  for (int ilayer = 1; ilayer <= posmap_->lastLayer; ++ilayer) {
42  sumCellsLayer = 0;
43  posmap_->zLayer[ilayer - 1] = static_cast<float>(ddd_->waferZ(ilayer, true)); //originally a double
44  posmap_->zLayer[ilayer - 1 + nlayers] = static_cast<float>(ddd_->waferZ(ilayer, true)); //originally a double
45 
46  for (int iwaferU = posmap_->waferMin; iwaferU < posmap_->waferMax; ++iwaferU) {
47  sumCellsWaferUChunk = 0;
48 
49  for (int iwaferV = posmap_->waferMin; iwaferV < posmap_->waferMax; ++iwaferV) {
50  //0: fine; 1: coarseThin; 2: coarseThick (as defined in DataFormats/ForwardDetId/interface/HGCSiliconDetId.h)
51  int type_ = ddd_->waferType(ilayer, iwaferU, iwaferV, false);
52 
53  int nCellsHexSide = ddd_->numberCellsHexagon(ilayer, iwaferU, iwaferV, false);
54  int nCellsHexTotal = ddd_->numberCellsHexagon(ilayer, iwaferU, iwaferV, true);
55  sumCellsLayer += nCellsHexTotal;
56  sumCellsWaferUChunk += nCellsHexTotal;
57  posmap_->nCellsHexagon.push_back(nCellsHexTotal);
58 
59  //left side of wafer
60  for (int cellUmax = nCellsHexSide, icellV = 0; cellUmax < 2 * nCellsHexSide and icellV < nCellsHexSide;
61  ++cellUmax, ++icellV) {
62  for (int icellU = 0; icellU <= cellUmax; ++icellU) {
63  HGCSiliconDetId detid_(DetId::HGCalHSi, 1, type_, ilayer, iwaferU, iwaferV, icellU, icellV);
64  posmap_->detid.push_back(detid_.rawId());
65  }
66  }
67  //right side of wafer
68  for (int cellUmin = 1, icellV = nCellsHexSide; cellUmin <= nCellsHexSide and icellV < 2 * nCellsHexSide;
69  ++cellUmin, ++icellV) {
70  for (int icellU = cellUmin; icellU < 2 * nCellsHexSide; ++icellU) {
71  HGCSiliconDetId detid_(DetId::HGCalHSi, 1, type_, ilayer, iwaferU, iwaferV, icellU, icellV);
72  posmap_->detid.push_back(detid_.rawId());
73  }
74  }
75  }
76  posmap_->nCellsWaferUChunk.push_back(sumCellsWaferUChunk);
77  }
78  posmap_->nCellsLayer.push_back(sumCellsLayer);
79  }
80 }
81 
82 std::unique_ptr<HeterogeneousHGCalHEFCellPositionsConditions> HeterogeneousHGCalHEFCellPositionsFiller::produce(
84  auto geom = iRecord.getTransientHandle(geometryToken_);
85  ddd_ = &(geom->topology().dddConstants());
87 
89 
90  std::unique_ptr<HeterogeneousHGCalHEFCellPositionsConditions> up =
91  std::make_unique<HeterogeneousHGCalHEFCellPositionsConditions>(posmap_);
92  return up;
93 }
94 
ESTransientHandle< ProductT > getTransientHandle(ESGetToken< ProductT, DepRecordT > const &iToken) const
double waferZ(int layer, bool reco) const
Definition: BitonicSort.h:7
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:165
const HGCalParameters * getParameter() const
int lastLayer(bool reco) const
int firstLayer() const
assert(be >=bs)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
hgcal_conditions::positions::HGCalPositionsMapping * posmap_
int numberCellsHexagon(int wafer) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int waferType(DetId const &id, bool fromFile) const
std::unique_ptr< HeterogeneousHGCalHEFCellPositionsConditions > produce(const HeterogeneousHGCalHEFCellPositionsConditionsRecord &)
edm::ESGetToken< HGCalGeometry, IdealGeometryRecord > geometryToken_