CMS 3D CMS Logo

Public Member Functions | Private Attributes

SiPixelRecHitSource Class Reference

#include <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 ()
virtual void beginRun (const edm::Run &, 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 firstRun
bool isPIB
bool ladOn
bool layOn
bool modOn
bool phiOn
std::map< uint32_t, int > rechit_count
bool reducedSet
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 54 of file SiPixelRecHitSource.h.


Constructor & Destructor Documentation

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

Definition at line 50 of file SiPixelRecHitSource.cc.

References cmsCodeRules::cppFunctionSkipper::operator, and theDMBE.

                                                                       :
  conf_(iConfig),
  src_( conf_.getParameter<edm::InputTag>( "src" ) ),
  saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
  isPIB( conf_.getUntrackedParameter<bool>("isPIB",false) ),
  slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
  modOn( conf_.getUntrackedParameter<bool>("modOn",true) ),
  twoDimOn( conf_.getUntrackedParameter<bool>("twoDimOn",true) ),
  reducedSet( conf_.getUntrackedParameter<bool>("reducedSet",false) ),
  ladOn( conf_.getUntrackedParameter<bool>("ladOn",false) ), 
  layOn( conf_.getUntrackedParameter<bool>("layOn",false) ), 
  phiOn( conf_.getUntrackedParameter<bool>("phiOn",false) ), 
  ringOn( conf_.getUntrackedParameter<bool>("ringOn",false) ), 
  bladeOn( conf_.getUntrackedParameter<bool>("bladeOn",false) ), 
  diskOn( conf_.getUntrackedParameter<bool>("diskOn",false) )
{
   theDMBE = edm::Service<DQMStore>().operator->();
   LogInfo ("PixelDQM") << "SiPixelRecHitSource::SiPixelRecHitSource: Got DQM BackEnd interface"<<endl;
}
SiPixelRecHitSource::~SiPixelRecHitSource ( )

Definition at line 71 of file SiPixelRecHitSource.cc.

References thePixelStructure.

{
   // do anything here that needs to be done at desctruction time
  // (e.g. close files, deallocate resources etc.)
  LogInfo ("PixelDQM") << "SiPixelRecHitSource::~SiPixelRecHitSource: Destructor"<<endl;
  std::map<uint32_t,SiPixelRecHitModule*>::iterator struct_iter;
  for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++){
    delete struct_iter->second;
    struct_iter->second = 0;
  }
}

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 123 of file SiPixelRecHitSource.cc.

References edmNew::DetSet< T >::begin(), bladeOn, diskOn, edmNew::DetSet< T >::end(), eventNo, edm::Event::getByLabel(), ladOn, layOn, match(), modOn, phiOn, rechit_count, reducedSet, ringOn, slowDown, mathSSE::sqrt(), src_, thePixelStructure, twoDimOn, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().

{
  eventNo++;
  //cout << eventNo << endl;
  // get input data
  edm::Handle<SiPixelRecHitCollection>  recHitColl;
  iEvent.getByLabel( src_, recHitColl );

  std::map<uint32_t,SiPixelRecHitModule*>::iterator struct_iter;
  for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++) {
    uint32_t TheID = (*struct_iter).first;

    SiPixelRecHitCollection::const_iterator match = recHitColl->find(TheID);
    
      // if( pixelrechitRangeIteratorBegin == pixelrechitRangeIteratorEnd) {cout << "oops" << endl;}
      float rechit_x = 0;
      float rechit_y = 0;
      int rechit_count = 0;

      if (match != recHitColl->end()) {
        SiPixelRecHitCollection::DetSet pixelrechitRange = *match;
        SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
        SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.end();
        SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;

       for ( ; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) 
         {
          

          rechit_count++;
          //cout << TheID << endl;
          SiPixelRecHit::ClusterRef const& clust = pixeliter->cluster();
          int sizeX = (*clust).sizeX();
          //cout << sizeX << endl;
          int sizeY = (*clust).sizeY();
          //cout << sizeY << endl;
          LocalPoint lp = pixeliter->localPosition();
          rechit_x = lp.x();
          rechit_y = lp.y();
          
          LocalError lerr = pixeliter->localPositionError();
          float lerr_x = sqrt(lerr.xx());
          float lerr_y = sqrt(lerr.yy());
          //std::cout << "errors " << lerr_x << " " << lerr_y << std::endl;
          //cout << "hh" << endl;
          (*struct_iter).second->fill(rechit_x, rechit_y, sizeX, sizeY, lerr_x, lerr_y, 
                                      modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn, 
                                      twoDimOn, reducedSet);
          //cout << "ii" << endl;
        
        }
      }
      if(rechit_count > 0) (*struct_iter).second->nfill(rechit_count, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
    
  }

  // slow down...
  if(slowDown) usleep(10000);
  
}
void SiPixelRecHitSource::beginJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 84 of file SiPixelRecHitSource.cc.

