CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
SiPixelDetectorStatus Class Reference

#include <SiPixelDetectorStatus.h>

Public Member Functions

void addModule (int detid, int nrocs)
 
void addModule (int detid, SiPixelModuleStatus a)
 
std::map< int, SiPixelModuleStatus >::iterator begin ()
 
unsigned long int digiOccDET ()
 
void dumpToFile (std::string filename)
 
std::map< int, SiPixelModuleStatus >::iterator end ()
 
void fillDIGI (int detid, int roc)
 
void fillFEDerror25 (int detid, PixelFEDChannel ch)
 
bool findModule (int detid)
 
std::map< int, SiPixelModuleStatusgetDetectorStatus ()
 
std::map< int, std::vector< int > > getFEDerror25Rocs ()
 
std::pair< int, int > getLSRange ()
 
SiPixelModuleStatusgetModule (int detid)
 
unsigned long int getNevents ()
 
std::pair< int, int > getRunRange ()
 
std::map< int, SiPixelModuleStatus >::iterator next ()
 
int nmodules ()
 
double perRocDigiOcc ()
 
double perRocDigiOccVar ()
 
void readFromFile (std::string filename)
 
void resetDetectorStatus ()
 
void setLSRange (int ls0, int ls1)
 
void setNevents (unsigned long int N)
 
void setRunRange (int run0, int run1)
 
 SiPixelDetectorStatus ()
 
void updateDetectorStatus (SiPixelDetectorStatus newData)
 
 ~SiPixelDetectorStatus ()
 

Private Attributes

unsigned long int fDetHits
 
int fLS0
 
int fLS1
 
std::map< int, SiPixelModuleStatusfModules
 
unsigned long int fNevents
 
int fRun0
 
int fRun1
 

Detailed Description

Definition at line 11 of file SiPixelDetectorStatus.h.

Constructor & Destructor Documentation

SiPixelDetectorStatus::SiPixelDetectorStatus ( )

Definition at line 14 of file SiPixelDetectorStatus.cc.

References fDetHits, and fNevents.

14  : fLS0(99999999), fLS1(0), fRun0(99999999), fRun1(0) {
15 
16  fDetHits = 0;
17  fNevents = 0;
18 
19 }
unsigned long int fDetHits
unsigned long int fNevents
SiPixelDetectorStatus::~SiPixelDetectorStatus ( )

Definition at line 22 of file SiPixelDetectorStatus.cc.

22  {
23 
24 }

Member Function Documentation

void SiPixelDetectorStatus::addModule ( int  detid,
int  nrocs 
)

Definition at line 114 of file SiPixelDetectorStatus.cc.

References a, and fModules.

Referenced by readFromFile().

114  {
115 
116  SiPixelModuleStatus a(detid, nrocs);
117  fModules.insert(std::make_pair(detid, a));
118 
119 }
std::map< int, SiPixelModuleStatus > fModules
double a
Definition: hdecay.h:121
void SiPixelDetectorStatus::addModule ( int  detid,
SiPixelModuleStatus  a 
)

Definition at line 122 of file SiPixelDetectorStatus.cc.

References fModules.

122  {
123 
124  fModules.insert(std::make_pair(detid, a));
125 
126 }
std::map< int, SiPixelModuleStatus > fModules
std::map< int, SiPixelModuleStatus >::iterator SiPixelDetectorStatus::begin ( )

Definition at line 172 of file SiPixelDetectorStatus.cc.

References fModules.

Referenced by digiOccDET(), dumpToFile(), getFEDerror25Rocs(), perRocDigiOcc(), perRocDigiOccVar(), and updateDetectorStatus().

172  {
173 
174  return fModules.begin();
175 
176 }
std::map< int, SiPixelModuleStatus > fModules
unsigned long int SiPixelDetectorStatus::digiOccDET ( )
inline

Definition at line 38 of file SiPixelDetectorStatus.h.

References begin(), end(), fDetHits, findModule(), getModule(), next(), and nmodules().

Referenced by SiPixelStatusManager::readLumi(), and updateDetectorStatus().

38 { return fDetHits; }
unsigned long int fDetHits
void SiPixelDetectorStatus::dumpToFile ( std::string  filename)

Definition at line 92 of file SiPixelDetectorStatus.cc.

References begin(), end(), fDetHits, fLS0, fLS1, fRun0, and fRun1.

