CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelClusterSource.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelMonitorCluster
4 // Class: SiPixelClusterSource
5 //
13 //
14 // Original Author: Vincenzo Chiochia & Andrew York
15 // Created:
16 //
17 //
18 // Updated by: Lukas Wehrli
19 // for pixel offline DQM
21 // Framework
24 // DQM Framework
27 // Geometry
32 // DataFormats
37 //
38 #include <string>
39 #include <stdlib.h>
40 
41 using namespace std;
42 using namespace edm;
43 
45  conf_(iConfig),
46  src_( conf_.getParameter<edm::InputTag>( "src" ) ),
47  saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
48  isPIB( conf_.getUntrackedParameter<bool>("isPIB",false) ),
49  slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
50  modOn( conf_.getUntrackedParameter<bool>("modOn",true) ),
51  twoDimOn( conf_.getUntrackedParameter<bool>("twoDimOn",true) ),
52  reducedSet( conf_.getUntrackedParameter<bool>("reducedSet",false) ),
53  ladOn( conf_.getUntrackedParameter<bool>("ladOn",false) ),
54  layOn( conf_.getUntrackedParameter<bool>("layOn",false) ),
55  phiOn( conf_.getUntrackedParameter<bool>("phiOn",false) ),
56  ringOn( conf_.getUntrackedParameter<bool>("ringOn",false) ),
57  bladeOn( conf_.getUntrackedParameter<bool>("bladeOn",false) ),
58  diskOn( conf_.getUntrackedParameter<bool>("diskOn",false) ),
59  smileyOn(conf_.getUntrackedParameter<bool>("smileyOn",false) ),
60  bigEventSize( conf_.getUntrackedParameter<int>("bigEventSize",100) )
61 {
63  LogInfo ("PixelDQM") << "SiPixelClusterSource::SiPixelClusterSource: Got DQM BackEnd interface"<<endl;
64 
65  //set Token(-s)
66  srcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(conf_.getParameter<edm::InputTag>("src"));
67 }
68 
69 
71 {
72  // do anything here that needs to be done at desctruction time
73  // (e.g. close files, deallocate resources etc.)
74  LogInfo ("PixelDQM") << "SiPixelClusterSource::~SiPixelClusterSource: Destructor"<<endl;
75 
76  std::map<uint32_t,SiPixelClusterModule*>::iterator struct_iter;
77  for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++){
78  delete struct_iter->second;
79  struct_iter->second = 0;
80  }
81 }
82 
83 
85  firstRun = true;
86 }
87 
89 
90  LogInfo ("PixelDQM") << " SiPixelClusterSource::beginJob - Initialisation ... " << std::endl;
91  LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/"
92  << layOn << "/" << phiOn << std::endl;
93  LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/"
94  << ringOn << std::endl;
95  LogInfo ("PixelDQM") << "2DIM IS " << twoDimOn << "\n";
96  LogInfo ("PixelDQM") << "Smiley (Cluster sizeY vs. Cluster eta) is " << smileyOn << "\n";
97 
98  if(firstRun){
99  eventNo = 0;
100  lumSec = 0;
101  nLumiSecs = 0;
102  nBigEvents = 0;
103  // Build map
104  buildStructure(iSetup);
105  // Book Monitoring Elements
106  bookMEs();
107  // Book occupancy maps in global coordinates for all clusters:
108  theDMBE->setCurrentFolder("Pixel/Clusters/OffTrack");
109  //bpix
110  meClPosLayer1 = theDMBE->book2D("position_siPixelClusters_Layer_1","Clusters Layer1;Global Z (cm);Global #phi",200,-30.,30.,128,-3.2,3.2);
111  meClPosLayer2 = theDMBE->book2D("position_siPixelClusters_Layer_2","Clusters Layer2;Global Z (cm);Global #phi",200,-30.,30.,128,-3.2,3.2);
112  meClPosLayer3 = theDMBE->book2D("position_siPixelClusters_Layer_3","Clusters Layer3;Global Z (cm);Global #phi",200,-30.,30.,128,-3.2,3.2);
113  //fpix
114  meClPosDisk1pz = theDMBE->book2D("position_siPixelClusters_pz_Disk_1","Clusters +Z Disk1;Global X (cm);Global Y (cm)",80,-20.,20.,80,-20.,20.);
115  meClPosDisk2pz = theDMBE->book2D("position_siPixelClusters_pz_Disk_2","Clusters +Z Disk2;Global X (cm);Global Y (cm)",80,-20.,20.,80,-20.,20.);
116  meClPosDisk1mz = theDMBE->book2D("position_siPixelClusters_mz_Disk_1","Clusters -Z Disk1;Global X (cm);Global Y (cm)",80,-20.,20.,80,-20.,20.);
117  meClPosDisk2mz = theDMBE->book2D("position_siPixelClusters_mz_Disk_2","Clusters -Z Disk2;Global X (cm);Global Y (cm)",80,-20.,20.,80,-20.,20.);
118 
119  firstRun = false;
120  }
121 }
122 
123 
125  if(saveFile){
126  LogInfo ("PixelDQM") << " SiPixelClusterSource::endJob - Saving Root File " << std::endl;
128  theDMBE->save( outputFile.c_str() );
129  }
130 }
131 
132 //------------------------------------------------------------------
133 // Method called for every event
134 //------------------------------------------------------------------
136 {
137  eventNo++;
138 
139  if(modOn){
140  MonitorElement* meReset = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_Layer_1");
141  MonitorElement* meReset1 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_Layer_2");
142  MonitorElement* meReset2 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_Layer_3");
143  MonitorElement* meReset3 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_mz_Disk_1");
144  MonitorElement* meReset4 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_mz_Disk_2");
145  MonitorElement* meReset5 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_pz_Disk_1");
146  MonitorElement* meReset6 = theDMBE->get("Pixel/Clusters/OffTrack/position_siPixelClusters_pz_Disk_2");
147  if(meReset && meReset->getEntries()>150000){
148  meReset->Reset();
149  meReset1->Reset();
150  meReset2->Reset();
151  meReset3->Reset();
152  meReset4->Reset();
153  meReset5->Reset();
154  meReset6->Reset();
155  }
156  }
157 
158  // get input data
160  iEvent.getByToken(srcToken_, input);
161 
163  iSetup.get<TrackerDigiGeometryRecord> ().get (pDD);
164  const TrackerGeometry* tracker = &(* pDD);
165 // const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(detId) );
166 
167  //float iOrbitSec = iEvent.orbitNumber()/11223.;
168  //int bx = iEvent.bunchCrossing();
169  //long long tbx = (long long)iEvent.orbitNumber() * 3564 + bx;
170  int lumiSection = (int)iEvent.luminosityBlock();
171  int nEventFpixClusters = 0;
172 
173 
174  std::map<uint32_t,SiPixelClusterModule*>::iterator struct_iter;
175  for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++) {
176 
177  int numberOfFpixClusters = (*struct_iter).second->fill(*input, tracker, modOn,
178  ladOn, layOn, phiOn,
179  bladeOn, diskOn, ringOn,
181  nEventFpixClusters = nEventFpixClusters + numberOfFpixClusters;
182 
183  }
184 
185 // if(lumiSection>lumSec){ lumSec = lumiSection; nLumiSecs++; }
186 // if(nEventFpixClusters>bigEventSize) nBigEvents++;
187 // if(nLumiSecs%5==0){
188 
189  if(nEventFpixClusters>bigEventSize){
190  MonitorElement* me = theDMBE->get("Pixel/bigFpixClusterEventRate");
191  if(me){
192  me->Fill(lumiSection,1./23.);
193  }
194  }
195  //std::cout<<"nEventFpixClusters: "<<nEventFpixClusters<<" , nLumiSecs: "<<nLumiSecs<<" , nBigEvents: "<<nBigEvents<<std::endl;
196 
197  // slow down...
198  if(slowDown) usleep(10000);
199 
200 }
201 
202 //------------------------------------------------------------------
203 // Build data structure
204 //------------------------------------------------------------------
206 
207  LogInfo ("PixelDQM") <<" SiPixelClusterSource::buildStructure" ;
209  iSetup.get<TrackerDigiGeometryRecord>().get( pDD );
210 
211  LogVerbatim ("PixelDQM") << " *** Geometry node for TrackerGeom is "<<&(*pDD)<<std::endl;
212  LogVerbatim ("PixelDQM") << " *** I have " << pDD->dets().size() <<" detectors"<<std::endl;
213  LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<" types"<<std::endl;
214 
215  for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
216 
217  if(dynamic_cast<PixelGeomDetUnit*>((*it))!=0){
218 
219  DetId detId = (*it)->geographicalId();
220  const GeomDetUnit * geoUnit = pDD->idToDetUnit( detId );
221  const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
222  int nrows = (pixDet->specificTopology()).nrows();
223  int ncols = (pixDet->specificTopology()).ncolumns();
224 
225 
226  if((detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) ||
227  (detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))){
228  uint32_t id = detId();
229  SiPixelClusterModule* theModule = new SiPixelClusterModule(id, ncols, nrows);
230  if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
231  if(isPIB) continue;
232  LogDebug ("PixelDQM") << " ---> Adding Barrel Module " << detId.rawId() << endl;
233  thePixelStructure.insert(pair<uint32_t,SiPixelClusterModule*> (id,theModule));
234  }else if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
235  LogDebug ("PixelDQM") << " ---> Adding Endcap Module " << detId.rawId() << endl;
237  int disk = PixelEndcapName(DetId(id)).diskName();
238  int blade = PixelEndcapName(DetId(id)).bladeName();
239  int panel = PixelEndcapName(DetId(id)).pannelName();
241  char sside[80]; sprintf(sside, "HalfCylinder_%i",side);
242  char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
243  char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
244  char spanel[80]; sprintf(spanel, "Panel_%i",panel);
245  char smodule[80];sprintf(smodule,"Module_%i",module);
246  std::string side_str = sside;
247  std::string disk_str = sdisk;
248  bool mask = side_str.find("HalfCylinder_1")!=string::npos||
249  side_str.find("HalfCylinder_2")!=string::npos||
250  side_str.find("HalfCylinder_4")!=string::npos||
251  disk_str.find("Disk_2")!=string::npos;
252  // clutch to take all of FPIX, but no BPIX:
253  mask = false;
254  if(isPIB && mask) continue;
255  thePixelStructure.insert(pair<uint32_t,SiPixelClusterModule*> (id,theModule));
256  }
257  }
258  }
259  }
260  LogInfo ("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
261 }
262 //------------------------------------------------------------------
263 // Book MEs
264 //------------------------------------------------------------------
266 
267  // Get DQM interface
269  theDMBE->setCurrentFolder("Pixel");
270  char title[256]; snprintf(title, 256, "Rate of events with >%i FPIX clusters;LumiSection;Rate of large FPIX events per LS [Hz]",bigEventSize);
271  bigFpixClusterEventRate = theDMBE->book1D("bigFpixClusterEventRate",title,5000,0.,5000.);
272 
273 
274  std::map<uint32_t,SiPixelClusterModule*>::iterator struct_iter;
275 
276  SiPixelFolderOrganizer theSiPixelFolder;
277 
278  for(struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++){
279 
281  if(modOn){
282  if(theSiPixelFolder.setModuleFolder((*struct_iter).first)){
283  (*struct_iter).second->book( conf_,0,twoDimOn,reducedSet);
284  } else {
285 
286  if(!isPIB) throw cms::Exception("LogicError")
287  << "[SiPixelClusterSource::bookMEs] Creation of DQM folder failed";
288  }
289  }
290  if(ladOn){
291  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,1)){
292  (*struct_iter).second->book( conf_,1,twoDimOn,reducedSet);
293  } else {
294  LogDebug ("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
295  }
296  }
297  if(layOn){
298  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,2)){
299  (*struct_iter).second->book( conf_,2,twoDimOn,reducedSet);
300  } else {
301  LogDebug ("PixelDQM") << "PROBLEM WITH LAYER-FOLDER\n";
302  }
303  }
304  if(phiOn){
305  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,3)){
306  (*struct_iter).second->book( conf_,3,twoDimOn,reducedSet);
307  } else {
308  LogDebug ("PixelDQM") << "PROBLEM WITH PHI-FOLDER\n";
309  }
310  }
311  if(bladeOn){
312  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,4)){
313  (*struct_iter).second->book( conf_,4,twoDimOn,reducedSet);
314  } else {
315  LogDebug ("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
316  }
317  }
318  if(diskOn){
319  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,5)){
320  (*struct_iter).second->book( conf_,5,twoDimOn,reducedSet);
321  } else {
322  LogDebug ("PixelDQM") << "PROBLEM WITH DISK-FOLDER\n";
323  }
324  }
325  if(ringOn){
326  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,6)){
327  (*struct_iter).second->book( conf_,6,twoDimOn,reducedSet);
328  } else {
329  LogDebug ("PixelDQM") << "PROBLEM WITH RING-FOLDER\n";
330  }
331  }
332  if(smileyOn){
333  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,7)){
334  (*struct_iter).second->book( conf_,7,twoDimOn,reducedSet);
335  } else {
336  LogDebug ("PixelDQM") << "PROBLEM WITH BARREL-FOLDER\n";
337  }
338  }
339 
340  }
341 
342 }
343 
344 //define this as a plug-in
#define LogDebug(id)
int plaquetteName() const
plaquetteId (in pannel)
virtual void buildStructure(edm::EventSetup const &)
T getParameter(std::string const &) const
MonitorElement * meClPosDisk2pz
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
double getEntries(void) const
get # of entries
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > srcToken_
SiPixelClusterSource(const edm::ParameterSet &conf)
static std::string const input
Definition: EdmProvDump.cc:44
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
MonitorElement * meClPosLayer3
int bladeName() const
blade id
int iEvent
Definition: GenABIO.cc:243
virtual void beginRun(const edm::Run &, edm::EventSetup const &)
MonitorElement * meClPosLayer1
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2296
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1623
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
edm::ParameterSet conf_
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:55
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
int pannelName() const
pannel id
MonitorElement * meClPosDisk2mz
int diskName() const
disk id
std::map< uint32_t, SiPixelClusterModule * > thePixelStructure
volatile std::atomic< bool > shutdown_flag false
MonitorElement * meClPosLayer2
HalfCylinder halfCylinder() const
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1000
void Reset(void)
reset ME (ie. contents, errors, etc)
Definition: vlib.h:208
MonitorElement * bigFpixClusterEventRate
MonitorElement * meClPosDisk1mz
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
MonitorElement * meClPosDisk1pz
Definition: Run.h:41
bool setModuleFolder(const uint32_t &rawdetid=0, int type=0)
Set folder name for a module or plaquette.