References firstRun.

                                  {
  firstRun = true;
}
void SiPixelRecHitSource::beginRun ( const edm::Run r,
edm::EventSetup const &  iSetup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 89 of file SiPixelRecHitSource.cc.

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

                                                                              {

  LogInfo ("PixelDQM") << " SiPixelRecHitSource::beginJob - Initialisation ... " << std::endl;
  LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/" 
            << layOn << "/" << phiOn << std::endl;
  LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/" 
            << ringOn << std::endl;
  LogInfo ("PixelDQM") << "2DIM IS " << twoDimOn << "\n";
  
  if(firstRun){
    eventNo = 0;
    // Build map
    buildStructure(iSetup);
    // Book Monitoring Elements
    bookMEs();
    firstRun = false;
  }
}
void SiPixelRecHitSource::bookMEs ( ) [virtual]

Create folder tree and book histograms

Definition at line 251 of file SiPixelRecHitSource.cc.

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

Referenced by beginRun().

                                 {
  
  std::map<uint32_t,SiPixelRecHitModule*>::iterator struct_iter;
    
  SiPixelFolderOrganizer theSiPixelFolder;
  
  for(struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++){
    
    if(modOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first)){
        (*struct_iter).second->book( conf_,0,twoDimOn, reducedSet);
      } else {
        if(!isPIB) throw cms::Exception("LogicError")
          << "[SiPixelDigiSource::bookMEs] Creation of DQM folder failed";
      }
    }
    if(ladOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,1)){
        (*struct_iter).second->book( conf_,1,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
      }
    }
    if(layOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,2)){
        (*struct_iter).second->book( conf_,2,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH LAYER-FOLDER\n";
      }
    }
    if(phiOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,3)){
        (*struct_iter).second->book( conf_,3,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH PHI-FOLDER\n";
      }
    }
    if(bladeOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,4)){
        (*struct_iter).second->book( conf_,4,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
      }
    }
    if(diskOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,5)){
        (*struct_iter).second->book( conf_,5,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH DISK-FOLDER\n";
      }
    }
    if(ringOn){
      if(theSiPixelFolder.setModuleFolder((*struct_iter).first,6)){
        (*struct_iter).second->book( conf_,6,twoDimOn, reducedSet);
        } else {
        LogDebug ("PixelDQM") << "PROBLEM WITH RING-FOLDER\n";
      }
    }

  }

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

Definition at line 187 of file SiPixelRecHitSource.cc.

References PixelEndcapName::bladeName(), PixelEndcapName::diskName(), edm::EventSetup::get(), PixelEndcapName::halfCylinder(), isPIB, LogDebug, PixelEndcapName::pannelName(), GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, PixelEndcapName::plaquetteName(), DetId::rawId(), DetId::subdetId(), and thePixelStructure.

