CMS 3D CMS Logo

SiPixelRecHitSource Class Reference

Description: header file for Pixel Monitor Rec Hits. More...

#include <DQM/SiPixelMonitorRecHit/interface/SiPixelRecHitSource.h>

Inheritance diagram for SiPixelRecHitSource:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (edm::EventSetup const &)
virtual void bookMEs ()
virtual void buildStructure (edm::EventSetup const &)
virtual void endJob ()
 SiPixelRecHitSource (const edm::ParameterSet &conf)
 ~SiPixelRecHitSource ()

Private Attributes

bool bladeOn
edm::ParameterSet conf_
bool diskOn
int eventNo
bool isPIB
bool ladOn
bool layOn
bool modOn
bool phiOn
std::map< uint32_t, intrechit_count
bool ringOn
bool saveFile
bool slowDown
edm::InputTag src_
DQMStoretheDMBE
std::map< uint32_t,
SiPixelRecHitModule * > 
thePixelStructure
bool twoDimOn


Detailed Description

Description: header file for Pixel Monitor Rec Hits.

Usage: see description

Definition at line 55 of file SiPixelRecHitSource.h.


Constructor & Destructor Documentation

SiPixelRecHitSource::SiPixelRecHitSource ( const edm::ParameterSet conf  )  [explicit]

Definition at line 50 of file SiPixelRecHitSource.cc.

References lat::endl(), and theDMBE.

00050                                                                        :
00051   conf_(iConfig),
00052   src_( conf_.getParameter<edm::InputTag>( "src" ) ),
00053   saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
00054   isPIB( conf_.getUntrackedParameter<bool>("isPIB",false) ),
00055   slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
00056   modOn( conf_.getUntrackedParameter<bool>("modOn",true) ),
00057   twoDimOn( conf_.getUntrackedParameter<bool>("twoDimOn",true) ),
00058   ladOn( conf_.getUntrackedParameter<bool>("ladOn",false) ), 
00059   layOn( conf_.getUntrackedParameter<bool>("layOn",false) ), 
00060   phiOn( conf_.getUntrackedParameter<bool>("phiOn",false) ), 
00061   ringOn( conf_.getUntrackedParameter<bool>("ringOn",false) ), 
00062   bladeOn( conf_.getUntrackedParameter<bool>("bladeOn",false) ), 
00063   diskOn( conf_.getUntrackedParameter<bool>("diskOn",false) )
00064 {
00065    theDMBE = edm::Service<DQMStore>().operator->();
00066    LogInfo ("PixelDQM") << "SiPixelRecHitSource::SiPixelRecHitSource: Got DQM BackEnd interface"<<endl;
00067 }

SiPixelRecHitSource::~SiPixelRecHitSource (  ) 

Definition at line 70 of file SiPixelRecHitSource.cc.

References lat::endl().

00071 {
00072    // do anything here that needs to be done at desctruction time
00073   // (e.g. close files, deallocate resources etc.)
00074   LogInfo ("PixelDQM") << "SiPixelRecHitSource::~SiPixelRecHitSource: Destructor"<<endl;
00075 }


Member Function Documentation

void SiPixelRecHitSource::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 111 of file SiPixelRecHitSource.cc.

References bladeOn, diskOn, eventNo, edm::Event::getByLabel(), ladOn, layOn, lp, modOn, phiOn, edm::Handle< T >::product(), rechit_count, ringOn, slowDown, src_, thePixelStructure, twoDimOn, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

