00001 #include "FWCore/Framework/interface/Event.h"
00002 #include "FWCore/Framework/interface/EventSetup.h"
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00006
00007
00008 #include "CondFormats/SiStripObjects/interface/SiStripLatency.h"
00009 #include "CondFormats/DataRecord/interface/SiStripCondDataRecords.h"
00010 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012
00013
00014 #include "DPGAnalysis/SiStripTools/interface/EventWithHistoryFilter.h"
00015
00016
00017 EventWithHistoryFilter::EventWithHistoryFilter():
00018 _historyProduct(),
00019 _partition(),
00020 _APVPhaseLabel(),
00021 _apvmodes(),
00022 _dbxrange(), _dbxrangelat(),
00023 _bxrange(), _bxrangelat(),
00024 _bxcyclerange(), _bxcyclerangelat(),
00025 _dbxcyclerange(), _dbxcyclerangelat(),
00026 _dbxtrpltrange(),
00027 _noAPVPhase(true)
00028 {
00029 printConfig();
00030 }
00031
00032 EventWithHistoryFilter::EventWithHistoryFilter(const edm::ParameterSet& iConfig):
00033 _historyProduct(iConfig.getUntrackedParameter<edm::InputTag>("historyProduct",edm::InputTag("consecutiveHEs"))),
00034 _partition(iConfig.getUntrackedParameter<std::string>("partitionName","Any")),
00035 _APVPhaseLabel(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel","APVPhases")),
00036 _apvmodes(iConfig.getUntrackedParameter<std::vector<int> >("apvModes",std::vector<int>())),
00037 _dbxrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxRange",std::vector<int>())),
00038 _dbxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware",std::vector<int>())),
00039 _bxrange(iConfig.getUntrackedParameter<std::vector<int> >("absBXRange",std::vector<int>())),
00040 _bxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware",std::vector<int>())),
00041 _bxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange",std::vector<int>())),
00042 _bxcyclerangelat(iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware",std::vector<int>())),
00043 _dbxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange",std::vector<int>())),
00044 _dbxcyclerangelat(iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware",std::vector<int>())),
00045 _dbxtrpltrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange",std::vector<int>()))
00046 {
00047 _noAPVPhase = isAPVPhaseNotNeeded();
00048 printConfig();
00049 }
00050
00051 void EventWithHistoryFilter::set(const edm::ParameterSet& iConfig) {
00052
00053
00054 _historyProduct = iConfig.getUntrackedParameter<edm::InputTag>("historyProduct",edm::InputTag("consecutiveHEs"));
00055 _partition = iConfig.getUntrackedParameter<std::string>("partitionName","Any");
00056 _APVPhaseLabel = iConfig.getUntrackedParameter<std::string>("APVPhaseLabel","APVPhases");
00057 _dbxrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxRange",std::vector<int>());
00058 _dbxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware",std::vector<int>());
00059 _bxrange = iConfig.getUntrackedParameter<std::vector<int> >("absBXRange",std::vector<int>());
00060 _bxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware",std::vector<int>());
00061 _bxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange",std::vector<int>());
00062 _bxcyclerangelat = iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware",std::vector<int>());
00063 _dbxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange",std::vector<int>());
00064 _dbxcyclerangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware",std::vector<int>());
00065 _dbxtrpltrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange",std::vector<int>());
00066
00067 _noAPVPhase = isAPVPhaseNotNeeded();
00068 printConfig();
00069
00070 }
00071
00072 const bool EventWithHistoryFilter::selected(const EventWithHistory& he, const edm::EventSetup& iSetup) const {
00073
00074 const std::vector<int> dummy;
00075 return is_selected(he,iSetup,dummy);
00076
00077 }
00078
00079 const bool EventWithHistoryFilter::selected(const EventWithHistory& he, const edm::Event& iEvent, const edm::EventSetup& iSetup) const {
00080
00081 const std::vector<int> apvphases = getAPVPhase(iEvent);
00082 return is_selected(he,iSetup,apvphases);
00083
00084 }
00085
00086 const bool EventWithHistoryFilter::selected(const edm::Event& event, const edm::EventSetup& iSetup) const {
00087
00088 const std::vector<int> apvphases = getAPVPhase(event);
00089
00090 edm::Handle<EventWithHistory> hEvent;
00091 event.getByLabel(_historyProduct,hEvent);
00092
00093 return is_selected(*hEvent,iSetup,apvphases);
00094
00095 }
00096
00097
00098 const bool EventWithHistoryFilter::is_selected(const EventWithHistory& he, const edm::EventSetup& iSetup, const std::vector<int> apvphases) const {
00099
00100
00101 const int latency = getAPVLatency(iSetup);
00102
00103
00104 bool selected = true;
00105
00106 if(!isAPVModeNotNeeded()) {
00107 const int apvmode = getAPVMode(iSetup);
00108 bool modeok = false;
00109 for(std::vector<int>::const_iterator wantedmode =_apvmodes.begin();wantedmode!=_apvmodes.end();++wantedmode) {
00110 modeok = modeok || (apvmode == *wantedmode);
00111 }
00112 if(!modeok) return false;
00113 }
00114
00115 selected = selected && (isCutInactive(_dbxrange) || isInRange(he.deltaBX(),_dbxrange,he.depth()!=0));
00116
00117 selected = selected && (isCutInactive(_dbxrangelat) ||
00118 isInRange(he.deltaBX()-latency,_dbxrangelat,he.depth()!=0 && latency>=0));
00119
00120 selected = selected && (isCutInactive(_bxrange) || isInRange(he.absoluteBX()%70,_bxrange,1));
00121
00122 selected = selected && (isCutInactive(_bxrangelat) ||
00123 isInRange((he.absoluteBX()-latency)%70,_bxrangelat,latency>=0));
00124
00125
00126
00127
00128 bool phaseselected;
00129
00130 phaseselected = isCutInactive(_bxcyclerange);
00131 for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) {
00132 phaseselected = phaseselected || isInRange(he.absoluteBXinCycle(*phase)%70,_bxcyclerange,*phase>0);
00133 }
00134 selected = selected && phaseselected;
00135
00136 phaseselected = isCutInactive(_bxcyclerangelat);
00137 for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) {
00138 phaseselected = phaseselected || isInRange((he.absoluteBXinCycle(*phase)-latency)%70,_bxcyclerangelat,
00139 *phase>=0 && latency>=0);
00140 }
00141 selected = selected && phaseselected;
00142
00143 phaseselected = isCutInactive(_dbxcyclerange);
00144 for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) {
00145 phaseselected = phaseselected || isInRange(he.deltaBXinCycle(*phase),_dbxcyclerange,he.depth()!=0 && *phase>=0);
00146 }
00147 selected = selected && phaseselected;
00148
00149 phaseselected = isCutInactive(_dbxcyclerangelat);
00150 for(std::vector<int>::const_iterator phase=apvphases.begin();phase!=apvphases.end();++phase) {
00151 phaseselected = phaseselected || isInRange(he.deltaBXinCycle(*phase)-latency,_dbxcyclerangelat,
00152 he.depth()!=0 && *phase>=0 && latency>=0);
00153 }
00154 selected = selected && phaseselected;
00155
00156
00157
00158 selected = selected && (isCutInactive(_dbxtrpltrange) ||
00159 isInRange(he.deltaBX(1,2),_dbxtrpltrange,he.depth()>1));
00160
00161 return selected;
00162
00163 }
00164
00165 const int EventWithHistoryFilter::getAPVLatency(const edm::EventSetup& iSetup) const {
00166
00167 if(isAPVLatencyNotNeeded()) return -1;
00168
00169 edm::ESHandle<SiStripLatency> apvlat;
00170 iSetup.get<SiStripLatencyRcd>().get(apvlat);
00171 const int latency = apvlat->singleLatency()!=255 ? apvlat->singleLatency(): -1;
00172
00173
00174
00175
00176
00177
00178
00179 return latency;
00180
00181 }
00182
00183 const int EventWithHistoryFilter::getAPVMode(const edm::EventSetup& iSetup) const {
00184
00185 if(isAPVModeNotNeeded()) return -1;
00186
00187 edm::ESHandle<SiStripLatency> apvlat;
00188 iSetup.get<SiStripLatencyRcd>().get(apvlat);
00189 int mode = -1;
00190 if(apvlat->singleReadOutMode()==1) mode = 47;
00191 if(apvlat->singleReadOutMode()==0) mode = 37;
00192
00193
00194
00195
00196
00197
00198
00199 return mode;
00200
00201 }
00202
00203 const std::vector<int> EventWithHistoryFilter::getAPVPhase(const edm::Event& iEvent) const {
00204
00205 if(_noAPVPhase) {
00206 const std::vector<int> dummy;
00207 return dummy;
00208 }
00209
00210 edm::Handle<APVCyclePhaseCollection> apvPhases;
00211 iEvent.getByLabel(_APVPhaseLabel,apvPhases);
00212
00213 const std::vector<int> phases = apvPhases->getPhases(_partition.c_str());
00214
00215
00216
00217
00218
00219
00220
00221
00222 return phases;
00223 }
00224
00225 const bool EventWithHistoryFilter::isAPVLatencyNotNeeded() const {
00226
00227 return
00228 isCutInactive(_bxrangelat) &&
00229 isCutInactive(_dbxrangelat) &&
00230 isCutInactive(_bxcyclerangelat) &&
00231 isCutInactive(_dbxcyclerangelat);
00232
00233 }
00234
00235 const bool EventWithHistoryFilter::isAPVPhaseNotNeeded() const {
00236
00237 return
00238 isCutInactive(_bxcyclerange) &&
00239 isCutInactive(_dbxcyclerange) &&
00240 isCutInactive(_bxcyclerangelat) &&
00241 isCutInactive(_dbxcyclerangelat);
00242
00243 }
00244
00245 const bool EventWithHistoryFilter::isAPVModeNotNeeded() const {
00246
00247 return (_apvmodes.size()==0) ;
00248
00249 }
00250
00251 const bool EventWithHistoryFilter::isCutInactive(const std::vector<int>& range) const {
00252
00253 return ((range.size()==0 ||
00254 (range.size()==1 && range[0]<0) ||
00255 (range.size()==2 && range[0]<0 && range[1]<0)));
00256
00257 }
00258
00259 const bool EventWithHistoryFilter::isInRange(const long long bx, const std::vector<int>& range, const bool extra) const {
00260
00261 bool cut1 = range.size()<1 || range[0]<0 || (extra && bx >= range[0]);
00262 bool cut2 = range.size()<2 || range[1]<0 || (extra && bx <= range[1]);
00263
00264 if(range.size()>=2 && range[0]>=0 && range[1]>=0 && (range[0] > range[1])) {
00265 return cut1 || cut2;
00266 }
00267 else {
00268 return cut1 && cut2;
00269 }
00270
00271 }
00272
00273 void EventWithHistoryFilter::printConfig() const {
00274
00275 std::string msgcategory = "EventWithHistoryFilterConfiguration";
00276
00277 if(!(
00278 isCutInactive(_bxrange) &&
00279 isCutInactive(_bxrangelat) &&
00280 isCutInactive(_bxcyclerange) &&
00281 isCutInactive(_bxcyclerangelat) &&
00282 isCutInactive(_dbxrange) &&
00283 isCutInactive(_dbxrangelat) &&
00284 isCutInactive(_dbxcyclerange) &&
00285 isCutInactive(_dbxcyclerangelat) &&
00286 isCutInactive(_dbxtrpltrange)
00287 )) {
00288
00289 edm::LogInfo(msgcategory.c_str()) << "historyProduct: "
00290 << _historyProduct.label() << " "
00291 << _historyProduct.instance() << " "
00292 << _historyProduct.process() << " "
00293 << " APVCyclePhase: "
00294 << _APVPhaseLabel;
00295
00296 edm::LogVerbatim(msgcategory.c_str()) << "-----------------------";
00297 edm::LogVerbatim(msgcategory.c_str()) << "List of active cuts:";
00298 if(!isCutInactive(_bxrange)) {
00299 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00300 if(_bxrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBX lower limit " << _bxrange[0];
00301 if(_bxrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBX upper limit " << _bxrange[1];
00302 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00303 }
00304 if(!isCutInactive(_bxrangelat)) {
00305 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00306 if(_bxrangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXLtcyAware lower limit "
00307 << _bxrangelat[0];
00308 if(_bxrangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXLtcyAware upper limit "
00309 << _bxrangelat[1];
00310 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00311 }
00312 if(!isCutInactive(_bxcyclerange)) {
00313 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00314 edm::LogVerbatim(msgcategory.c_str()) <<"absoluteBXinCycle partition: " << _partition;
00315 if(_bxcyclerange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycle lower limit "
00316 << _bxcyclerange[0];
00317 if(_bxcyclerange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycle upper limit "
00318 << _bxcyclerange[1];
00319 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00320 }
00321 if(!isCutInactive(_bxcyclerangelat)) {
00322 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00323 edm::LogVerbatim(msgcategory.c_str()) <<"absoluteBXinCycleLtcyAware partition: " << _partition;
00324 if(_bxcyclerangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycleLtcyAware lower limit "
00325 << _bxcyclerangelat[0];
00326 if(_bxcyclerangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "absoluteBXinCycleLtcyAware upper limit "
00327 << _bxcyclerangelat[1];
00328 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00329 }
00330 if(!isCutInactive(_dbxrange)) {
00331 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00332 if(_dbxrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBX lower limit " << _dbxrange[0];
00333 if(_dbxrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBX upper limit " << _dbxrange[1];
00334 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00335 }
00336 if(!isCutInactive(_dbxrangelat)) {
00337 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00338 if(_dbxrangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXLtcyAware lower limit "
00339 << _dbxrangelat[0];
00340 if(_dbxrangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXLtcyAware upper limit "
00341 << _dbxrangelat[1];
00342 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00343 }
00344 if(!isCutInactive(_dbxcyclerange)) {
00345 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00346 edm::LogVerbatim(msgcategory.c_str()) <<"deltaBXinCycle partition: " << _partition;
00347 if(_dbxcyclerange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycle lower limit "
00348 << _dbxcyclerange[0];
00349 if(_dbxcyclerange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycle upper limit "
00350 << _dbxcyclerange[1];
00351 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00352 }
00353 if(!isCutInactive(_dbxcyclerangelat)) {
00354 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00355 edm::LogVerbatim(msgcategory.c_str()) <<"deltaBXinCycleLtcyAware partition: " << _partition;
00356 if(_dbxcyclerangelat.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycleLtcyAware lower limit "
00357 << _dbxcyclerangelat[0];
00358 if(_dbxcyclerangelat.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "deltaBXinCycleLtcyAware upper limit "
00359 << _dbxcyclerangelat[1];
00360 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00361 }
00362 if(!isCutInactive(_dbxtrpltrange)) {
00363 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00364 if(_dbxtrpltrange.size()>=1) edm::LogVerbatim(msgcategory.c_str()) << "TripletIsolation lower limit "
00365 << _dbxtrpltrange[0];
00366 if(_dbxtrpltrange.size()>=2) edm::LogVerbatim(msgcategory.c_str()) << "TripletIsolation upper limit "
00367 << _dbxtrpltrange[1];
00368 edm::LogVerbatim(msgcategory.c_str()) << "......................";
00369 }
00370 edm::LogVerbatim(msgcategory.c_str()) << "-----------------------";
00371 }
00372 }