CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/DQM/DTMonitorModule/src/DTTestPulsesTask.cc

Go to the documentation of this file.
00001 /*
00002  * \file DTTestPulsesTask.cc
00003  * 
00004  * $Date: 2010/01/05 10:14:40 $
00005  * $Revision: 1.16 $
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(){
00062 
00063   cout<<"[DTTestPulsesTask]: BeginJob"<<endl;
00064 
00065   nevents = 0;
00066 
00067 }
00068 
00069 void DTTestPulsesTask::beginRun(const edm::Run& run, const edm::EventSetup& context) {
00070 
00071   // Get the geometry
00072   context.get<MuonGeometryRecord>().get(muonGeom);
00073 
00074   // Get the pedestals tTrig (always get it, even if the TPRange is taken from conf)
00075   //context.get<DTRangeT0Rcd>().get(t0RangeMap);
00076 
00077 }
00078 
00079 void DTTestPulsesTask::bookHistos(const DTLayerId& dtLayer, string folder, string histoTag) {
00080 
00081 
00082   stringstream wheel; wheel << dtLayer.wheel(); 
00083   stringstream station; station << dtLayer.station();   
00084   stringstream sector; sector << dtLayer.sector();      
00085   stringstream superLayer; superLayer << dtLayer.superlayer();  
00086   stringstream layer; layer << dtLayer.layer(); 
00087 
00088   cout<<"[DTTestPulseTask]: booking"<<endl;
00089 
00090   // TP Profiles  
00091   if ( folder == "TPProfile" ) {
00092 
00093     const int nWires = muonGeom->layer(DTLayerId(dtLayer.wheel(),
00094                                                  dtLayer.station(),
00095                                                  dtLayer.sector(),
00096                                                  dtLayer.superlayer(),
00097                                                  dtLayer.layer()))->specificTopology().channels();
00098 
00099     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00100                           "/Station" + station.str() +
00101                           "/Sector" + sector.str() + 
00102                           "/SuperLayer" + superLayer.str() + 
00103                           "/" +folder);
00104     
00105     string histoName = histoTag 
00106       + "_W" + wheel.str() 
00107       + "_St" + station.str() 
00108       + "_Sec" + sector.str() 
00109       + "_SL" + superLayer.str() 
00110       + "_L" + layer.str();
00111 
00112     // Setting the range 
00113     if ( parameters.getUntrackedParameter<bool>("readDB", false) ) {
00114       t0RangeMap->slRangeT0( dtLayer.superlayerId() , t0sPeakRange.first, t0sPeakRange.second);
00115     }
00116     
00117     
00118     cout<<"t0sRangeLowerBound "<<t0sPeakRange.first<<"; "
00119         <<"t0sRangeUpperBound "<<t0sPeakRange.second<<endl;
00120     
00121     
00122     testPulsesProfiles[int(DTLayerId(dtLayer.wheel(),
00123                                    dtLayer.station(),
00124                                    dtLayer.sector(),
00125                                    dtLayer.superlayer(),
00126                                    dtLayer.layer()).rawId())] =
00127       dbe->bookProfile(histoName,histoName,
00128                        nWires, 0, nWires, // Xaxis: channels
00129                        t0sPeakRange.first - t0sPeakRange.second, t0sPeakRange.first, t0sPeakRange.second); // Yaxis: times
00130   }
00131 
00132   // TP Occupancies
00133   else if ( folder == "TPOccupancy" ) {
00134 
00135     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00136                           "/Station" + station.str() +
00137                           "/Sector" + sector.str() + 
00138                           "/SuperLayer" + superLayer.str() + 
00139                           "/" +folder);
00140     
00141     string histoName = histoTag 
00142       + "_W" + wheel.str() 
00143       + "_St" + station.str() 
00144       + "_Sec" + sector.str() 
00145       + "_SL" + superLayer.str() 
00146       + "_L" + layer.str();
00147     
00148     const int nWires = muonGeom->layer(DTLayerId(dtLayer.wheel(),
00149                                                  dtLayer.station(),
00150                                                  dtLayer.sector(),
00151                                                  dtLayer.superlayer(),
00152                                                  dtLayer.layer()))->specificTopology().channels();
00153     
00154     testPulsesOccupancies[int(DTLayerId(dtLayer.wheel(),
00155                                         dtLayer.station(),
00156                                         dtLayer.sector(),
00157                                         dtLayer.superlayer(),
00158                                         dtLayer.layer()).rawId())] =
00159       dbe->book1D(histoName, histoName, nWires, 0, nWires); 
00160   }
00161 
00162   // Time Box per Chamber
00163   else if ( folder == "TPTimeBox" ) {
00164 
00165     dbe->setCurrentFolder("DT/DTTestPulsesTask/Wheel" + wheel.str() +
00166                           "/Station" + station.str() +
00167                           "/Sector" + sector.str() + 
00168                           "/" +folder);
00169     
00170     string histoName = histoTag 
00171       + "_W" + wheel.str() 
00172       + "_St" + station.str() 
00173       + "_Sec" + sector.str(); 
00174     
00175     testPulsesTimeBoxes[int( DTLayerId(dtLayer.wheel(),
00176                                        dtLayer.station(),
00177                                        dtLayer.sector(),
00178                                        dtLayer.superlayer(),
00179                                        dtLayer.layer()).chamberId().rawId())] = 
00180       dbe->book1D(histoName, histoName, 10000, 0, 10000); // Overview of the TP (and noise) times
00181   }
00182 
00183 }
00184 
00185 
00186 void DTTestPulsesTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00187 
00188   nevents++;
00189   
00190   edm::Handle<DTDigiCollection> dtdigis;
00191   e.getByLabel("dtunpacker", dtdigis);
00192   
00193   DTDigiCollection::DigiRangeIterator dtLayerId_It;
00194   for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It){
00195     
00196     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
00197          digiIt!=((*dtLayerId_It).second).second; ++digiIt){
00198       
00199       // for clearness..
00200       int layerIndex = ((*dtLayerId_It).first).rawId();
00201       int chIndex = ((*dtLayerId_It).first).chamberId().rawId();
00202 
00203 
00204       if ((int)(*digiIt).countsTDC() > t0sPeakRange.first &&
00205           (int)(*digiIt).countsTDC() < t0sPeakRange.second ) {
00206 
00207         // Occupancies
00208         if (testPulsesOccupancies.find(layerIndex) != testPulsesOccupancies.end())
00209           testPulsesOccupancies.find(layerIndex)->second->Fill((*digiIt).wire());
00210         else {
00211           bookHistos( (*dtLayerId_It).first , string("TPOccupancy"), string("TestPulses") );
00212           testPulsesOccupancies.find(layerIndex)->second->Fill((*digiIt).wire());
00213         }
00214 
00215         // Profiles
00216         if (testPulsesProfiles.find(layerIndex) != testPulsesProfiles.end())
00217           testPulsesProfiles.find(layerIndex)->second->Fill((*digiIt).wire(),(*digiIt).countsTDC());
00218         else {
00219           bookHistos( (*dtLayerId_It).first , string("TPProfile"), string("TestPulses2D") );
00220           testPulsesProfiles.find(layerIndex)->second->Fill((*digiIt).wire(),(*digiIt).countsTDC());
00221         }
00222       }
00223 
00224       // Time Box
00225       if (testPulsesTimeBoxes.find(chIndex) != testPulsesTimeBoxes.end())
00226         testPulsesTimeBoxes.find(chIndex)->second->Fill((*digiIt).countsTDC());
00227       else {
00228         bookHistos( (*dtLayerId_It).first , string("TPTimeBox"), string("TestPulsesTB") );
00229         testPulsesTimeBoxes.find(chIndex)->second->Fill((*digiIt).countsTDC());
00230       }
00231     }
00232   }
00233 
00234 }
00235