CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CastorTTRecord.cc
Go to the documentation of this file.
2 
13 
15 {
16  CastorDigiColl_ = consumes<CastorDigiCollection>(ps.getParameter<edm::InputTag>("CastorDigiCollection")) ;
17  CastorSignalTS_ = ps.getParameter< unsigned int >("CastorSignalTS") ;
18 
19  ttpBits_ = ps.getParameter< std::vector<unsigned int> >("ttpBits");
20  TrigNames_ = ps.getParameter< std::vector<std::string> >("TriggerBitNames");
21  TrigThresholds_ = ps.getParameter< std::vector<double> >("TriggerThresholds");
22 
23  reweighted_gain = 1.0;
24 
25  produces<L1GtTechnicalTriggerRecord>();
26 }
27 
28 
30 }
31 
33 
34  // std::cerr << "**** RUNNING THROUGH CastorTTRecord::produce" << std::endl;
35 
36  std::vector<L1GtTechnicalTrigger> vecTT(ttpBits_.size()) ;
37 
38  // Get Inputs
39  edm::Handle<CastorDigiCollection> CastorDigiColl ;
40  e.getByToken(CastorDigiColl_,CastorDigiColl) ;
41 
42  if ( !CastorDigiColl.failedToGet() ) {
43 
44  double cas_efC[16][14];
45  getEnergy_fC(cas_efC,CastorDigiColl,e,eventSetup);
46 
47 
48  std::vector<bool> decision(ttpBits_.size());
49 
50  getTriggerDecisions(decision,cas_efC);
51 
52  for(unsigned int i=0; i<ttpBits_.size(); i++) {
53  // if( decision.at(i) ) std::cerr << "**** Something Triggered" << std::endl;
54  // std::cout << "Run CastorTTRecord::produce. TriggerBit = " << ttpBits_.at(i) << "; TriggerName = " << TrigNames_.at(i) << "; Decision = " << decision[i] << std::endl;
55  vecTT.at(i) = L1GtTechnicalTrigger(TrigNames_.at(i), ttpBits_.at(i), 0, decision.at(i)) ;
56  }
57 
58  } else {
59  vecTT.clear() ;
60  }
61 
62  // Put output into event
63  std::auto_ptr<L1GtTechnicalTriggerRecord> output(new L1GtTechnicalTriggerRecord()) ;
64  output->setGtTechnicalTrigger(vecTT) ;
65  e.put(output) ;
66 }
67 
68 
70  edm::Event& e, const edm::EventSetup& eventSetup)
71 {
72  // std::cerr << "**** RUNNING THROUGH CastorTTRecord::getEnergy_fC" << std::endl;
73 
74  // Get Conditions
76  eventSetup.get<CastorDbRecord>().get(conditions) ;
77  const CastorQIEShape* shape = conditions->getCastorShape () ; // this one is generic
78 
79  for(int isec=0; isec<16; isec++) for(int imod=0; imod<14; imod++) energy[isec][imod] = 0;
80 
81  // Loop over digis
83  for (idigi=CastorDigiColl->begin(); idigi!=CastorDigiColl->end(); idigi++) {
84  const CastorDataFrame & digi = (*idigi) ;
85  HcalCastorDetId cell = digi.id() ;
86 
87  // Get Castor Coder
88  const CastorQIECoder* channelCoder = conditions->getCastorCoder(cell);
89  CastorCoderDb coder (*channelCoder, *shape);
90 
91  // Get Castor Calibration
92  const CastorCalibrations& calibrations=conditions->getCastorCalibrations(cell);
93 
94  // convert adc to fC
95  CaloSamples tool ;
96  coder.adc2fC(digi,tool) ;
97 
98  // pedestal substraction
99  int capid=digi[CastorSignalTS_].capid();
100  double fC = tool[CastorSignalTS_] - calibrations.pedestal(capid);
101 
102  // to correct threshold levels in fC for different gains
103  reweighted_gain = calibrations.gain(capid) / 0.015;
104 
105  energy[digi.id().sector()-1][digi.id().module()-1] = fC;
106  }
107 }
108 
109 void CastorTTRecord::getTriggerDecisions(std::vector<bool>& decision, double energy[16][14]) const
110 {
111  // std::cerr << "**** RUNNING THROUGH CastorTTRecord::getTriggerDecisions" << std::endl;
112 
113  // check if number of bits is at least four
114  if( decision.size() < 4 ) return;
115 
116  std::vector<bool> tdpo[8]; // TriggerDecisionsPerOctant
117  getTriggerDecisionsPerOctant(tdpo,energy);
118 
119 
120  // preset trigger decisions
121  decision.at(0) = true;
122  decision.at(1) = false;
123  decision.at(2) = false;
124  decision.at(3) = false;
125 
126  // canceld for low pt jet
127  // bool EM_decision = false;
128  // bool HAD_decision = false;
129  // loop over castor octants
130  for(int ioct=0; ioct<8; ioct++) {
131  int next_oct = (ioct+1)%8;
132  int prev_oct = (ioct+8-1)%8;
133 
134  // gap Trigger
135  if( !tdpo[ioct].at(0) ) decision.at(0) = false;
136  if( !tdpo[ioct].at(1) ) decision.at(0) = false;
137 
138  // jet Trigger
139  if( tdpo[ioct].at(2) ) decision.at(1) = true;
140 
141  // electron
142  // canceld for low pt jet
143  // if( tdpo[ioct].at(3) ) EM_decision = true;
144  // if( tdpo[ioct].at(4) ) HAD_decision = true;
145 
146  // iso muon
147  if( tdpo[ioct].at(5) ) {
148  // was one of the other sectors
149  // in the octant empty ?
150  if( tdpo[ioct].at(0) ) {
151  if( tdpo[prev_oct].at(1) &&
152  tdpo[next_oct].at(0) &&
153  tdpo[next_oct].at(1) )
154  decision.at(3) = true;
155  }
156  else if( tdpo[ioct].at(1) ) {
157  if( tdpo[prev_oct].at(0) &&
158  tdpo[prev_oct].at(1) &&
159  tdpo[next_oct].at(0) )
160  decision.at(3) = true;
161  }
162  // when not no iso muon
163  }
164 
165  // low pt jet Trigger
166  if( tdpo[ioct].at(6) ) decision.at(2) = true;
167  }
168 
169  // for EM Trigger whole castor not hadronic and somewhere EM
170  // canceld for low pt jet
171  // decision.at(2) = EM_decision && !HAD_decision;
172 }
173 
174 void CastorTTRecord::getTriggerDecisionsPerOctant(std::vector<bool> tdpo[8], double energy[16][14]) const
175 {
176  // std::cerr << "**** RUNNING THROUGH CastorTTRecord::getTriggerDecisionsPerOctant" << std::endl;
177 
178  // loop over octatants
179  for(int ioct=0; ioct<8; ioct++)
180  {
181  // six bits from HTR card
182  // 0. first sector empty
183  // 1. second sector empty
184  // 2. jet any sector
185  // 3. EM any sector
186  // 4. HAD any sector
187  // 5. muon any sector
188  // add instead of EM Trigger (not bit 6 in real)
189  // 6. low pt jet any sector
190  tdpo[ioct].resize(7);
191 
192  for(int ibit=0; ibit<7; ibit++)
193  tdpo[ioct].at(ibit) = false;
194 
195  // loop over castor sectors in octant
196  for(int ioctsec=0; ioctsec<2; ioctsec++)
197  {
198  // absolute sector number
199  int isec = 2*ioct+ioctsec;
200 
201  // init module sums for every sector
202  double fCsum_mod = 0;
203  double fCsum_em = 0, fCsum_ha = 0;
204  double fCsum_jet_had = 0;
205  double fCsum_col[3] = { 0, 0, 0 };
206 
207  // loop over modules
208  for(int imod=0; imod<14; imod++) {
209  // total sum
210  fCsum_mod += energy[isec][imod];
211 
212  // EM & HAD sum
213  if( imod < 2 ) fCsum_em += energy[isec][imod];
214  if( imod > 2 && imod < 12 ) fCsum_ha += energy[isec][imod];
215 
216  // sum over three sector parts
217  if( imod < 4 ) fCsum_col[0] += energy[isec][imod];
218  else if( imod < 8 ) fCsum_col[1] += energy[isec][imod];
219  else if( imod < 12 ) fCsum_col[2] += energy[isec][imod];
220 
221  // HAD sum for jet trigger v2
222  if( imod > 1 && imod < 5 ) fCsum_jet_had += energy[isec][imod];
223  }
224 
225  // gap Trigger
226  if( fCsum_mod < TrigThresholds_.at(0) ) {
227  if( ioctsec == 0 ) tdpo[ioct].at(0) = true;
228  else if( ioctsec == 1 ) tdpo[ioct].at(1) = true;
229  }
230 
231  // jet Trigger
232  // with gain correction
233  /* old version of jet trigger ( deprecated because of saturation )
234  if( fCsum_mod > TrigThresholds_.at(1) / reweighted_gain )
235  tdpo[ioct].at(2) = true;
236  */
237  if( fCsum_jet_had > TrigThresholds_.at(1) / reweighted_gain )
238  // additional high threshold near saturation for EM part
239  if( energy[isec][0] > 26000 / reweighted_gain && energy[isec][1] > 26000 / reweighted_gain )
240  tdpo[ioct].at(2) = true;
241 
242  // low pt jet Trigger
243  if( fCsum_mod > TrigThresholds_.at(5) / reweighted_gain )
244  tdpo[ioct].at(6) = true;
245 
246  // egamma Trigger
247  // with gain correction only in the EM threshold
248  if( fCsum_em > TrigThresholds_.at(2) / reweighted_gain )
249  tdpo[ioct].at(3) = true;
250  if( fCsum_ha > TrigThresholds_.at(3) )
251  tdpo[ioct].at(4) = true;
252 
253  // muon Trigger
254  int countColumns = 0;
255  for( int icol=0; icol<3; icol++ )
256  if( fCsum_col[icol] > TrigThresholds_.at(4) )
257  countColumns++;
258  if( countColumns >= 2 )
259  tdpo[ioct].at(5) = true;
260  }
261  }
262 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
int sector() const
get the sector (1-16)
std::vector< unsigned int > ttpBits_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
CastorTTRecord(const edm::ParameterSet &ps)
unsigned int CastorSignalTS_
std::vector< std::string > TrigNames_
std::vector< CastorDataFrame >::const_iterator const_iterator
int module() const
get the module (1-2 for EM, 1-12 for HAD)
std::vector< double > TrigThresholds_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
virtual void produce(edm::Event &e, const edm::EventSetup &c)
void getTriggerDecisionsPerOctant(std::vector< bool > tdps[16], double energy[16][14]) const
double pedestal(int fCapId) const
get pedestal for capid=0..3
double gain(int fCapId) const
get gain for capid=0..3
void getTriggerDecisions(std::vector< bool > &decision, double energy[16][14]) const
const T & get() const
Definition: EventSetup.h:56
virtual void adc2fC(const CastorDataFrame &df, CaloSamples &lf) const
edm::EDGetTokenT< CastorDigiCollection > CastorDigiColl_
const HcalCastorDetId & id() const
void getEnergy_fC(double energy[16][14], edm::Handle< CastorDigiCollection > &CastorDigiColl, edm::Event &e, const edm::EventSetup &c)
double reweighted_gain
virtual ~CastorTTRecord()