CMS 3D CMS Logo

DTTestPulsesTask.cc

Go to the documentation of this file.
00001 /*
00002  * \file DTTestPulsesTask.cc
00003  * 
00004  * $Date: 2008/05/06 13:26:47 $
00005  * $Revision: 1.14 $
00006  * \author M. Zanetti - INFN Padova
00007  *
00008 */
00009 
00010 #include <DQM/DTMonitorModule/interface/DTTestPulsesTask.h>
00011 
00012 // Framework
00013 #include <FWCore/Framework/interface/EventSetup.h>
00014 
00015 // Digis
00016 #include <DataFormats/DTDigi/interface/DTDigi.h>
00017 #include <DataFormats/DTDigi/interface/DTDigiCollection.h>
00018 #include <DataFormats/MuonDetId/interface/DTLayerId.h>
00019 
00020 // Geometry
00021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00022 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00023 #include "Geometry/DTGeometry/interface/DTLayer.h"
00024 #include "Geometry/DTGeometry/interface/DTTopology.h"
00025 
00026 // Pedestals
00027 #include <CondFormats/DTObjects/interface/DTTtrig.h>
00028 #include <CondFormats/DataRecord/interface/DTTtrigRcd.h>
00029 #include <CondFormats/DTObjects/interface/DTRangeT0.h>
00030 #include <CondFormats/DataRecord/interface/DTRangeT0Rcd.h>
00031 #include "DQMServices/Core/interface/DQMStore.h"
00032 
00033 
00034 using namespace edm;
00035 using namespace std;
00036 
00037 
00038 DTTestPulsesTask::DTTestPulsesTask(const edm::ParameterSet& ps){
00039 
00040 
00041   cout<<"[DTTestPulseTask]: Constructor"<<endl;
00042 
00043   parameters = ps;
00044 
00045 
00046   t0sPeakRange = make_pair( parameters.getUntrackedParameter<int>("t0sRangeLowerBound", -100), 
00047                             parameters.getUntrackedParameter<int>("t0sRangeUpperBound", 100));
00048 
00049   
00050   dbe = edm::Service<DQMStore>().operator->();
00051 
00052 }
00053 
00054 DTTestPulsesTask::~DTTestPulsesTask(){
00055 
00056   cout <<"[DTTestPulsesTask]: analyzed " << nevents << " events" << endl;
00057  
00058 }
00059 
00060 
00061 void DTTestPulsesTask::beginJob(const edm::EventSetup& context){
00062 
00063   cout<<"[DTTestPulsesTask]: BeginJob"<<endl;
00064 
00065   nevents = 0;
00066 
00067   // Get the geometry
00068   context.get<MuonGeometryRecord>().get(muonGeom);
00069 
00070   // Get the pedestals tTrig (always get it, even if the TPRange is taken from conf)
00071   //context.get<DTRangeT0Rcd>().get(t0RangeMap);
00072 
00073 }
00074 
00075 
00076 void DTTestPulsesTask::bookHistos(const DTLayerId& dtLayer, string folder, string histoTag) {
00077 
00078 
00079   stringstream wheel; wheel << dtLayer.wheel(); 
00080   stringstream station; station << dtLayer.station();   
00081   stringstream sector; sector << dtLayer.sector();      
00082   stringstream superLayer; superLayer << dtLayer.superlayer();  
00083   stringstream layer; layer << dtLayer.layer(); 
00084 
00085   cout<<"[DTTestPulseTask]: booking"<<endl;
00086 
00087   // TP Profiles  
00088   if ( folder == "TPProfile" ) {
00089 
00090     const int nWires = muonGeom->layer(DTLayerId(dtLayer.wheel(),
00091                                                  dtLayer.station(),
00092                                                  dtLayer.sector(),
00093                                                  dtLayer.superlayer(),
00094                                                  dtLayer.layer()))->specificTopology().channels();
00095 
00096     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00097                           "/Station" + station.str() +
00098                           "/Sector" + sector.str() + 
00099                           "/SuperLayer" + superLayer.str() + 
00100                           "/" +folder);
00101     
00102     string histoName = histoTag 
00103       + "_W" + wheel.str() 
00104       + "_St" + station.str() 
00105       + "_Sec" + sector.str() 
00106       + "_SL" + superLayer.str() 
00107       + "_L" + layer.str();
00108 
00109     // Setting the range 
00110     if ( parameters.getUntrackedParameter<bool>("readDB", false) ) {
00111       t0RangeMap->slRangeT0( dtLayer.superlayerId() , t0sPeakRange.first, t0sPeakRange.second);
00112     }
00113     
00114     
00115     cout<<"t0sRangeLowerBound "<<t0sPeakRange.first<<"; "
00116         <<"t0sRangeUpperBound "<<t0sPeakRange.second<<endl;
00117     
00118     
00119     testPulsesProfiles[int(DTLayerId(dtLayer.wheel(),
00120                                    dtLayer.station(),
00121                                    dtLayer.sector(),
00122                                    dtLayer.superlayer(),
00123                                    dtLayer.layer()).rawId())] =
00124       dbe->bookProfile(histoName,histoName,
00125                        nWires, 0, nWires, // Xaxis: channels
00126                        t0sPeakRange.first - t0sPeakRange.second, t0sPeakRange.first, t0sPeakRange.second); // Yaxis: times
00127   }
00128 
00129   // TP Occupancies
00130   else if ( folder == "TPOccupancy" ) {
00131 
00132     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00133                           "/Station" + station.str() +
00134                           "/Sector" + sector.str() + 
00135                           "/SuperLayer" + superLayer.str() + 
00136                           "/" +folder);
00137     
00138     string histoName = histoTag 
00139       + "_W" + wheel.str() 
00140       + "_St" + station.str() 
00141       + "_Sec" + sector.str() 
00142       + "_SL" + superLayer.str() 
00143       + "_L" + layer.str();
00144     
00145     const int nWires = muonGeom->layer(DTLayerId(dtLayer.wheel(),
00146                                                  dtLayer.station(),
00147                                                  dtLayer.sector(),
00148                                                  dtLayer.superlayer(),
00149                                                  dtLayer.layer()))->specificTopology().channels();
00150     
00151     testPulsesOccupancies[int(DTLayerId(dtLayer.wheel(),
00152                                         dtLayer.station(),
00153                                         dtLayer.sector(),
00154                                         dtLayer.superlayer(),
00155                                         dtLayer.layer()).rawId())] =
00156       dbe->book1D(histoName, histoName, nWires, 0, nWires); 
00157   }
00158 
00159   // Time Box per Chamber
00160   else if ( folder == "TPTimeBox" ) {
00161 
00162     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00163                           "/Station" + station.str() +
00164                           "/Sector" + sector.str() + 
00165                           "/" +folder);
00166     
00167     string histoName = histoTag 
00168       + "_W" + wheel.str() 
00169       + "_St" + station.str() 
00170       + "_Sec" + sector.str(); 
00171     
00172     testPulsesTimeBoxes[int( DTLayerId(dtLayer.wheel(),
00173                                        dtLayer.station(),
00174                                        dtLayer.sector(),
00175                                        dtLayer.superlayer(),
00176                                        dtLayer.layer()).chamberId().rawId())] = 
00177       dbe->book1D(histoName, histoName, 10000, 0, 10000); // Overview of the TP (and noise) times
00178   }
00179 
00180 }
00181 
00182 
00183 void DTTestPulsesTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00184 
00185   nevents++;
00186   
00187   edm::Handle<DTDigiCollection> dtdigis;
00188   e.getByLabel("dtunpacker", dtdigis);
00189   
00190   DTDigiCollection::DigiRangeIterator dtLayerId_It;
00191   for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It){
00192     
00193     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
00194          digiIt!=((*dtLayerId_It).second).second; ++digiIt){
00195       
00196       // for clearness..
00197       int layerIndex = ((*dtLayerId_It).first).rawId();
00198       int chIndex = ((*dtLayerId_It).first).chamberId().rawId();
00199 
00200 
00201       if ((*digiIt).countsTDC() > t0sPeakRange.first &&
00202           (*digiIt).countsTDC() < t0sPeakRange.second ) {
00203 
00204         // Occupancies
00205         if (testPulsesOccupancies.find(layerIndex) != testPulsesOccupancies.end())
00206           testPulsesOccupancies.find(layerIndex)->second->Fill((*digiIt).wire());
00207         else {
00208           bookHistos( (*dtLayerId_It).first , string("TPOccupancy"), string("TestPulses") );
00209           testPulsesOccupancies.find(layerIndex)->second->Fill((*digiIt).wire());
00210         }
00211 
00212         // Profiles
00213         if (testPulsesProfiles.find(layerIndex) != testPulsesProfiles.end())
00214           testPulsesProfiles.find(layerIndex)->second->Fill((*digiIt).wire(),(*digiIt).countsTDC());
00215         else {
00216           bookHistos( (*dtLayerId_It).first , string("TPProfile"), string("TestPulses2D") );
00217           testPulsesProfiles.find(layerIndex)->second->Fill((*digiIt).wire(),(*digiIt).countsTDC());
00218         }
00219       }
00220 
00221       // Time Box
00222       if (testPulsesTimeBoxes.find(chIndex) != testPulsesTimeBoxes.end())
00223         testPulsesTimeBoxes.find(chIndex)->second->Fill((*digiIt).countsTDC());
00224       else {
00225         bookHistos( (*dtLayerId_It).first , string("TPTimeBox"), string("TestPulsesTB") );
00226         testPulsesTimeBoxes.find(chIndex)->second->Fill((*digiIt).countsTDC());
00227       }
00228     }
00229   }
00230 
00231 }
00232 

Generated on Tue Jun 9 17:32:41 2009 for CMSSW by  doxygen 1.5.4