CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TkMeasurementDetSet.cc
Go to the documentation of this file.
1 #include "TkMeasurementDetSet.h"
4 
6 
7 
8 void StMeasurementDetSet::init(std::vector<TkStripMeasurementDet> & stripDets) {
9  // assume vector is full and ordered!
10  int size = stripDets.size();
11 
12  empty_.resize(size,true);
13  activeThisEvent_.resize(size,true);
14  activeThisPeriod_.resize(size,true);
15  id_.resize(size);
16  subId_.resize(size);
17  totalStrips_.resize(size);
18 
19  bad128Strip_.resize(size*6);
20  hasAny128StripBad_.resize(size);
21  badStripBlocks_.resize(size);
22 
23  if (isRegional()) {
24  clusterI_.resize(2*size);
25  } else {
26  detSet_.resize(size);
27  }
28 
29  for (int i=0; i!=size; ++i) {
30  auto & mdet = stripDets[i];
31  mdet.setIndex(i);
32  //intialize the detId !
33  id_[i] = mdet.specificGeomDet().geographicalId().rawId();
35  //initalize the total number of strips
36  totalStrips_[i] = mdet.specificGeomDet().specificTopology().nstrips();
37  }
38 }
39 
40 
41 void StMeasurementDetSet::set128StripStatus(int i, bool good, int idx) {
42  int offset = nbad128*i;
43  if (idx == -1) {
44  std::fill(bad128Strip_.begin()+offset, bad128Strip_.begin()+offset+6, !good);
45  hasAny128StripBad_[i] = !good;
46  } else {
47  bad128Strip_[offset+idx] = !good;
48  if (good == false) {
49  hasAny128StripBad_[i] = false;
50  } else { // this should not happen, as usually you turn on all fibers
51  // and then turn off the bad ones, and not vice-versa,
52  // so I don't care if it's not optimized
53  hasAny128StripBad_[i] = true;
54  for (int j = 0; i < (totalStrips_[j] >> 7); j++) {
55  if (bad128Strip_[j+offset] == false) hasAny128StripBad_[i] = false; break;
56  }
57  }
58  }
59  }
60 
61 
62 void StMeasurementDetSet::initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags, const edm::ParameterSet& cutPset) {
63  if (qualityFlags & BadStrips) {
68  }
69  setMaskBad128StripBlocks((qualityFlags & MaskBad128StripBlocks) != 0);
70 
71 
72  if ((quality != 0) && (qualityFlags != 0)) {
73  edm::LogInfo("MeasurementTracker") << "qualityFlags = " << qualityFlags;
74  unsigned int on = 0, tot = 0;
75  unsigned int foff = 0, ftot = 0, aoff = 0, atot = 0;
76  for (int i=0; i!=nDet(); i++) {
77  uint32_t detid = id(i);
78  if (qualityFlags & BadModules) {
79  bool isOn = quality->IsModuleUsable(detid);
80  setActive(i,isOn);
81  tot++; on += (unsigned int) isOn;
82  if (qualityDebugFlags & BadModules) {
83  edm::LogInfo("MeasurementTracker")<< "MeasurementTrackerImpl::initializeStripStatus : detid " << detid << " is " << (isOn ? "on" : "off");
84  }
85  } else {
86  setActive(i,true);
87  }
88  // first turn all APVs and fibers ON
89  set128StripStatus(i,true);
90  if (qualityFlags & BadAPVFibers) {
91  short badApvs = quality->getBadApvs(detid);
92  short badFibers = quality->getBadFibers(detid);
93  for (int j = 0; j < 6; j++) {
94  atot++;
95  if (badApvs & (1 << j)) {
96  set128StripStatus(i,false, j);
97  aoff++;
98  }
99  }
100  for (int j = 0; j < 3; j++) {
101  ftot++;
102  if (badFibers & (1 << j)) {
103  set128StripStatus(i,false, 2*j);
104  set128StripStatus(i,false, 2*j+1);
105  foff++;
106  }
107  }
108  }
109  auto & badStrips = getBadStripBlocks(i);
110  badStrips.clear();
111  if (qualityFlags & BadStrips) {
112  SiStripBadStrip::Range range = quality->getRange(detid);
113  for (SiStripBadStrip::ContainerIterator bit = range.first; bit != range.second; ++bit) {
114  badStrips.push_back(quality->decode(*bit));
115  }
116  }
117  }
118  if (qualityDebugFlags & BadModules) {
119  edm::LogInfo("MeasurementTracker StripModuleStatus") <<
120  " Total modules: " << tot << ", active " << on <<", inactive " << (tot - on);
121  }
122  if (qualityDebugFlags & BadAPVFibers) {
123  edm::LogInfo("MeasurementTracker StripAPVStatus") <<
124  " Total APVs: " << atot << ", active " << (atot-aoff) <<", inactive " << (aoff);
125  edm::LogInfo("MeasurementTracker StripFiberStatus") <<
126  " Total Fibers: " << ftot << ", active " << (ftot-foff) <<", inactive " << (foff);
127  }
128  } else {
129  for (int i=0; i!=nDet(); i++) {
130  setActive(i,true); // module ON
131  set128StripStatus(i,true); // all APVs and fibers ON
132  }
133  }
134 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
short getBadFibers(const uint32_t &detid) const
string fill
Definition: lumiContext.py:319
std::vector< unsigned int >::const_iterator ContainerIterator
std::vector< bool > activeThisPeriod_
void initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags, const edm::ParameterSet &cutPset)
std::vector< BadStripBlock > & getBadStripBlocks(int i)
std::vector< unsigned int > id_
void set128StripStatus(int i, bool good, int idx=-1)
std::vector< bool > bad128Strip_
unsigned int id(int i) const
BadStripCuts badStripCuts_[4]
std::vector< bool > hasAny128StripBad_
std::vector< std::vector< BadStripBlock > > badStripBlocks_
short getBadApvs(const uint32_t &detid) const
bool IsModuleUsable(const uint32_t &detid) const
std::vector< unsigned int > clusterI_
std::vector< unsigned char > subId_
int j
Definition: DBlmapReader.cc:9
unsigned int offset(bool)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
void setMaskBad128StripBlocks(bool maskThem)
std::vector< bool > empty_
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const Range getRange(const uint32_t &detID) const
std::vector< int > totalStrips_
void setActive(int i, bool active)
Turn on/off the module for reconstruction, for the full run or lumi (using info from DB...
static const int nbad128
std::vector< StripDetset > detSet_
std::pair< ContainerIterator, ContainerIterator > Range
std::vector< bool > activeThisEvent_
void init(std::vector< TkStripMeasurementDet > &stripDets)
tuple size
Write out results.
data decode(const unsigned int &value) const