00112 {
00113   eventNo++;
00114   //cout << eventNo << endl;
00115   // get input data
00116   edm::Handle<SiPixelRecHitCollection>  recHitColl;
00117   iEvent.getByLabel( src_, recHitColl );
00118 
00119   std::map<uint32_t,SiPixelRecHitModule*>::iterator struct_iter;
00120   for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++) {
00121     uint32_t TheID = (*struct_iter).first;
00122         
00123     SiPixelRecHitCollection::range pixelrechitRange = (recHitColl.product())->get(TheID);
00124     SiPixelRecHitCollection::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.first;
00125     
00126     SiPixelRecHitCollection::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.second;
00127     SiPixelRecHitCollection::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
00128 
00129       // if( pixelrechitRangeIteratorBegin == pixelrechitRangeIteratorEnd) {cout << "oops" << endl;}
00130       float rechit_x = 0;
00131       float rechit_y = 0;
00132       float rechit_count = 0;
00133   for ( ; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) 
00134         {
00135           
00136 
00137           rechit_count++;
00138           //cout << TheID << endl;
00139           SiPixelRecHit::ClusterRef const& clust = pixeliter->cluster();
00140           int sizeX = (*clust).sizeX();
00141           //cout << sizeX << endl;
00142           int sizeY = (*clust).sizeY();
00143           //cout << sizeY << endl;
00144           LocalPoint lp = pixeliter->localPosition();
00145           rechit_x = lp.x();
00146           rechit_y = lp.y();
00147           
00148           LocalError lerr = pixeliter->localPositionError();
00149           //  float lerr_x = sqrt(lerr.xx());
00150           //  float lerr_y = sqrt(lerr.yy());
00151           //cout << "hh" << endl;
00152           (*struct_iter).second->fill(rechit_x, rechit_y, sizeX, sizeY,modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn, twoDimOn);
00153           //cout << "ii" << endl;
00154         
00155         }
00156             (*struct_iter).second->nfill(rechit_count, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
00157     
00158   }
00159 
00160   // slow down...
00161   if(slowDown) usleep(10000);
00162   
00163 }

void SiPixelRecHitSource::beginJob ( edm::EventSetup const &  iSetup  )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 78 of file SiPixelRecHitSource.cc.

References bladeOn, bookMEs(), buildStructure(), diskOn, lat::endl(), eventNo, ladOn, layOn, modOn, phiOn, ringOn, and twoDimOn.

00078                                                              {
00079 
00080   LogInfo ("PixelDQM") << " SiPixelRecHitSource::beginJob - Initialisation ... " << std::endl;
00081   LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/" 
00082             << layOn << "/" << phiOn << std::endl;
00083   LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/" 
00084             << ringOn << std::endl;
00085   LogInfo ("PixelDQM") << "2DIM IS " << twoDimOn << "\n";
00086   eventNo = 0;
00087         
00088   // Build map
00089   buildStructure(iSetup);
00090         
00091   // Book Monitoring Elements
00092   bookMEs();
00093 
00094 }

void SiPixelRecHitSource::bookMEs (  )  [virtual]

Create folder tree and book histograms

Definition at line 233 of file SiPixelRecHitSource.cc.

References bladeOn, conf_, diskOn, Exception, isPIB, ladOn, layOn, LogDebug, modOn, phiOn, ringOn, SiPixelFolderOrganizer::setModuleFolder(), DQMStore::setVerbose(), theDMBE, thePixelStructure, and twoDimOn.

Referenced by beginJob().

