CMS 3D CMS Logo

HGCalVFEProcessorSums.cc
Go to the documentation of this file.
2 #include <limits>
3 
6 
8 
11  "HGCalVFEProcessorSums");
12 
13 
16  vfeLinearizationImpl_(conf),
17  vfeSummationImpl_(conf),
18  calibration_( conf.getParameterSet("calib_parameters") )
19 {
20 }
21 
22 void
24  l1t::HGCalTriggerCellBxCollection& triggerCellColl,
25  const edm::EventSetup& es)
26 {
28 
29  std::vector<HGCalDataFrame> dataframes;
30  std::vector<std::pair<DetId, uint32_t >> linearized_dataframes;
31  std::map<HGCalDetId, uint32_t> payload;
32 
33  // convert ee and fh hit collections into the same object
34  for(const auto& digiData : digiColl)
35  {
36  if(DetId(digiData.id()).det()==DetId::Hcal && HcalDetId(digiData.id()).subdetId()!=HcalEndcap) continue;
37  uint32_t module = geometry_->getModuleFromCell(digiData.id());
38  if(geometry_->disconnectedModule(module)) continue;
39  dataframes.emplace_back(digiData.id());
40  for(int i=0; i<digiData.size(); i++)
41  {
42  dataframes.back().setSample(i, digiData.sample(i));
43  }
44  }
45 
46  vfeLinearizationImpl_.linearize(dataframes, linearized_dataframes);
47  vfeSummationImpl_.triggerCellSums(*geometry_, linearized_dataframes, payload);
48 
49  // Transform map to trigger cell vector vector<HGCalTriggerCell>
50  for(const auto& id_value : payload)
51  {
52  if (id_value.second>0){
53  l1t::HGCalTriggerCell triggerCell(reco::LeafCandidate::LorentzVector(), id_value.second, 0, 0, 0, id_value.first.rawId());
54  GlobalPoint point = geometry_->getTriggerCellPosition(id_value.first.rawId());
55 
56  // 'value' is hardware, so p4 is meaningless, except for eta and phi
57  math::PtEtaPhiMLorentzVector p4((double)id_value.second/cosh(point.eta()), point.eta(), point.phi(), 0.);
58  triggerCell.setP4(p4);
59  triggerCell.setPosition(point);
60 
61  // calibration part ---------------------------
62  if( triggerCell.hwPt() > 0 )
63  {
64  l1t::HGCalTriggerCell calibratedtriggercell( triggerCell );
65  calibration_.calibrateInGeV( calibratedtriggercell);
66  triggerCellColl.push_back(0, calibratedtriggercell);
67  }
68  }
69  }
70 }
HGCalVFESummationImpl vfeSummationImpl_
void run(const HGCalDigiCollection &digiColl, l1t::HGCalTriggerCellBxCollection &triggerCellColl, const edm::EventSetup &es) override
const HGCalTriggerGeometryBase * geometry_
ParameterSet const & getParameterSet(ParameterSetID const &id)
HGCalVFEProcessorSums(const edm::ParameterSet &conf)
void eventSetup(const edm::EventSetup &es)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
double p4[4]
Definition: TauolaWrapper.h:92
HGCalVFELinearizationImpl vfeLinearizationImpl_
void linearize(const std::vector< HGCDataFrame< DetId, HGCSample >> &, std::vector< std::pair< DetId, uint32_t > > &)
virtual bool disconnectedModule(const unsigned module_id) const =0
virtual unsigned getModuleFromCell(const unsigned cell_det_id) const =0
Definition: DetId.h:18
HGCalTriggerCellCalibration calibration_
#define DEFINE_EDM_PLUGIN(factory, type, name)
void triggerCellSums(const HGCalTriggerGeometryBase &, const std::vector< std::pair< DetId, uint32_t > > &, std::map< HGCalDetId, uint32_t > &payload)
void calibrateInGeV(l1t::HGCalTriggerCell &)
Definition: vlib.h:208
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
virtual GlobalPoint getTriggerCellPosition(const unsigned trigger_cell_det_id) const =0
void push_back(int bx, T object)