92  {
93 
94  std::ofstream OD(filename.c_str());
95  OD << "# SiPixelDetectorStatus START" << std::endl;
96  OD << "# SiPixelDetectorStatus for LS " << fLS0 << " .. " << fLS1 << std::endl;
97  OD << "# SiPixelDetectorStatus for run " << fRun0 << " .. " << fRun1 << std::endl;
98  OD << "# SiPixelDetectorStatus total hits = " << fDetHits << std::endl;
99 
100  for (std::map<int, SiPixelModuleStatus>::iterator it = SiPixelDetectorStatus::begin(); it != SiPixelDetectorStatus::end(); ++it) {
101  for (int iroc = 0; iroc < it->second.nrocs(); ++iroc) {
102  for (int idc = 0; idc < 26; ++idc) {
103  OD << Form("%10d %2d %3d", it->first, iroc, int(it->second.getRoc(iroc)->digiOccROC())) << std::endl;
104  }
105  }
106  }
107  OD << "# SiPixelDetectorStatus END" << std::endl;
108  OD.close();
109 
110 }
unsigned long int fDetHits
std::map< int, SiPixelModuleStatus >::iterator end()
std::map< int, SiPixelModuleStatus >::iterator begin()
std::map< int, SiPixelModuleStatus >::iterator SiPixelDetectorStatus::end ( )

Definition at line 184 of file SiPixelDetectorStatus.cc.

References fModules.

Referenced by Types.LuminosityBlockRange::cppID(), Types.EventRange::cppID(), digiOccDET(), dumpToFile(), getFEDerror25Rocs(), perRocDigiOcc(), perRocDigiOccVar(), and updateDetectorStatus().

184  {
185 
186  return fModules.end();
187 
188 }
std::map< int, SiPixelModuleStatus > fModules
void SiPixelDetectorStatus::fillDIGI ( int  detid,
int  roc 
)

Definition at line 130 of file SiPixelDetectorStatus.cc.

References fDetHits, and fModules.

130  {
131 
132  ++fDetHits;
133  fModules[detid].fillDIGI(roc);
134 
135 }
std::map< int, SiPixelModuleStatus > fModules
unsigned long int fDetHits
void SiPixelDetectorStatus::fillFEDerror25 ( int  detid,
PixelFEDChannel  ch 
)

Definition at line 138 of file SiPixelDetectorStatus.cc.

References fModules.

138  {
139 
140  if (fModules.find(detid) != fModules.end()){
141  fModules[detid].fillFEDerror25(ch);
142  }
143 
144 }
std::map< int, SiPixelModuleStatus > fModules
bool SiPixelDetectorStatus::findModule ( int  detid)

Definition at line 207 of file SiPixelDetectorStatus.cc.

References fModules.

Referenced by digiOccDET().

207  {
208 
209  if (fModules.find(detid) == fModules.end())
210  return false;
211  else
212  return true;
213 
214 }
std::map< int, SiPixelModuleStatus > fModules
std::map<int, SiPixelModuleStatus> SiPixelDetectorStatus::getDetectorStatus ( )
inline

Definition at line 70 of file SiPixelDetectorStatus.h.

References fModules.

Referenced by SiPixelStatusHarvester::endRunProduce().

70 { return fModules; }
std::map< int, SiPixelModuleStatus > fModules
std::map< int, std::vector< int > > SiPixelDetectorStatus::getFEDerror25Rocs ( )

Definition at line 147 of file SiPixelDetectorStatus.cc.

References begin(), end(), SiPixelModuleStatus::getRoc(), SiPixelRocStatus::isFEDerror25(), list(), and SiPixelModuleStatus::nrocs().

Referenced by SiPixelStatusManager::createFEDerror25().

147  {
148 
149  std::map<int, std::vector<int>> badRocLists_;
150 
151  for(std::map<int, SiPixelModuleStatus>::iterator itMod = SiPixelDetectorStatus::begin(); itMod != SiPixelDetectorStatus::end(); ++itMod)
152  {
153  int detid = itMod->first;
154  // FEDerror25 effected ROCs in a given module
155  std::vector<int> list;
156  SiPixelModuleStatus modStatus = itMod->second;
157  for (int iroc = 0; iroc < modStatus.nrocs(); ++iroc) {
158 
159  SiPixelRocStatus* roc = modStatus.getRoc(iroc);
160  if(roc->isFEDerror25()){
161  list.push_back(iroc);
162  }
163  badRocLists_[detid]=list;
164  }
165 
166  }
167 
168  return badRocLists_;
169 }
std::map< int, SiPixelModuleStatus >::iterator end()
std::map< int, SiPixelModuleStatus >::iterator begin()
SiPixelRocStatus * getRoc(int i)
get a ROC
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
std::pair<int,int> SiPixelDetectorStatus::getLSRange ( )
inline