00233                                  {
00234   
00235   std::map<uint32_t,SiPixelRecHitModule*>::iterator struct_iter;
00236   theDMBE->setVerbose(0);
00237     
00238   SiPixelFolderOrganizer theSiPixelFolder;
00239   
00240   for(struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++){
00241     
00243     if(modOn){
00244       if(theSiPixelFolder.setModuleFolder((*struct_iter).first)){
00245         (*struct_iter).second->book( conf_,0,twoDimOn);
00246       } else {
00247         if(!isPIB) throw cms::Exception("LogicError")
00248           << "[SiPixelDigiSource::bookMEs] Creation of DQM folder failed";
00249       }
00250     }
00251     if(ladOn){
00252       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,1)){
00253         (*struct_iter).second->book( conf_,1,twoDimOn);
00254         } else {
00255         LogDebug ("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
00256       }
00257     }
00258     if(layOn){
00259       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,2)){
00260         (*struct_iter).second->book( conf_,2,twoDimOn);
00261         } else {
00262         LogDebug ("PixelDQM") << "PROBLEM WITH LAYER-FOLDER\n";
00263       }
00264     }
00265     if(phiOn){
00266       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,3)){
00267         (*struct_iter).second->book( conf_,3,twoDimOn);
00268         } else {
00269         LogDebug ("PixelDQM") << "PROBLEM WITH PHI-FOLDER\n";
00270       }
00271     }
00272     if(bladeOn){
00273       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,4)){
00274         (*struct_iter).second->book( conf_,4,twoDimOn);
00275         } else {
00276         LogDebug ("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
00277       }
00278     }
00279     if(diskOn){
00280       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,5)){
00281         (*struct_iter).second->book( conf_,5,twoDimOn);
00282         } else {
00283         LogDebug ("PixelDQM") << "PROBLEM WITH DISK-FOLDER\n";
00284       }
00285     }
00286     if(ringOn){
00287       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,6)){
00288         (*struct_iter).second->book( conf_,6,twoDimOn);
00289         } else {
00290         LogDebug ("PixelDQM") << "PROBLEM WITH RING-FOLDER\n";
00291       }
00292     }
00293 
00294   }
00295 
00296 }

void SiPixelRecHitSource::buildStructure ( edm::EventSetup const &  iSetup  )  [virtual]

Definition at line 168 of file SiPixelRecHitSource.cc.

References DetId::DetId(), detId, muonGeometry::disk, lat::endl(), edm::EventSetup::get(), isPIB, it, LogDebug, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, PixelEndcapName::PixelEndcapName(), DetId::rawId(), DetId::subdetId(), and thePixelStructure.

Referenced by beginJob().

00168                                                                    {
00169 
00170   LogInfo ("PixelDQM") <<" SiPixelRecHitSource::buildStructure" ;
00171   edm::ESHandle<TrackerGeometry> pDD;
00172   iSetup.get<TrackerDigiGeometryRecord>().get( pDD );
00173 
00174   LogVerbatim ("PixelDQM") << " *** Geometry node for TrackerGeom is  "<<&(*pDD)<<std::endl;
00175   LogVerbatim ("PixelDQM") << " *** I have " << pDD->dets().size() <<" detectors"<<std::endl;
00176   LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<" types"<<std::endl;
00177   
00178   for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
00179     
00180     if(dynamic_cast<PixelGeomDetUnit*>((*it))!=0){
00181 
00182       DetId detId = (*it)->geographicalId();
00183       // const GeomDetUnit      * geoUnit = pDD->idToDetUnit( detId );
00184       //const PixelGeomDetUnit * pixDet  = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
00185 
00186           
00187           
00188               // SiPixelRecHitModule *theModule = new SiPixelRecHitModule(id, rechit_x, rechit_y, x_res, y_res, x_pull, y_pull);
00189         
00190         
00191               if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
00192                 if(isPIB) continue;
00193                 LogDebug ("PixelDQM") << " ---> Adding Barrel Module " <<  detId.rawId() << endl;
00194                 uint32_t id = detId();
00195                 SiPixelRecHitModule* theModule = new SiPixelRecHitModule(id);
00196                 thePixelStructure.insert(pair<uint32_t,SiPixelRecHitModule*> (id,theModule));
00197                 
00198               } else if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
00199                 LogDebug ("PixelDQM") << " ---> Adding Endcap Module " <<  detId.rawId() << endl;
00200                 uint32_t id = detId();
00201                 SiPixelRecHitModule* theModule = new SiPixelRecHitModule(id);
00202 
00203                 PixelEndcapName::HalfCylinder side = PixelEndcapName::PixelEndcapName(DetId::DetId(id)).halfCylinder();
00204                 int disk   = PixelEndcapName::PixelEndcapName(DetId::DetId(id)).diskName();
00205                 int blade  = PixelEndcapName::PixelEndcapName(DetId::DetId(id)).bladeName();
00206                 int panel  = PixelEndcapName::PixelEndcapName(DetId::DetId(id)).pannelName();
00207                 int module = PixelEndcapName::PixelEndcapName(DetId::DetId(id)).plaquetteName();
00208 
00209                 char sside[80];  sprintf(sside,  "HalfCylinder_%i",side);
00210                 char sdisk[80];  sprintf(sdisk,  "Disk_%i",disk);
00211                 char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
00212                 char spanel[80]; sprintf(spanel, "Panel_%i",panel);
00213                 char smodule[80];sprintf(smodule,"Module_%i",module);
00214                 std::string side_str = sside;
00215                 std::string disk_str = sdisk;
00216                 bool mask = side_str.find("HalfCylinder_1")!=string::npos||
00217                             side_str.find("HalfCylinder_2")!=string::npos||
00218                             side_str.find("HalfCylinder_4")!=string::npos||
00219                             disk_str.find("Disk_2")!=string::npos;
00220                 if(isPIB && mask) continue;
00221         
00222                 thePixelStructure.insert(pair<uint32_t,SiPixelRecHitModule*> (id,theModule));
00223               }
00224       
00225         }           
00226   }
00227 
00228   LogInfo ("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
00229 }

