CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ESDaqInfoTask.cc
Go to the documentation of this file.
1 #include <iostream>
2 
7 
10 
12 
16 
19 
21 
23 
25 
26 using namespace cms;
27 using namespace edm;
28 using namespace std;
29 
31 
32  dqmStore_ = Service<DQMStore>().operator->();
33 
34  prefixME_ = ps.getUntrackedParameter<string>("prefixME", "");
35 
36  enableCleanup_ = ps.getUntrackedParameter<bool>("enableCleanup", false);
37 
38  mergeRuns_ = ps.getUntrackedParameter<bool>("mergeRuns", false);
39 
40  ESFedRangeMin_ = ps.getUntrackedParameter<int>("ESFedRangeMin", 520);
41  ESFedRangeMax_ = ps.getUntrackedParameter<int>("ESFedRangeMax", 575);
42 
43  meESDaqFraction_ = 0;
44  meESDaqActiveMap_ = 0;
45  meESDaqError_ = 0;
46 
47  for (int i = 0; i < 56; i++) {
48  meESDaqActive_[i] = 0;
49  }
50 
51  if (ps.exists("esMapping")){
52  edm::ParameterSet esMap=ps.getParameter<edm::ParameterSet>("esMapping");
53  es_mapping_ = new ESElectronicsMapper(esMap);
54  }else{
55  edm::LogError("ESDaqInfoTask")<<"preshower mapping pointer not initialized. Temporary.";
56  es_mapping_=0;
57  }
58 
59 
60 
61 }
62 
64  delete es_mapping_;
65 }
66 
68 
69  char histo[200];
70 
71  if ( dqmStore_ ) {
72 
73  dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
74 
75  sprintf(histo, "DAQSummary");
76  meESDaqFraction_ = dqmStore_->bookFloat(histo);
77  meESDaqFraction_->Fill(0.0);
78 
79  sprintf(histo, "DAQSummaryMap");
80  meESDaqActiveMap_ = dqmStore_->book2D(histo,histo, 80, 0.5, 80.5, 80, 0.5, 80.5);
81  meESDaqActiveMap_->setAxisTitle("Si X", 1);
82  meESDaqActiveMap_->setAxisTitle("Si Y", 2);
83 
84  dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/DAQContents");
85 
86  for (int i = 0; i < 56; i++) {
87  sprintf(histo, "EcalPreshower_%d", ESFedRangeMin_+i);
88  meESDaqActive_[i] = dqmStore_->bookFloat(histo);
89  meESDaqActive_[i]->Fill(0.0);
90 
91  ESOnFed_[i] = false;
92  for ( int x = 0; x < 80; x++ ) {
93  for ( int y = 0; y < 80; y++ ) {
94  int iz = (x<40)? 1:2;
95  int ip = (y>=40)? 1:2;
96  int ix = (x<40)? x:x-40;
97  int iy = (y<40)? y:y-40;
98  int ifed = (*es_mapping_).getFED( iz, ip, ix, iy);
99  if(ifed == ESFedRangeMin_+i){
100  ESOnFed_[i] = true;
101  break;
102  }
103  }
104  if(ESOnFed_[i] == true) break;
105  }
106  }
107 
108  dqmStore_->setCurrentFolder(prefixME_ + "/ESIntegrityTask");
109  sprintf(histo, "DAQError");
110  meESDaqError_ = dqmStore_->book1D(histo, histo, 56, ESFedRangeMin_-0.5, ESFedRangeMax_+0.5);
111  meESDaqError_->setAxisTitle("FedID", 1);
112 
113  }
114 
115 }
116 
118 
119  if ( enableCleanup_ ) this->cleanup();
120 
121 }
122 
124 
125  this->reset();
126 
127  for (int x = 0; x < 80; ++x) {
128  for (int y = 0; y < 80; ++y) {
129  int iz = (x<40) ? 1:2;
130  int ip = (y>=40) ? 1:2;
131  int ix = (x<40) ? x:x-40;
132  int iy = (y<40) ? y:y-40;
133  int ifed = (*es_mapping_).getFED( iz, ip, ix+1, iy+1);
134  if( ifed > 0 ) meESDaqActiveMap_->setBinContent( x+1, y+1, 0.0 );
135  else meESDaqActiveMap_->setBinContent( x+1, y+1, -1.0 );
136  }
137  }
138 
139  for (int i = 0; i < 56; i++) {
140  if ( meESDaqError_ ) meESDaqError_->setBinContent(i, 0.0);
141  }
142 
144 
145  if( iSetup.find( recordKey ) ) {
146 
147  edm::ESHandle<RunInfo> sumFED;
148  iSetup.get<RunInfoRcd>().get(sumFED);
149 
150  std::vector<int> FedsInIds= sumFED->m_fed_in;
151 
152  float ESFedCount = 0.;
153 
154  for( unsigned int fedItr=0; fedItr<FedsInIds.size(); ++fedItr ) {
155 
156  int fedID=FedsInIds[fedItr];
157 
158  if ( fedID >= ESFedRangeMin_ && fedID <= ESFedRangeMax_ ) {
159 
160  if( ESOnFed_[fedID - ESFedRangeMin_] ) ESFedCount++;
161 
162  if ( meESDaqActive_[fedID-ESFedRangeMin_] ) meESDaqActive_[fedID-ESFedRangeMin_]->Fill(1.0);
163 
164  if( meESDaqActiveMap_ ) {
165 
166  for( int x = 0; x < 80; x++ ) {
167  for( int y = 0; y < 80; y++ ) {
168  int iz = (x<40)? 1:2;
169  int ip = (y>=40)? 1:2;
170  int ix = (x<40) ? x:x-40;
171  int iy = (x<40) ? y:y-40;
172  int ifed = es_mapping_->getFED(iz, ip, ix, iy);
173  if( fedID==ifed ) meESDaqActiveMap_->setBinContent( x+1, y+1, 1.0 );
174  }
175  }
176 
177  }
178 
179  if( meESDaqFraction_ ) meESDaqFraction_->Fill( ESFedCount/40. );
180 
181  if( meESDaqError_ ){
182  for( int i = 0; i < 56; i++){
183  if( ESOnFed_[fedID-ESFedRangeMin_] ) meESDaqError_->setBinContent(i+1, 1.0);
184  else meESDaqError_->setBinContent(i+1, 2.0);
185  }
186  }
187 
188  }
189 
190  }
191 
192  } else {
193 
194  LogWarning("ESDaqInfoTask") << "Cannot find any RunInfoRcd" << endl;
195 
196  }
197 
198 }
199 
201 
202 }
203 
205 
206  if ( meESDaqFraction_ ) meESDaqFraction_->Reset();
207 
208  for (int i = 0; i < 56; i++) {
209  if ( meESDaqActive_[i] ) meESDaqActive_[i]->Reset();
210  }
211 
212  if ( meESDaqActiveMap_ ) meESDaqActiveMap_->Reset();
213 
214  if ( meESDaqError_ ) meESDaqError_->Reset();
215 
216 }
217 
218 
220 
221  if ( dqmStore_ ) {
222 
223  dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
224 
225  if ( meESDaqFraction_ ) dqmStore_->removeElement( meESDaqFraction_->getName() );
226 
227  if ( meESDaqActiveMap_ ) dqmStore_->removeElement( meESDaqActiveMap_->getName() );
228 
229  if ( meESDaqError_ ) dqmStore_->removeElement( meESDaqError_->getName() );
230 
231  dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/DAQContents");
232 
233  for (int i = 0; i < 56; i++) {
234  if ( meESDaqActive_[i] ) dqmStore_->removeElement( meESDaqActive_[i]->getName() );
235  }
236 
237  }
238 
239 }
240 
241 void ESDaqInfoTask::analyze(const Event& e, const EventSetup& c){
242 
243 }
244 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void beginJob(void)
BeginJob.
void cleanup(void)
Cleanup.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
static void cleanup(const Factory::MakerMap::value_type &v)
Definition: Factory.cc:12
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:90
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
void reset(void)
Reset.
std::string getName(Reflex::Type &cc)
Definition: ClassFiller.cc:18
const T & get() const
Definition: EventSetup.h:55
void endJob(void)
EndJob.
x
Definition: VDTMath.h:216
void beginLuminosityBlock(const edm::LuminosityBlock &lumiBlock, const edm::EventSetup &iSetup)
BeginLuminosityBlock.
void reset(double vett[256])
Definition: TPedValues.cc:11
ESDaqInfoTask(const edm::ParameterSet &ps)
Constructor.
void endLuminosityBlock(const edm::LuminosityBlock &lumiBlock, const edm::EventSetup &iSetup)
EndLuminosityBlock.
static HCTypeTag findType(char const *iTypeName)
find a type based on the types name, if not found will return default HCTypeTag
Definition: HCTypeTag.cc:129
virtual ~ESDaqInfoTask()
Destructor.