Definition at line 56 of file SiPixelDetectorStatus.h.

References fLS0, and fLS1.

Referenced by SiPixelStatusManager::rankByLumi().

SiPixelModuleStatus * SiPixelDetectorStatus::getModule ( int  detid)

Definition at line 198 of file SiPixelDetectorStatus.cc.

References fModules.

Referenced by digiOccDET(), readFromFile(), and updateDetectorStatus().

198  {
199 
200  if (fModules.find(detid) == fModules.end()) {
201  return nullptr;
202  }
203  return &(fModules[detid]);
204 
205 }
std::map< int, SiPixelModuleStatus > fModules
unsigned long int SiPixelDetectorStatus::getNevents ( )
inline

Definition at line 60 of file SiPixelDetectorStatus.h.

References fNevents.

Referenced by updateDetectorStatus().

60 { return fNevents; }
unsigned long int fNevents
std::pair<int,int> SiPixelDetectorStatus::getRunRange ( )
inline

Definition at line 54 of file SiPixelDetectorStatus.h.

References fRun0, and fRun1.

std::map<int, SiPixelModuleStatus>::iterator SiPixelDetectorStatus::next ( )

Referenced by digiOccDET().

int SiPixelDetectorStatus::nmodules ( )

Definition at line 191 of file SiPixelDetectorStatus.cc.

References fModules.

Referenced by digiOccDET().

191  {
192 
193  return static_cast<int>(fModules.size());
194 
195 }
std::map< int, SiPixelModuleStatus > fModules
double SiPixelDetectorStatus::perRocDigiOcc ( )

Definition at line 217 of file SiPixelDetectorStatus.cc.

References begin(), and end().

Referenced by SiPixelStatusManager::createBadComponents(), SiPixelStatusHarvester::endRunProduce(), and perRocDigiOccVar().

217  {
218 
219  unsigned long int ave(0);
220  int nrocs(0);
221  for (std::map<int, SiPixelModuleStatus>::iterator it = SiPixelDetectorStatus::begin(); it != SiPixelDetectorStatus::end(); ++it) {
222  unsigned long int inc = it->second.digiOccMOD();
223  ave += inc;
224  nrocs += it->second.nrocs();
225  }
226  return (1.0*ave)/nrocs;
227 
228 }
std::map< int, SiPixelModuleStatus >::iterator end()
std::map< int, SiPixelModuleStatus >::iterator begin()
double SiPixelDetectorStatus::perRocDigiOccVar ( )

Definition at line 230 of file SiPixelDetectorStatus.cc.

References begin(), end(), and perRocDigiOcc().

230  {
231 
232  double fDetAverage = SiPixelDetectorStatus::perRocDigiOcc();
233 
234  double sig = 0.0;
235  int nrocs(0);
236  for(std::map<int, SiPixelModuleStatus>::iterator it = SiPixelDetectorStatus::begin(); it != SiPixelDetectorStatus::end(); ++it) {
237  unsigned long int inc = it->second.digiOccMOD();
238  sig += (fDetAverage - inc) * (fDetAverage - inc);
239  nrocs += it->second.nrocs();
240  }
241 
242  double fDetSigma = sig/(nrocs - 1);
243  return TMath::Sqrt(fDetSigma);
244 }
std::map< int, SiPixelModuleStatus >::iterator end()
std::map< int, SiPixelModuleStatus >::iterator begin()
void SiPixelDetectorStatus::readFromFile ( std::string  filename)

Definition at line 28 of file SiPixelDetectorStatus.cc.

References addModule(), fDetHits, fLS0, fLS1, fRun0, fRun1, getModule(), hfClusterShapes_cfi::hits, SiPixelModuleStatus::setNrocs(), AlCaHLTBitMon_QueryRunRegistry::string, and SiPixelModuleStatus::updateModuleDIGI().