void SiPixelRecHitSource::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 97 of file SiPixelRecHitSource.cc.

References conf_, lat::endl(), edm::ParameterSet::getParameter(), DQMStore::save(), saveFile, and theDMBE.

00097                                     {
00098 
00099 
00100   if(saveFile){
00101     LogInfo ("PixelDQM") << " SiPixelRecHitSource::endJob - Saving Root File " << std::endl;
00102     std::string outputFile = conf_.getParameter<std::string>("outputFile");
00103     theDMBE->save( outputFile );
00104   }
00105 
00106 }


Member Data Documentation

bool SiPixelRecHitSource::bladeOn [private]

Definition at line 84 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

edm::ParameterSet SiPixelRecHitSource::conf_ [private]

Definition at line 70 of file SiPixelRecHitSource.h.

Referenced by bookMEs(), and endJob().

bool SiPixelRecHitSource::diskOn [private]

Definition at line 84 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

int SiPixelRecHitSource::eventNo [private]

Definition at line 75 of file SiPixelRecHitSource.h.

Referenced by analyze(), and beginJob().

bool SiPixelRecHitSource::isPIB [private]

Definition at line 73 of file SiPixelRecHitSource.h.

Referenced by bookMEs(), and buildStructure().

bool SiPixelRecHitSource::ladOn [private]

Definition at line 82 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

bool SiPixelRecHitSource::layOn [private]

Definition at line 82 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

bool SiPixelRecHitSource::modOn [private]

Definition at line 79 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

bool SiPixelRecHitSource::phiOn [private]

Definition at line 82 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

std::map<uint32_t,int> SiPixelRecHitSource::rechit_count [private]

Definition at line 78 of file SiPixelRecHitSource.h.

Referenced by analyze().

bool SiPixelRecHitSource::ringOn [private]

Definition at line 84 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().

bool SiPixelRecHitSource::saveFile [private]

Definition at line 72 of file SiPixelRecHitSource.h.

Referenced by endJob().

bool SiPixelRecHitSource::slowDown [private]

Definition at line 74 of file SiPixelRecHitSource.h.

Referenced by analyze().

edm::InputTag SiPixelRecHitSource::src_ [private]

Definition at line 71 of file SiPixelRecHitSource.h.

Referenced by analyze().

DQMStore* SiPixelRecHitSource::theDMBE [private]

Definition at line 76 of file SiPixelRecHitSource.h.

Referenced by bookMEs(), endJob(), and SiPixelRecHitSource().

std::map<uint32_t,SiPixelRecHitModule*> SiPixelRecHitSource::thePixelStructure [private]

Definition at line 77 of file SiPixelRecHitSource.h.

Referenced by analyze(), bookMEs(), and buildStructure().

bool SiPixelRecHitSource::twoDimOn [private]

Definition at line 80 of file SiPixelRecHitSource.h.

Referenced by analyze(), beginJob(), and bookMEs().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:32:04 2009 for CMSSW by  doxygen 1.5.4