CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DetStatus.cc
Go to the documentation of this file.
5 #include <iostream>
6 
7 using namespace std;
8 
9 //
10 // -- Constructor
11 //
13  verbose_ = pset.getUntrackedParameter<bool>("DebugOn", false);
14  AndOr_ = pset.getParameter<bool>("AndOr");
15  applyfilter_ = pset.getParameter<bool>("ApplyFilter");
16  DetNames_ = pset.getParameter<std::vector<std::string> >("DetectorType");
17 
18  // build a map
19  DetMap_ = 0;
20  for (unsigned int detreq = 0; detreq < DetNames_.size(); detreq++) {
21  for (unsigned int detlist = 0; detlist < DcsStatus::nPartitions; detlist++) {
22  if (DetNames_[detreq] == DcsStatus::partitionName[detlist]) {
23  DetMap_ |= (1 << DcsStatus::partitionList[detlist]);
24  if (verbose_)
25  std::cout << "DCSStatus filter: asked partition " << DcsStatus::partitionName[detlist] << " bit "
26  << DcsStatus::partitionList[detlist] << std::endl;
27  }
28  }
29  }
30  edm::InputTag scalersTag("scalersRawToDigi");
31  scalersToken_ = consumes<DcsStatusCollection>(scalersTag);
32 }
33 
34 //
35 // -- Destructor
36 //
38 
40  bool accepted = false;
41 
43  evt.getByToken(scalersToken_, dcsStatus);
44 
45  if (dcsStatus.isValid()) {
46  if (dcsStatus->empty()) {
47  // Only complain in case it's real data. Accept in any case.
48  if (evt.eventAuxiliary().isRealData())
49  edm::LogError("DetStatus") << "Error! dcsStatus has size 0, accept in any case";
50  accepted = true;
51  } else {
52  unsigned int curr_dcs = (*dcsStatus)[0].ready();
53  if (verbose_)
54  std::cout << "curr_dcs = " << curr_dcs << std::endl;
55  if (AndOr_)
56  accepted = ((DetMap_ & curr_dcs) == DetMap_);
57  else
58  accepted = ((DetMap_ & curr_dcs) != 0);
59 
60  if (verbose_)
61 
62  {
63  std::cout << "DCSStatus filter: requested map: " << DetMap_ << " dcs in event: " << curr_dcs
64  << " filter: " << accepted << std::endl;
65  std::cout << "Partitions ON: ";
66  for (unsigned int detlist = 0; detlist < DcsStatus::nPartitions; detlist++) {
67  if ((*dcsStatus)[0].ready(DcsStatus::partitionList[detlist])) {
68  std::cout << " " << DcsStatus::partitionName[detlist];
69  }
70  }
71  std::cout << std::endl;
72  }
73  }
74  } else {
75  edm::LogError("DetStatus") << "Error! can't get the product: scalersRawToDigi, accept in any case";
76  accepted = true;
77  }
78 
79  if (!applyfilter_)
80  accepted = true;
81  return accepted;
82 }
83 
DetStatus(const edm::ParameterSet &)
Definition: DetStatus.cc:12
~DetStatus() override
Definition: DetStatus.cc:37
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool isRealData() const
EventAuxiliary const & eventAuxiliary() const override
Definition: Event.h:93
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool filter(edm::Event &, edm::EventSetup const &) override
Definition: DetStatus.cc:39
static const char *const partitionName[]
Definition: DcsStatus.h:31
bool isValid() const
Definition: HandleBase.h:70
bool accepted(std::vector< std::string_view > const &, std::string_view)
static const int partitionList[]
Definition: DcsStatus.h:30