CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalTTPDigiProducer.cc
Go to the documentation of this file.
2 
10 #include <cstdio>
11 
12 // DO NOT MODIFY: Mapping between iphi (array index) and TTP input (value) for HF
13 const int HcalTTPDigiProducer::inputs_[] = { 30,66,4,44,4,44,0,68,
14  0,68,16,48,16,48,6,46,
15  6,46,2,70,2,70,18,50,
16  18,50,12,40,12,40,8,52,
17  8,52,20,36,20,36,14,42,
18  14,42,10,54,10,54,22,38,
19  22,38,24,56,24,56,32,60,
20  32,60,28,64,28,64,26,58,
21  26,58,34,62,34,62,30,66 } ;
22 
24 {
25  tok_hf_ = consumes<HFDigiCollection>(ps.getParameter<edm::InputTag>("HFDigiCollection"));
26  maskedChannels_ = ps.getParameter< std::vector<unsigned int> >("maskedChannels") ;
27  bit_[0] = ps.getParameter<std::string>("defTT8") ;
28  bit_[1] = ps.getParameter<std::string>("defTT9") ;
29  bit_[2] = ps.getParameter<std::string>("defTT10") ;
30  bit_[3] = ps.getParameter<std::string>("defTTLocal") ;
31 
32  for (int i=0; i<4; i++) {
33  nHits_[i] = -1 ; nHFp_[i] = -1 ; nHFm_[i] = -1 ;
34  pReq_[i] = ' ' ; mReq_[i] = ' ' ; pmLogic_[i] = ' ' ;
35  calc_[i] = sscanf(bit_[i].c_str(),"hits>=%d:hfp%c=%d%chfm%c=%d",
36  &(nHits_[i]),&(pReq_[i]),&(nHFp_[i]),
37  &(pmLogic_[i]),&(mReq_[i]),&(nHFm_[i])) ;
38  if ( calc_[i] == 1 ) {
39  if ( nHits_[i] < 0 )
40  throw cms::Exception("HcalTTPDigiProducer")
41  << "Unable to read logic for technical trigger" ;
42  } else if ( calc_[i] == 6 ) {
43  if ( nHits_[i] < 0 || nHFp_[i] < 0 || nHFm_[i] < 0 )
44  throw cms::Exception("HcalTTPDigiProducer")
45  << "Unable to read logic for technical trigger" ;
46  if ( (pReq_[i] != '>' && pReq_[i] != '<') ||
47  (mReq_[i] != '>' && mReq_[i] != '<') ||
48  (pmLogic_[i] != ':' && pmLogic_[i] != '|') )
49  throw cms::Exception("HcalTTPDigiProducer")
50  << "Technical Trigger logic must obey the following format:\n"
51  "\"hits>=[A1]:hfp[B1]=[A2][C]hfm[B2]=[A3]\",\n"
52  "or \"hits>=[A1]\",\n"
53  "with A# >= 0, B# = (</>) and C = (:/|)" ;
54  } else {
55  throw cms::Exception("HcalTTPDigiProducer")
56  << "Unable to read logic for technical trigger" ;
57  }
58  }
59 
60  id_ = ps.getUntrackedParameter<int>("id",-1) ;
61  samples_ = ps.getParameter<int>("samples") ;
62  presamples_ = ps.getParameter<int>("presamples") ;
63  iEtaMin_ = ps.getParameter<int>("iEtaMin") ;
64  iEtaMax_ = ps.getParameter<int>("iEtaMax") ;
65  threshold_ = ps.getParameter<unsigned int>("threshold") ;
66  fwAlgo_ = ps.getParameter<int>("fwAlgorithm") ;
67 
68  SoI_ = ps.getParameter<int>("HFSoI") ;
69 
70  if ( samples_ > 8 ) {
71  samples_ = 8 ;
72  edm::LogWarning("HcalTTPDigiProducer") << "Samples forced to maximum value of 8" ;
73  }
74  if ( presamples_ - SoI_ > 0 ) { // Too many presamples
75  presamples_ = SoI_ ;
76  edm::LogWarning("HcalTTPDigiProducer") << "Presamples reset to HF SoI value" ;
77  }
78 
79  produces<HcalTTPDigiCollection>();
80 }
81 
82 
84 }
85 
87 
88  for ( unsigned int i=0; i<maskedChannels_.size(); i++ )
89  if ( id.rawId() == maskedChannels_.at(i) ) return true ;
90  return false ;
91 }
92 
93 bool HcalTTPDigiProducer::decision(int nP, int nM, int bit) {
94 
95  bool pOK = false ; bool mOK = false ;
96  if ( (nP + nM) < nHits_[bit] ) return false ;
97  if ( calc_[bit] == 1 ) return ( (nP + nM) >= nHits_[bit] ) ;
98 
99  if ( pReq_[bit] == '>' ) pOK = ( nP >= nHFp_[bit] ) ;
100  else if ( pReq_[bit] == '<' ) pOK = ( nP <= nHFp_[bit] ) ;
101 
102  if ( mReq_[bit] == '>' ) mOK = ( nM >= nHFm_[bit] ) ;
103  else if ( mReq_[bit] == '<' ) mOK = ( nM <= nHFm_[bit] ) ;
104 
105  if ( pmLogic_[bit] == ':' ) return ( pOK && mOK ) ;
106  else if ( pmLogic_[bit] == '|' ) return ( pOK || mOK ) ;
107 
108  // Should not ever get here...need to create a warning message
109  edm::LogWarning("HcalTTPDigiProducer") << "Trigger logic exhausted. Returning false" ;
110  return false ;
111 }
112 
114 
115  // Step A: Get Inputs
116  edm::Handle<HFDigiCollection> hfDigiCollection ;
117  e.getByToken(tok_hf_,hfDigiCollection) ;
118  edm::ESHandle<HcalTPGCoder> inputCoder ;
119  eventSetup.get<HcalTPGRecord>().get(inputCoder) ;
120 
121  // Step B: Create empty output
122  std::auto_ptr<HcalTTPDigiCollection> ttpResult(new HcalTTPDigiCollection()) ;
123 
124  // Step C: Compute TTP inputs
125  uint16_t trigInputs[40] ;
126  int nP[8] ; int nM[8] ;
127  for (int i=0; i<8; i++) {
128  nP[i] = 0 ; nM[i] = 0 ;
129  for (int j=0; j<5; j++) trigInputs[j*8+i] = 0 ;
130  }
131  for (HFDigiCollection::const_iterator theDigi=hfDigiCollection->begin();
132  theDigi!=hfDigiCollection->end(); theDigi++) {
133  HcalDetId id = HcalDetId(theDigi->id()) ;
134  if ( isMasked(id) ) continue ;
135  if ( id.ietaAbs() < iEtaMin_ || id.ietaAbs() > iEtaMax_ ) continue ;
136 
137  IntegerCaloSamples samples(id,theDigi->size()) ;
138  inputCoder->adc2Linear(*theDigi,samples) ;
139 
140  for (int relSample=-presamples_; relSample<(samples_-presamples_); relSample++) {
141  if ( samples[SoI_+relSample] >= threshold_ ) {
142  int linSample = presamples_ + relSample ;
143  int offset = (-1+id.zside())/2 ;
144  int shift = inputs_[id.iphi()+offset] ;
145  int group = 0 ;
146  while ( shift >= 16 ) { shift -= 16 ; group++ ; }
147  if ( !(trigInputs[(linSample*8)+group]&(1<<shift)) )
148  ( id.ieta() > 0 ) ? ( nP[linSample]++) : ( nM[linSample]++ ) ;
149  trigInputs[(linSample*8)+group] |= (1<<shift) ;
150  }
151  }
152  }
153 
154  // Step D: Compute trigger decision and fill TTP digi
155  uint8_t trigOutput[8] ;
156  uint32_t algoDepBits[8] ;
158  for (int linSample=0; linSample<8; linSample++) {
159  trigOutput[linSample] = 0 ; algoDepBits[linSample] = 0 ;
160  if ( linSample<samples_) {
161  for (int j=0; j<4; j++)
162  trigOutput[linSample] |= (decision(nP[linSample],nM[linSample],j)<<j) ;
163  int nT = nP[linSample] + nM[linSample] ;
164 
165  // Algorithm Dependent bits for FW flavor = 1
166  // NOTE: this disagrees with the fw var. names that implies (LSB) T,M,P (MSB)
167  if ( fwAlgo_ == 1 ) algoDepBits[linSample] = (nT&0x7F) | ((nP[linSample]&0x3F)<<7) | ((nM[linSample]&0x3F)<<13) ;
168  ttpDigi.setSample((linSample-presamples_),&trigInputs[linSample*8],algoDepBits[linSample],trigOutput[linSample]) ;
169  }
170  }
171  ttpResult->push_back( ttpDigi ) ;
172 
173  // Step E: Put outputs into event
174  e.put(ttpResult);
175 }
176 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void produce(edm::Event &e, const edm::EventSetup &c)
std::vector< HFDataFrame >::const_iterator const_iterator
std::vector< unsigned int > maskedChannels_
bool isMasked(HcalDetId id)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
static const int inputs_[]
int j
Definition: DBlmapReader.cc:9
void setSample(int relativeSample, const uint16_t *triggerInputs, const uint32_t algodep, const uint8_t outputTrigger)
Definition: HcalTTPDigi.cc:24
bool decision(int nP, int nM, int bit)
unsigned int offset(bool)
tuple group
Definition: watchdog.py:82
const T & get() const
Definition: EventSetup.h:55
edm::SortedCollection< HcalTTPDigi > HcalTTPDigiCollection
edm::EDGetTokenT< HFDigiCollection > tok_hf_
static unsigned int const shift
volatile std::atomic< bool > shutdown_flag false
HcalTTPDigiProducer(const edm::ParameterSet &ps)