Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include <memory>
00022
00023
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/EDProducer.h"
00026
00027 #include "FWCore/Framework/interface/Event.h"
00028 #include "FWCore/Framework/interface/Run.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 #include "FWCore/Utilities/interface/Exception.h"
00035
00036 #include <map>
00037 #include <string>
00038
00039 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00040
00041
00042
00043
00044
00045 class ConfigurableAPVCyclePhaseProducer : public edm::EDProducer {
00046 public:
00047 explicit ConfigurableAPVCyclePhaseProducer(const edm::ParameterSet&);
00048 ~ConfigurableAPVCyclePhaseProducer();
00049
00050 private:
00051 virtual void beginJob() ;
00052 virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
00053 virtual void produce(edm::Event&, const edm::EventSetup&) override;
00054 virtual void endJob() ;
00055
00056
00057
00058 const std::vector<std::string> _defpartnames;
00059 const std::vector<int> _defphases;
00060
00061 std::map<int,std::vector<std::string> > _runpartnames;
00062 std::map<int,std::vector<int> > _runphases;
00063
00064 APVCyclePhaseCollection _currapvphases;
00065 };
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 ConfigurableAPVCyclePhaseProducer::ConfigurableAPVCyclePhaseProducer(const edm::ParameterSet& iConfig):
00080 _defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")),
00081 _defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")),
00082 _currapvphases()
00083 {
00084
00085 produces<APVCyclePhaseCollection,edm::InEvent>();
00086
00087
00088
00089 if(_defphases.size() < _defpartnames.size() ) {
00090
00091 throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent default phases/partitions vector sizes: "
00092 << _defphases.size() << " "
00093 << _defpartnames.size();
00094 }
00095
00096 std::vector<edm::ParameterSet> vps(iConfig.getParameter<std::vector<edm::ParameterSet> >("runPhases"));
00097
00098 for(std::vector<edm::ParameterSet>::const_iterator ps = vps.begin();ps!=vps.end();ps++) {
00099 _runphases[ps->getParameter<int>("runNumber")] = ps->getUntrackedParameter<std::vector<int> >("phases",_defphases);
00100 _runpartnames[ps->getParameter<int>("runNumber")] = ps->getUntrackedParameter<std::vector<std::string> >("partitions",_defpartnames);
00101
00102 if(_runphases[ps->getParameter<int>("runNumber")].size() < _runpartnames[ps->getParameter<int>("runNumber")].size() ) {
00103
00104 throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent run " << ps->getParameter<int>("runNumber")
00105 << " phases/partitions vector sizes: "
00106 << _runphases[ps->getParameter<int>("runNumber")].size() << " "
00107 << _runpartnames[ps->getParameter<int>("runNumber")].size();
00108 }
00109
00110 }
00111
00112 }
00113
00114
00115 ConfigurableAPVCyclePhaseProducer::~ConfigurableAPVCyclePhaseProducer()
00116 {
00117
00118
00119
00120
00121 }
00122
00123
00124
00125
00126
00127
00128
00129 void
00130 ConfigurableAPVCyclePhaseProducer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
00131 {
00132
00133 using namespace edm;
00134
00135 _currapvphases.get().clear();
00136
00137
00138
00139 const std::map<int,std::vector<std::string> >& _crunpartnames = _runpartnames;
00140 const std::map<int,std::vector<int> >& _crunphases = _runphases;
00141
00142 std::map<int,std::vector<int> >::const_iterator trphases = _crunphases.find(iRun.run());
00143 std::map<int,std::vector<std::string> >::const_iterator trpartnames = _crunpartnames.find(iRun.run());
00144
00145 std::vector<int> phases = _defphases;
00146 std::vector<std::string> partnames = _defpartnames;
00147
00148 if(trphases != _crunphases.end()) {
00149 phases = trphases->second;
00150 }
00151 if(trpartnames != _crunpartnames.end()) {
00152 partnames = trpartnames->second;
00153 }
00154
00155 if(phases.size() < partnames.size() ) {
00156
00157 throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: "
00158 << phases.size() << " "
00159 << partnames.size();
00160 }
00161
00162 for(unsigned int ipart=0;ipart<partnames.size();++ipart) {
00163 if(phases[ipart]>=0) {
00164 _currapvphases.get()[partnames[ipart]] = phases[ipart];
00165 }
00166 }
00167
00168
00169 for(std::map<std::string,int>::const_iterator it=_currapvphases.get().begin(); it!=_currapvphases.get().end();it++) {
00170
00171 edm::LogInfo("APVCyclePhaseProducerDebug") << "partition " << it->first << " phase " << it->second;
00172
00173 }
00174
00175 }
00176
00177 void
00178 ConfigurableAPVCyclePhaseProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00179 {
00180
00181 using namespace edm;
00182
00183 std::auto_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection(_currapvphases) );
00184
00185 iEvent.put(apvphases);
00186
00187 }
00188
00189
00190 void
00191 ConfigurableAPVCyclePhaseProducer::beginJob()
00192 {
00193 }
00194
00195
00196 void
00197 ConfigurableAPVCyclePhaseProducer::endJob() {
00198 }
00199
00200
00201 DEFINE_FWK_MODULE(ConfigurableAPVCyclePhaseProducer);