Referenced by beginRun().

                                                                   {

  LogInfo ("PixelDQM") <<" SiPixelRecHitSource::buildStructure" ;
  edm::ESHandle<TrackerGeometry> pDD;
  iSetup.get<TrackerDigiGeometryRecord>().get( pDD );

  LogVerbatim ("PixelDQM") << " *** Geometry node for TrackerGeom is  "<<&(*pDD)<<std::endl;
  LogVerbatim ("PixelDQM") << " *** I have " << pDD->dets().size() <<" detectors"<<std::endl;
  LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<" types"<<std::endl;
  
  for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
    
    if(dynamic_cast<PixelGeomDetUnit*>((*it))!=0){

      DetId detId = (*it)->geographicalId();
      // const GeomDetUnit      * geoUnit = pDD->idToDetUnit( detId );
      //const PixelGeomDetUnit * pixDet  = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);

          
          
              // SiPixelRecHitModule *theModule = new SiPixelRecHitModule(id, rechit_x, rechit_y, x_res, y_res, x_pull, y_pull);
        
        
            if((detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) ||
               (detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))){ 
              uint32_t id = detId();
              SiPixelRecHitModule* theModule = new SiPixelRecHitModule(id);
              if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
                if(isPIB) continue;
                LogDebug ("PixelDQM") << " ---> Adding Barrel Module " <<  detId.rawId() << endl;
                thePixelStructure.insert(pair<uint32_t,SiPixelRecHitModule*> (id,theModule));
                
              } else if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
                LogDebug ("PixelDQM") << " ---> Adding Endcap Module " <<  detId.rawId() << endl;
                PixelEndcapName::HalfCylinder side = PixelEndcapName(DetId(id)).halfCylinder();
                int disk   = PixelEndcapName(DetId(id)).diskName();
                int blade  = PixelEndcapName(DetId(id)).bladeName();
                int panel  = PixelEndcapName(DetId(id)).pannelName();
                int module = PixelEndcapName(DetId(id)).plaquetteName();

                char sside[80];  sprintf(sside,  "HalfCylinder_%i",side);
                char sdisk[80];  sprintf(sdisk,  "Disk_%i",disk);
                char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
                char spanel[80]; sprintf(spanel, "Panel_%i",panel);
                char smodule[80];sprintf(smodule,"Module_%i",module);
                std::string side_str = sside;
                std::string disk_str = sdisk;
                bool mask = side_str.find("HalfCylinder_1")!=string::npos||
                            side_str.find("HalfCylinder_2")!=string::npos||
                            side_str.find("HalfCylinder_4")!=string::npos||
                            disk_str.find("Disk_2")!=string::npos;
                if(isPIB && mask) continue;
        
                thePixelStructure.insert(pair<uint32_t,SiPixelRecHitModule*> (id,theModule));
              }
            }
        }           
  }

  LogInfo ("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
}
void SiPixelRecHitSource::endJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 109 of file SiPixelRecHitSource.cc.

References conf_, edm::ParameterSet::getParameter(), download_sqlite_cfg::outputFile, DQMStore::save(), saveFile, and theDMBE.

                                    {


  if(saveFile){
    LogInfo ("PixelDQM") << " SiPixelRecHitSource::endJob - Saving Root File " << std::endl;
    std::string outputFile = conf_.getParameter<std::string>("outputFile");
    theDMBE->save( outputFile );
  }

}

Member Data Documentation

Definition at line 85 of file SiPixelRecHitSource.h.

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

Definition at line 70 of file SiPixelRecHitSource.h.

Referenced by bookMEs(), and endJob().

Definition at line 85 of file SiPixelRecHitSource.h.

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

Definition at line 75 of file SiPixelRecHitSource.h.

Referenced by analyze(), and beginRun().

Definition at line 87 of file SiPixelRecHitSource.h.

Referenced by beginJob(), and beginRun().

Definition at line 73 of file SiPixelRecHitSource.h.

Referenced by bookMEs(), and buildStructure().

Definition at line 83 of file SiPixelRecHitSource.h.

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

Definition at line 83 of file SiPixelRecHitSource.h.

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

Definition at line 79 of file SiPixelRecHitSource.h.

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

Definition at line 83 of file SiPixelRecHitSource.h.

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

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

Definition at line 78 of file SiPixelRecHitSource.h.

Referenced by analyze().

Definition at line 81 of file SiPixelRecHitSource.h.

Referenced by analyze(), and bookMEs().

Definition at line 85 of file SiPixelRecHitSource.h.

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

Definition at line 72 of file SiPixelRecHitSource.h.

Referenced by endJob().

Definition at line 74 of file SiPixelRecHitSource.h.

Referenced by analyze().

Definition at line 71 of file SiPixelRecHitSource.h.

Referenced by analyze().

Definition at line 76 of file SiPixelRecHitSource.h.

Referenced by endJob(), and SiPixelRecHitSource().

Definition at line 77 of file SiPixelRecHitSource.h.

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

Definition at line 80 of file SiPixelRecHitSource.h.

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