Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:54

0001 #include "RecoLocalCalo/HGCalRecProducers/plugins/HeterogeneousHGCalHEFCellPositionsFiller.h"
0002 
0003 HeterogeneousHGCalHEFCellPositionsFiller::HeterogeneousHGCalHEFCellPositionsFiller(const edm::ParameterSet& ps) {
0004   geometryToken_ = setWhatProduced(this).consumesFrom<HGCalGeometry, IdealGeometryRecord>(
0005       edm::ESInputTag{"", "HGCalHESiliconSensitive"});
0006   posmap_ = new hgcal_conditions::positions::HGCalPositionsMapping();
0007 }
0008 
0009 HeterogeneousHGCalHEFCellPositionsFiller::~HeterogeneousHGCalHEFCellPositionsFiller() { delete posmap_; }
0010 
0011 //the geometry is not required if the layer offset is hardcoded (potential speed-up)
0012 void HeterogeneousHGCalHEFCellPositionsFiller::set_conditions_() {
0013   //fill the CPU position structure from the geometry
0014   posmap_->zLayer.clear();
0015   posmap_->nCellsLayer.clear();
0016   posmap_->nCellsWaferUChunk.clear();
0017   posmap_->nCellsHexagon.clear();
0018   posmap_->detid.clear();
0019 
0020   int nlayers = ddd_->lastLayer(true) - ddd_->firstLayer() + 1;
0021   int upper_estimate_wafer_number_1D = 2 * nlayers * (ddd_->waferMax() - ddd_->waferMin());
0022   int upper_estimate_wafer_number_2D = upper_estimate_wafer_number_1D * (ddd_->waferMax() - ddd_->waferMin());
0023   int upper_estimate_cell_number = upper_estimate_wafer_number_2D * 3 * 12 * 12;
0024   posmap_->zLayer.resize(nlayers * 2);
0025   posmap_->nCellsLayer.reserve(nlayers * 2);
0026   posmap_->nCellsWaferUChunk.reserve(upper_estimate_wafer_number_1D);
0027   posmap_->nCellsHexagon.reserve(upper_estimate_wafer_number_2D);
0028   posmap_->detid.reserve(upper_estimate_cell_number);
0029   //set position-related variables
0030   posmap_->waferSize = static_cast<float>(params_->waferSize_);
0031   posmap_->sensorSeparation = static_cast<float>(params_->sensorSeparation_);
0032   posmap_->firstLayer = ddd_->firstLayer();
0033   assert(posmap_->firstLayer == 1);  //otherwise the loop over the layers has to be changed
0034   posmap_->lastLayer = ddd_->lastLayer(true);
0035   posmap_->waferMin = ddd_->waferMin();
0036   posmap_->waferMax = ddd_->waferMax();
0037 
0038   unsigned sumCellsLayer, sumCellsWaferUChunk;
0039 
0040   //store detids following a geometry ordering
0041   for (int ilayer = 1; ilayer <= posmap_->lastLayer; ++ilayer) {
0042     sumCellsLayer = 0;
0043     posmap_->zLayer[ilayer - 1] = static_cast<float>(ddd_->waferZ(ilayer, true));            //originally a double
0044     posmap_->zLayer[ilayer - 1 + nlayers] = static_cast<float>(ddd_->waferZ(ilayer, true));  //originally a double
0045 
0046     for (int iwaferU = posmap_->waferMin; iwaferU < posmap_->waferMax; ++iwaferU) {
0047       sumCellsWaferUChunk = 0;
0048 
0049       for (int iwaferV = posmap_->waferMin; iwaferV < posmap_->waferMax; ++iwaferV) {
0050         //0: fine; 1: coarseThin; 2: coarseThick (as defined in DataFormats/ForwardDetId/interface/HGCSiliconDetId.h)
0051         int type_ = ddd_->waferType(ilayer, iwaferU, iwaferV, false);
0052 
0053         int nCellsHexSide = ddd_->numberCellsHexagon(ilayer, iwaferU, iwaferV, false);
0054         int nCellsHexTotal = ddd_->numberCellsHexagon(ilayer, iwaferU, iwaferV, true);
0055         sumCellsLayer += nCellsHexTotal;
0056         sumCellsWaferUChunk += nCellsHexTotal;
0057         posmap_->nCellsHexagon.push_back(nCellsHexTotal);
0058 
0059         //left side of wafer
0060         for (int cellUmax = nCellsHexSide, icellV = 0; cellUmax < 2 * nCellsHexSide and icellV < nCellsHexSide;
0061              ++cellUmax, ++icellV) {
0062           for (int icellU = 0; icellU <= cellUmax; ++icellU) {
0063             HGCSiliconDetId detid_(DetId::HGCalHSi, 1, type_, ilayer, iwaferU, iwaferV, icellU, icellV);
0064             posmap_->detid.push_back(detid_.rawId());
0065           }
0066         }
0067         //right side of wafer
0068         for (int cellUmin = 1, icellV = nCellsHexSide; cellUmin <= nCellsHexSide and icellV < 2 * nCellsHexSide;
0069              ++cellUmin, ++icellV) {
0070           for (int icellU = cellUmin; icellU < 2 * nCellsHexSide; ++icellU) {
0071             HGCSiliconDetId detid_(DetId::HGCalHSi, 1, type_, ilayer, iwaferU, iwaferV, icellU, icellV);
0072             posmap_->detid.push_back(detid_.rawId());
0073           }
0074         }
0075       }
0076       posmap_->nCellsWaferUChunk.push_back(sumCellsWaferUChunk);
0077     }
0078     posmap_->nCellsLayer.push_back(sumCellsLayer);
0079   }
0080 }
0081 
0082 std::unique_ptr<HeterogeneousHGCalHEFCellPositionsConditions> HeterogeneousHGCalHEFCellPositionsFiller::produce(
0083     const HeterogeneousHGCalHEFCellPositionsConditionsRecord& iRecord) {
0084   auto geom = iRecord.getTransientHandle(geometryToken_);
0085   ddd_ = &(geom->topology().dddConstants());
0086   params_ = ddd_->getParameter();
0087 
0088   set_conditions_();
0089 
0090   std::unique_ptr<HeterogeneousHGCalHEFCellPositionsConditions> up =
0091       std::make_unique<HeterogeneousHGCalHEFCellPositionsConditions>(posmap_);
0092   return up;
0093 }
0094 
0095 #include "FWCore/Framework/interface/MakerMacros.h"
0096 #include "FWCore/Utilities/interface/typelookup.h"
0097 #include "FWCore/Framework/interface/eventsetuprecord_registration_macro.h"
0098 DEFINE_FWK_EVENTSETUP_MODULE(HeterogeneousHGCalHEFCellPositionsFiller);