CMS 3D CMS Logo

DTDigiTask.cc

Go to the documentation of this file.
00001  /*
00002  * \file DTDigiTask.cc
00003  * 
00004  * $Date: 2008/11/24 09:12:11 $
00005  * $Revision: 1.54 $
00006  * \author M. Zanetti - INFN Padova
00007  *
00008  */
00009 
00010 #include <DQM/DTMonitorModule/interface/DTDigiTask.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 #include <DataFormats/MuonDetId/interface/DTChamberId.h>
00020 
00021 // Geometry
00022 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00023 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00024 #include "Geometry/DTGeometry/interface/DTLayer.h"
00025 #include "Geometry/DTGeometry/interface/DTTopology.h"
00026 
00027 // T0s
00028 #include <CondFormats/DTObjects/interface/DTT0.h>
00029 #include <CondFormats/DataRecord/interface/DTT0Rcd.h>
00030 #include <CondFormats/DTObjects/interface/DTTtrig.h>
00031 #include <CondFormats/DataRecord/interface/DTTtrigRcd.h>
00032 
00033 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00034 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00035 
00036 #include "DQMServices/Core/interface/DQMStore.h"
00037 #include "DQMServices/Core/interface/MonitorElement.h"
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 
00040 #include <sstream>
00041 #include <math.h>
00042 
00043 using namespace edm;
00044 using namespace std;
00045 
00046 
00047 
00048 // Contructor
00049 DTDigiTask::DTDigiTask(const edm::ParameterSet& ps){
00050   // switch for the verbosity
00051   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") <<"[DTDigiTask]: Constructor"<<endl;
00052 
00053   // The label to retrieve the digis 
00054   dtDigiLabel = ps.getParameter<InputTag>("dtDigiLabel");
00055   // Read the configuration parameters
00056   maxTDCHits = ps.getUntrackedParameter<int>("maxTDCHitsPerChamber",30000);
00057   // Set to true to read the ttrig from DB (useful to determine in-time and out-of-time hits)
00058   readTTrigDB = ps.getUntrackedParameter<bool>("readDB", false);
00059   // Set to true to subtract t0 from test pulses
00060   subtractT0 = ps.getParameter<bool>("performPerWireT0Calibration");
00061   // Tmax value (TDC counts)
00062   defaultTmax = ps.getParameter<int>("defaultTmax");
00063   // Switch from static to dinamic histo booking
00064   doStaticBooking =  ps.getUntrackedParameter<bool>("staticBooking", true);
00065   // Switch for local/global runs
00066   isLocalRun = ps.getUntrackedParameter<bool>("localrun", true);
00067   // Setting for the reset of the ME after n (= ResetCycle) luminosity sections
00068   resetCycle = ps.getUntrackedParameter<int>("ResetCycle", 3);
00069   // Check the DB of noisy channels
00070   checkNoisyChannels = ps.getUntrackedParameter<bool>("checkNoisyChannels","false");
00071   // Default TTrig to be used when not reading the TTrig DB
00072   defaultTTrig = ps.getParameter<int>("defaultTtrig");
00073   inTimeHitsLowerBound = ps.getParameter<int>("inTimeHitsLowerBound");
00074   inTimeHitsUpperBound = ps.getParameter<int>("inTimeHitsUpperBound");
00075   timeBoxGranularity = ps.getUntrackedParameter<int>("timeBoxGranularity",4);
00076   tdcRescale = ps.getUntrackedParameter<int>("tdcRescale", 1);
00077 
00078   doAllHitsOccupancies = ps.getUntrackedParameter<bool>("doAllHitsOccupancies", true);
00079   doNoiseOccupancies = ps.getUntrackedParameter<bool>("doNoiseOccupancies", false);
00080   doInTimeOccupancies = ps.getUntrackedParameter<bool>("doInTimeOccupancies", false);
00081 
00082   // switch on the mode for running on test pulses (different top folder)
00083   tpMode = ps.getUntrackedParameter<bool>("testPulseMode", false);
00084   // switch on/off the filtering of synchronous noise events (cutting on the # of digis)
00085   // time-boxes and occupancy plots are not filled and summary plots are created to report the problem
00086   filterSyncNoise = ps.getUntrackedParameter<bool>("filterSyncNoise", false);
00087 
00088   dbe = edm::Service<DQMStore>().operator->();
00089 
00090   syncNumTot = 0;
00091   syncNum = 0;
00092 
00093 }
00094 
00095 
00096 
00097 // destructor
00098 DTDigiTask::~DTDigiTask(){
00099   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "DTDigiTask: analyzed " << nevents << " events" << endl;
00100 
00101 }
00102 
00103 
00104 
00105 
00106 void DTDigiTask::endJob(){
00107   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") <<"[DTDigiTask] endjob called!"<<endl;
00108   
00109 }
00110 
00111 
00112 
00113 
00114 void DTDigiTask::beginJob(const edm::EventSetup& context){
00115   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") <<"[DTDigiTask]: BeginJob"<<endl;
00116 
00117   nevents = 0;
00118 
00119   // Get the geometry
00120   context.get<MuonGeometryRecord>().get(muonGeom);
00121 }
00122 
00123 
00124 void DTDigiTask::beginRun(const edm::Run& run, const edm::EventSetup& context) {
00125   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") <<"[DTDigiTask]: begin run"<<endl;
00126 
00127   // tTrig 
00128   if (readTTrigDB) 
00129     context.get<DTTtrigRcd>().get(tTrigMap);
00130   // t0s 
00131   if (subtractT0) 
00132     context.get<DTT0Rcd>().get(t0Map);
00133   // FIXME: tMax (not yet from the DB)
00134   tMax = defaultTmax;
00135   
00136   // ----------------------------------------------------------------------
00137   if(doStaticBooking) {  // Static histo booking
00138     // book the event counter
00139     dbe->setCurrentFolder(topFolder());
00140     nEventMonitor = dbe->bookFloat("nProcessedEvents");
00141     for(int wh = -2; wh <= 2; ++wh) { // loop over wheels
00142       if(doAllHitsOccupancies) bookHistos(wh, string("Occupancies"), "OccupancyAllHits");
00143       if(doNoiseOccupancies) bookHistos(wh, string("Occupancies"), "OccupancyNoiseHits");
00144       if(doInTimeOccupancies) bookHistos(wh, string("Occupancies"), "OccupancyInTimeHits");
00145 
00146       if(filterSyncNoise) bookHistos(wh, string("SynchNoise"), "SyncNoiseEvents" );
00147       for(int st = 1; st <= 4; ++st) { // loop over stations
00148         for(int sect = 1; sect <= 14; ++sect) { // loop over sectors
00149           if((sect == 13 || sect == 14) && st != 4) continue;
00150           // Get the chamber ID
00151           const  DTChamberId dtChId(wh,st,sect);
00152 
00153           // Occupancies 
00154           if (doAllHitsOccupancies) 
00155             bookHistos(dtChId, string("Occupancies"), "OccupancyAllHits_perCh");
00156           if(doNoiseOccupancies) 
00157             bookHistos(dtChId, string("Occupancies"), "OccupancyNoise_perCh");
00158           if(doInTimeOccupancies)
00159             bookHistos(dtChId, string("Occupancies"), "OccupancyInTimeHits_perCh" );
00160 
00161 
00162 
00163 
00164           for(int sl = 1; sl <= 3; ++sl) { // Loop over SLs
00165             if(st == 4 && sl == 2) continue;
00166             const  DTSuperLayerId dtSLId(wh,st,sect,sl);
00167             if(isLocalRun) {
00168               bookHistos(dtSLId, string("TimeBoxes"), "TimeBox");
00169             } else {
00170               // TimeBoxes for different triggers
00171               bookHistos(dtSLId, string("TimeBoxes"), "TimeBoxDTonly");
00172               bookHistos(dtSLId, string("TimeBoxes"), "TimeBoxNoDT");
00173               bookHistos(dtSLId, string("TimeBoxes"), "TimeBoxDTalso");
00174             }
00175           }
00176         }
00177       }
00178     }
00179   }
00180 }
00181 
00182 
00183 
00184 
00185 void DTDigiTask::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00186   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") <<"[DTDigiTask]: Begin of LS transition"<<endl;
00187   
00188   // Reset the MonitorElements every n (= ResetCycle) Lumi Blocks
00189   if(lumiSeg.id().luminosityBlock() % resetCycle == 0) {
00190     LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask")
00191       <<"[DTDigiTask]: Reset at the LS transition : "<<lumiSeg.id().luminosityBlock()<<endl;
00192     // Loop over all ME
00193     for(map<string, map<uint32_t, MonitorElement*> > ::const_iterator histo = digiHistos.begin();
00194         histo != digiHistos.end(); histo++) {
00195       for(map<uint32_t, MonitorElement*> ::const_iterator ht = (*histo).second.begin();
00196           ht != (*histo).second.end(); ht++) {
00197         (*ht).second->Reset();
00198       }
00199     }
00200 
00201     // loop over wheel summaries
00202     for(map<string, map<int, MonitorElement*> > ::const_iterator histos = wheelHistos.begin();
00203         histos != wheelHistos.end(); ++histos) {
00204       for(map<int, MonitorElement*>::const_iterator histo = (*histos).second.begin();
00205           histo != (*histos).second.end(); ++histo) {
00206         (*histo).second->Reset();
00207       }
00208     }
00209   }
00210   
00211 }
00212 
00213 
00214 
00215 
00216 void DTDigiTask::bookHistos(const DTSuperLayerId& dtSL, string folder, string histoTag) {
00217   // set the folder
00218   stringstream wheel; wheel << dtSL.wheel();    
00219   stringstream station; station << dtSL.station();      
00220   stringstream sector; sector << dtSL.sector(); 
00221   stringstream superLayer; superLayer << dtSL.superlayer();
00222   dbe->setCurrentFolder(topFolder() + "Wheel" + wheel.str() +
00223                         "/Station" + station.str() +
00224                         "/Sector" + sector.str());
00225 
00226   // Build the histo name
00227   string histoName = histoTag 
00228     + "_W" + wheel.str() 
00229     + "_St" + station.str() 
00230     + "_Sec" + sector.str() 
00231     + "_SL" + superLayer.str(); 
00232 
00233   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
00234     << "[DTDigiTask]: booking SL histo:" << histoName
00235     << " (tag: " << histoTag
00236     << ") folder: " << topFolder() + "Wheel" + wheel.str() +
00237     "/Station" + station.str() +
00238     "/Sector" + sector.str() + "/" + folder << endl;
00239     
00240 
00241   // ttrig and rms are TDC counts
00242   if ( readTTrigDB ) 
00243     tTrigMap->get(dtSL, tTrig, tTrigRMS, DTTimeUnits::counts); 
00244   else tTrig = defaultTTrig;
00245   
00246 
00247   if ( folder == "TimeBoxes") {
00248     string histoTitle = histoName + " (TDC Counts)";
00249 
00250     if (!readTTrigDB) {
00251       int maxTDCCounts = 6400 * tdcRescale;
00252       (digiHistos[histoTag])[dtSL.rawId()] = 
00253         dbe->book1D(histoName,histoTitle, maxTDCCounts/timeBoxGranularity, 0, maxTDCCounts);
00254     }    
00255     else {
00256       (digiHistos[histoTag])[dtSL.rawId()] = 
00257         dbe->book1D(histoName,histoTitle, 3*tMax/timeBoxGranularity, tTrig-tMax, tTrig+2*tMax);
00258     }
00259   }
00260 
00261   if ( folder == "CathodPhotoPeaks" ) {
00262     dbe->setCurrentFolder(topFolder() + "Wheel" + wheel.str() +
00263                           "/Station" + station.str() +
00264                           "/Sector" + sector.str() + "/" + folder);
00265     (digiHistos[histoTag])[dtSL.rawId()] = dbe->book1D(histoName,histoName,500,0,1000);
00266   }
00267   
00268 }
00269 
00270 
00271 
00272 
00273 void DTDigiTask::bookHistos(const DTChamberId& dtCh, string folder, string histoTag) {
00274   // set the current folder
00275   stringstream wheel; wheel << dtCh.wheel();    
00276   stringstream station; station << dtCh.station();      
00277   stringstream sector; sector << dtCh.sector();
00278   dbe->setCurrentFolder(topFolder() + "Wheel" + wheel.str() +
00279                         "/Station" + station.str() +
00280                         "/Sector" + sector.str());
00281 
00282   // build the histo name
00283   string histoName = histoTag 
00284     + "_W" + wheel.str() 
00285     + "_St" + station.str() 
00286     + "_Sec" + sector.str(); 
00287 
00288 
00289   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
00290     << "[DTDigiTask]: booking chamber histo:" 
00291     << " (tag: " << histoTag
00292     << ") folder: " << topFolder() + "Wheel" + wheel.str() +
00293     "/Station" + station.str() +
00294     "/Sector" + sector.str() << endl;
00295     
00296   
00297   if (folder == "Occupancies")    {
00298     
00299     const DTChamber* dtchamber = muonGeom->chamber(dtCh);
00300     const std::vector<const DTSuperLayer*> dtSupLylist = dtchamber->superLayers();
00301     std::vector<const DTSuperLayer*>::const_iterator suly = dtSupLylist.begin();
00302     std::vector<const DTSuperLayer*>::const_iterator sulyend = dtSupLylist.end();
00303     
00304     int nWires = 0;
00305     int firstWire = 0;
00306     int nWires_max = 0;
00307     
00308     while(suly != sulyend) {
00309       const std::vector<const DTLayer*> dtLyList = (*suly)->layers();
00310       std::vector<const DTLayer*>::const_iterator ly = dtLyList.begin();
00311       std::vector<const DTLayer*>::const_iterator lyend = dtLyList.end();
00312       stringstream superLayer; superLayer << (*suly)->id().superlayer();
00313       
00314       while(ly != lyend) {
00315         nWires = muonGeom->layer((*ly)->id())->specificTopology().channels();
00316         firstWire = muonGeom->layer((*ly)->id())->specificTopology().firstChannel();
00317         stringstream layer; layer << (*ly)->id().layer();
00318         string histoName_layer = histoName + "_SL" + superLayer.str()  + "_L" + layer.str();
00319         if(histoTag == "OccupancyAllHits_perL" 
00320            || histoTag == "OccupancyNoise_perL"
00321            || histoTag == "OccupancyInTimeHits_perL")
00322           (digiHistos[histoTag])[(*ly)->id().rawId()] = dbe->book1D(histoName_layer,histoName_layer,nWires,firstWire,nWires+firstWire);
00323         ++ly;
00324         if((nWires+firstWire) > nWires_max) nWires_max = (nWires+firstWire);
00325         
00326       }
00327       ++suly;
00328     }
00329    
00330     if(histoTag != "OccupancyAllHits_perL" 
00331            && histoTag != "OccupancyNoise_perL"
00332            && histoTag != "OccupancyInTimeHits_perL"){
00333       (digiHistos[histoTag])[dtCh.rawId()] = dbe->book2D(histoName,histoName,nWires_max,1,nWires_max+1,12,0,12);
00334     
00335       for(int i=1;i<=12;i++) { 
00336         if(i<5){
00337           stringstream layer;
00338           string layer_name;
00339           layer<<i;
00340           layer>>layer_name;
00341           string label="SL1: L"+layer_name;
00342           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i,label,2);
00343         }
00344         else if(i>4 && i<9){
00345           stringstream layer;
00346           string layer_name;
00347           layer<<(i-4);
00348           layer>>layer_name;
00349           string label="SL2: L"+layer_name;
00350           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i,label,2);
00351         }
00352         else if(i>8 && i<13){
00353           stringstream layer;
00354           string layer_name;
00355           layer<<(i-8);
00356           layer>>layer_name;
00357           string label="SL3: L"+layer_name;
00358           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i,label,2);
00359         }
00360       }
00361     }
00362   }
00363 }
00364 
00365 
00366 
00367 
00368 void DTDigiTask::bookHistos(const int wheelId, string folder, string histoTag) {
00369   // Set the current folder
00370   stringstream wheel; wheel << wheelId; 
00371 
00372 
00373   // build the histo name
00374   string histoName = histoTag + "_W" + wheel.str(); 
00375 
00376 
00377   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
00378     << "[DTDigiTask]: booking wheel histo:" << histoName
00379     << " (tag: " << histoTag
00380     << ") folder: " << topFolder() + "Wheel" + wheel.str() + "/" <<endl;
00381 
00382   if(folder == "Occupancies") {
00383     dbe->setCurrentFolder(topFolder() + "Wheel" + wheel.str());
00384     string histoTitle = "# of digis per chamber WHEEL: "+wheel.str();
00385     (wheelHistos[histoTag])[wheelId] = dbe->book2D(histoName,histoTitle,12,1,13,4,1,5);
00386     (wheelHistos[histoTag])[wheelId]->setBinLabel(1,"MB1",2);
00387     (wheelHistos[histoTag])[wheelId]->setBinLabel(2,"MB2",2);
00388     (wheelHistos[histoTag])[wheelId]->setBinLabel(3,"MB3",2);
00389     (wheelHistos[histoTag])[wheelId]->setBinLabel(4,"MB4",2);
00390     (wheelHistos[histoTag])[wheelId]->setAxisTitle("sector",1);
00391   } else if(folder == "SynchNoise") {
00392     dbe->setCurrentFolder("DT/04-Noise/SynchNoise");
00393     string histoTitle = "Event rate of Syncronous Noisy events WHEEL: "+wheel.str();
00394     (wheelHistos[histoTag])[wheelId] = dbe->book2D(histoName,histoTitle,12,1,13,4,1,5);
00395     (wheelHistos[histoTag])[wheelId]->setBinLabel(1,"MB1",2);
00396     (wheelHistos[histoTag])[wheelId]->setBinLabel(2,"MB2",2);
00397     (wheelHistos[histoTag])[wheelId]->setBinLabel(3,"MB3",2);
00398     (wheelHistos[histoTag])[wheelId]->setBinLabel(4,"MB4",2);
00399     (wheelHistos[histoTag])[wheelId]->setAxisTitle("sector",1);
00400   }
00401 
00402 }
00403 
00404 
00405 
00406 // does the real job
00407 void DTDigiTask::analyze(const edm::Event& event, const edm::EventSetup& c) {
00408   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask] analyze" << endl;
00409   
00410   nevents++;
00411   nEventMonitor->Fill(nevents);
00412   if (nevents%1000 == 0) {
00413     LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask] Analyze #Run: " << event.id().run()
00414                                                  << " #Event: " << event.id().event() << endl;
00415   }
00416   
00417   // Get the ingredients from the event
00418   
00419   // Digi collection
00420   edm::Handle<DTDigiCollection> dtdigis;
00421   event.getByLabel(dtDigiLabel, dtdigis);
00422 
00423   // LTC digis
00424   if (!isLocalRun) event.getByType(ltcdigis);
00425 
00426   // Status map (for noisy channels)
00427   ESHandle<DTStatusFlag> statusMap;
00428   if(checkNoisyChannels) {
00429     LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "    get the map of noisy channels" << endl;
00430     // Get the map of noisy channels
00431     c.get<DTStatusFlagRcd>().get(statusMap);
00432   }
00433 
00434   string histoTag;
00435 
00436 
00437   // Check if the digi container is empty
00438   if(dtdigis->begin() == dtdigis->end()) {
00439     LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "Event " << nevents << " empty." << endl;
00440   }
00441 
00442   if (filterSyncNoise) { // dosync
00443     LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "     Filter Sync Noise" << endl;
00444     
00445     // Count the # of digis per chamber
00446     DTDigiCollection::DigiRangeIterator dtLayerId_It;
00447     for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); dtLayerId_It++) {
00448       DTChamberId chId = ((*dtLayerId_It).first).chamberId();
00449       if(hitMap.find(chId) == hitMap.end()) {// new chamber
00450         hitMap[chId] = 0;
00451       }
00452       hitMap[chId] += (((*dtLayerId_It).second).second - ((*dtLayerId_It).second).first);
00453     }
00454 
00455     
00456     // check chamber with # of digis above threshold and flag them as noisy
00457     for (map<DTChamberId,int>::iterator iter = hitMap.begin(); iter != hitMap.end(); iter++) {
00458       if((iter->second) > maxTDCHits) { 
00459         syncNoisyChambers.insert(iter->first);
00460         nSynchNoiseEvents[iter->first]++;// FIXME check and optimize
00461         // FIXME: log noisy event in this chamber
00462         wheelHistos["SyncNoiseEvents"][(*iter).first.wheel()]->setBinContent((*iter).first.sector(),(*iter).first.station(),
00463                                                                              (double)nSynchNoiseEvents[iter->first]/(double)nevents); 
00464       }
00465     }
00466     
00467     // clear the map of # of digis per chamber: not needed anymore
00468     hitMap.clear();
00469 
00470     // log
00471     if (syncNoisyChambers.size() != 0) {
00472       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask]** Synch Noise in event: " << nevents
00473                                                       << " noisy time-boxes and occupancy will not be filled!" << endl; 
00474       syncNumTot++;
00475       syncNum++;
00476     }
00477 
00478     if (nevents%1000 == 0) {
00479       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask") << (syncNumTot*100./nevents) << "% sync noise events since the beginning \n"
00480                                                       << (syncNum*0.1) << "% sync noise events in the last 1000 events " << endl;
00481       syncNum = 0;
00482     }
00483   }
00484 
00485   bool isSyncNoisy = false;
00486 
00487   DTDigiCollection::DigiRangeIterator dtLayerId_It;
00488   for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It) { // Loop over layers
00489     isSyncNoisy = false;
00490     // check if chamber labeled as synch noisy
00491     if (filterSyncNoise) {
00492       DTChamberId chId = ((*dtLayerId_It).first).chamberId();
00493       if(syncNoisyChambers.find(chId) != syncNoisyChambers.end()) {
00494         isSyncNoisy = true;
00495       }
00496     }
00497 
00498     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
00499          digiIt!=((*dtLayerId_It).second).second; ++digiIt) { // Loop over all digis
00500 
00501         bool isNoisy = false;
00502         bool isFEMasked = false;
00503         bool isTDCMasked = false;
00504         bool isTrigMask = false;
00505         bool isDead = false;
00506         bool isNohv = false;
00507         if(checkNoisyChannels) {
00508           const DTWireId wireId(((*dtLayerId_It).first), (*digiIt).wire());
00509           statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00510         }
00511         
00512         
00513 
00514         // Get the useful IDs
00515         const  DTSuperLayerId dtSLId = ((*dtLayerId_It).first).superlayerId();
00516         uint32_t indexSL = dtSLId.rawId();
00517         const  DTChamberId dtChId = dtSLId.chamberId(); 
00518         uint32_t indexCh = dtChId.rawId();
00519         int layer_number=((*dtLayerId_It).first).layer();
00520         int superlayer_number=dtSLId.superlayer();
00521         const  DTLayerId dtLId = (*dtLayerId_It).first;
00522         
00523         // Read the ttrig DB or set a rough value from config
00524         // ttrig and rms are TDC counts
00525         if (readTTrigDB)
00526           tTrigMap->get( ((*dtLayerId_It).first).superlayerId(),
00527                          tTrig, tTrigRMS, DTTimeUnits::counts); 
00528         else tTrig = defaultTTrig;
00529         
00530         int inTimeHitsLowerBoundCorr = int(round(tTrig)) - inTimeHitsLowerBound;
00531         int inTimeHitsUpperBoundCorr = int(round(tTrig)) + tMax + inTimeHitsUpperBound;
00532         
00533         float t0; float t0RMS;
00534         int tdcTime = (*digiIt).countsTDC();
00535         
00536         if (subtractT0) {
00537           const DTWireId dtWireId(((*dtLayerId_It).first), (*digiIt).wire());
00538           // t0s and rms are TDC counts
00539           t0Map->get(dtWireId, t0, t0RMS, DTTimeUnits::counts) ;
00540           tdcTime += int(round(t0));
00541         }
00542 
00543         
00544 
00545         // Fill Time-Boxes
00546         // NOTE: avoid to fill TB and PhotoPeak with noise. Occupancy are filled anyway
00547         if (( !isNoisy ) && (!isSyncNoisy)) { // Discard noisy channels
00548           // TimeBoxes per SL
00549           histoTag = "TimeBox" + triggerSource();
00550           if (digiHistos[histoTag].find(indexSL) == digiHistos[histoTag].end())
00551             bookHistos( dtSLId, string("TimeBoxes"), histoTag );
00552           (digiHistos.find(histoTag)->second).find(indexSL)->second->Fill(tdcTime);
00553           
00554           // FIXME: remove the time distribution for the after-pulses     
00555           // 2nd - 1st (CathodPhotoPeak) per SL
00556           //      if ( (*digiIt).number() == 1 ) {
00557             
00558           //        DTDigiCollection::const_iterator firstDigiIt = digiIt;
00559           //        firstDigiIt--;
00560             
00561           //        histoTag = "CathodPhotoPeak";
00562           //        if (digiHistos[histoTag].find(indexSL) == digiHistos[histoTag].end())
00563           //          bookHistos( dtSLId, string("CathodPhotoPeaks"), histoTag );
00564           //        (digiHistos.find(histoTag)->second).find(indexSL)->second->Fill((*digiIt).countsTDC()-
00565           //                                                                        (*firstDigiIt).countsTDC());
00566           //      }
00567         }
00568 
00569         // Fill Occupancies
00570         if (!isSyncNoisy) { // Discard synch noisy channels 
00571 
00572           if (doAllHitsOccupancies) { // fill occupancies for all hits
00573             //Occupancies per chamber & layer
00574             histoTag = "OccupancyAllHits_perCh";
00575             map<uint32_t, MonitorElement*>::const_iterator mappedHisto =
00576               digiHistos[histoTag].find(indexCh);
00577             if (mappedHisto == digiHistos[histoTag].end()) { // dynamic booking
00578               bookHistos(dtChId, string("Occupancies"), histoTag);
00579               mappedHisto = digiHistos[histoTag].find(indexCh);
00580             }
00581             mappedHisto->second->Fill((*digiIt).wire(),(layer_number+(superlayer_number-1)*4)-1);
00582 
00583             
00584             // Fill the chamber occupancy
00585             histoTag = "OccupancyAllHits";
00586             map<int, MonitorElement*>::const_iterator histoPerWheel =
00587               wheelHistos[histoTag].find(dtChId.wheel());
00588             if(histoPerWheel ==  wheelHistos[histoTag].end()) { // dynamic booking
00589               bookHistos(dtChId.wheel(), string("Occupancies"), histoTag);
00590               histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
00591             }
00592             histoPerWheel->second->Fill(dtChId.sector(),dtChId.station()); // FIXME: normalize to # of layers
00593            
00594             
00595           } 
00596 
00597           if(doNoiseOccupancies) { // fill occupancies for hits before the ttrig
00598             if (tdcTime < inTimeHitsLowerBoundCorr ) { 
00599               // FIXME: what about tdcTime > inTimeHitsUpperBoundCorr ???
00600 
00601               // Noise: Before tTrig
00602               //Occupancies Noise per chamber & layer
00603               histoTag = "OccupancyNoise_perCh";
00604               map<uint32_t, MonitorElement*>::const_iterator mappedHisto =
00605                 digiHistos[histoTag].find(indexCh);
00606               if(mappedHisto == digiHistos[histoTag].end()) {
00607                 bookHistos(dtChId, string("Occupancies"), histoTag);
00608                 mappedHisto = digiHistos[histoTag].find(indexCh);
00609               }
00610               mappedHisto->second->Fill((*digiIt).wire(),
00611                                         (layer_number+(superlayer_number-1)*4)-1);
00612 
00613               // Fill the chamber occupancy
00614               histoTag = "OccupancyNoise";
00615               map<int, MonitorElement*>::const_iterator histoPerWheel =
00616                 wheelHistos[histoTag].find(dtChId.wheel());
00617               if(histoPerWheel ==  wheelHistos[histoTag].end()) { // dynamic booking
00618                 bookHistos(dtChId.wheel(), string("Occupancies"), histoTag);
00619                 histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
00620               }
00621               histoPerWheel->second->Fill(dtChId.sector(),dtChId.station()); // FIXME: normalize to # of layers
00622 
00623             } 
00624           }
00625           
00626           if(doInTimeOccupancies) { // fill occpunacies for in-time hits only
00627             if (tdcTime > inTimeHitsLowerBoundCorr && tdcTime < inTimeHitsUpperBoundCorr) { 
00628               // Physical hits: within the time window  
00629 
00630               //Occupancies Signal per chamber & layer
00631               histoTag = "OccupancyInTimeHits_perCh";
00632               map<uint32_t, MonitorElement*>::const_iterator mappedHisto =
00633                 digiHistos[histoTag].find(indexCh);
00634               if(mappedHisto == digiHistos[histoTag].end()) {
00635                 bookHistos(dtChId, string("Occupancies"), histoTag);
00636                 mappedHisto = digiHistos[histoTag].find(indexCh);
00637               }
00638               mappedHisto->second->Fill((*digiIt).wire(),
00639                                         (layer_number+(superlayer_number-1)*4)-1);
00640 
00641               // Fill the chamber occupancy
00642               histoTag = "OccupancyInTimeHits";
00643               map<int, MonitorElement*>::const_iterator histoPerWheel =
00644                 wheelHistos[histoTag].find(dtChId.wheel());
00645               if(histoPerWheel ==  wheelHistos[histoTag].end()) { // dynamic booking
00646                 bookHistos(dtChId.wheel(), string("Occupancies"), histoTag);
00647                 histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
00648               }
00649               histoPerWheel->second->Fill(dtChId.sector(),dtChId.station()); // FIXME: normalize to # of layers
00650 
00651             }
00652           }
00653         }
00654     }
00655   }
00656   
00657   syncNoisyChambers.clear();
00658 }
00659 
00660 
00661 string DTDigiTask::triggerSource() {
00662 
00663   string l1ASource;
00664 
00665   if (!isLocalRun) {
00666     for (std::vector<LTCDigi>::const_iterator ltc_it = ltcdigis->begin(); ltc_it != ltcdigis->end(); ltc_it++){
00667       int otherTriggerSum=0;
00668       for (int i = 1; i < 6; i++)
00669         otherTriggerSum += int((*ltc_it).HasTriggered(i));
00670       
00671       if ((*ltc_it).HasTriggered(0) && otherTriggerSum == 0) 
00672         l1ASource = "DTonly";
00673       else if (!(*ltc_it).HasTriggered(0))
00674         l1ASource = "NoDT";
00675       else if ((*ltc_it).HasTriggered(0) && otherTriggerSum > 0)
00676         l1ASource = "DTalso";
00677     }
00678   }
00679 
00680   return l1ASource;
00681 
00682 }
00683 
00684 
00685 string DTDigiTask::topFolder() const {
00686   if(tpMode) return string("DT/10-TestPulses/");
00687   return string("DT/01-Digi/");
00688 }
00689 
00690 
00691 
00692 

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