28  {
29 
30  std::ifstream INS;
31  std::string sline;
32  INS.open(filename.c_str());
33 
34  int oldDetId(-1);
35  int detid(0), roc(0), hits(0), nroc(0);
36  SiPixelModuleStatus *pMod(nullptr);
37  bool readOK(false);
38  while (std::getline(INS, sline)) {
39 
40  if (std::string::npos != sline.find("# SiPixelDetectorStatus START")) {
41  readOK = true;
42  continue;
43  }
44  if (!readOK) continue;
45 
46  if (std::string::npos != sline.find("# SiPixelDetectorStatus END")) {
47  pMod->setNrocs(nroc+1);
48  break;
49  }
50 
51  if (sline.find("# SiPixelDetectorStatus for LS") != std::string::npos) {
52  std::sscanf(sline.c_str(), "# SiPixelDetectorStatus for LS %d .. %d", &fLS0, &fLS1);
53  continue;
54  }
55  if (sline.find("# SiPixelDetectorStatus for run") != std::string::npos) {
56  std::sscanf(sline.c_str(), "# SiPixelDetectorStatus for run %d .. %d", &fRun0, &fRun1);
57  continue;
58  }
59  if (sline.find("# SiPixelDetectorStatus total hits = ") != std::string::npos) {
60  std::sscanf(sline.c_str(), "# SiPixelDetectorStatus total hits = %ld", &fDetHits);
61  continue;
62  }
63 
64  std::sscanf(sline.c_str(), "%d %d %d", &detid, &roc, &hits);
65  if (roc > nroc) nroc = roc;
66  if (detid != oldDetId) {
67  if (pMod) {
68  pMod->setNrocs(nroc+1);
69  }
70 
71  oldDetId = detid;
72  if (getModule(detid) == nullptr) {
73  addModule(detid,nroc+1);
74  }
75 
76  pMod = getModule(detid);
77  nroc = 0;
78  }
79  if (pMod) {
80  fDetHits += hits;
81  pMod->updateModuleDIGI(roc, hits);
82  }
83 
84  }
85 
86  INS.close();
87 
88 }
void addModule(int detid, int nrocs)
unsigned long int fDetHits
SiPixelModuleStatus * getModule(int detid)
void SiPixelDetectorStatus::resetDetectorStatus ( )
inline

Definition at line 62 of file SiPixelDetectorStatus.h.

References fDetHits, fLS0, fLS1, fModules, fNevents, fRun0, fRun1, and updateDetectorStatus().

62  { fModules.clear(); fDetHits=0; fNevents=0;
63  fRun0 = 99999999; fRun1 = 0; fLS0 = 99999999; fLS1 = 0;
64  }
std::map< int, SiPixelModuleStatus > fModules
unsigned long int fDetHits
unsigned long int fNevents
void SiPixelDetectorStatus::setLSRange ( int  ls0,
int  ls1 
)
inline

Definition at line 55 of file SiPixelDetectorStatus.h.

References fLS0, and fLS1.

Referenced by SiPixelStatusManager::createBadComponents().

void SiPixelDetectorStatus::setNevents ( unsigned long int  N)
inline

Definition at line 59 of file SiPixelDetectorStatus.h.

References fNevents, and N.

59 { fNevents = N; }
unsigned long int fNevents
#define N
Definition: blowfish.cc:9
void SiPixelDetectorStatus::setRunRange ( int  run0,
int  run1 
)
inline
void SiPixelDetectorStatus::updateDetectorStatus ( SiPixelDetectorStatus  newData)

Definition at line 247 of file SiPixelDetectorStatus.cc.

References begin(), digiOccDET(), end(), fDetHits, fModules, fNevents, getModule(), and getNevents().

Referenced by SiPixelStatusManager::createBadComponents(), and resetDetectorStatus().

247  {
248 
249  // loop over new data status
250  for(std::map<int, SiPixelModuleStatus>::iterator it = newData.begin(); it != newData.end(); ++it) {
251 
252  int detid = it->first;
253  if(fModules.find(detid) != fModules.end()){// if the detid is in the module lists
254  fModules[detid].updateModuleStatus( *(newData.getModule(detid)) );
255  }
256  else{
257  fModules.insert(std::make_pair(detid, *(newData.getModule(detid))));
258  }
259 
260  }
261 
262  fDetHits = fDetHits + newData.digiOccDET();
263  fNevents = fNevents + newData.getNevents();
264 
265 }
unsigned long int digiOccDET()
std::map< int, SiPixelModuleStatus > fModules
unsigned long int fDetHits
unsigned long int fNevents
std::map< int, SiPixelModuleStatus >::iterator end()
std::map< int, SiPixelModuleStatus >::iterator begin()
unsigned long int getNevents()
SiPixelModuleStatus * getModule(int detid)

Member Data Documentation

unsigned long int SiPixelDetectorStatus::fDetHits
private
int SiPixelDetectorStatus::fLS0
private
int SiPixelDetectorStatus::fLS1
private
std::map<int, SiPixelModuleStatus> SiPixelDetectorStatus::fModules
private
unsigned long int SiPixelDetectorStatus::fNevents
private
int SiPixelDetectorStatus::fRun0
private
int SiPixelDetectorStatus::fRun1
private