CMS 3D CMS Logo

ScalersRawToDigi.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EventFilter/ScalersRawToDigi
00004 // Class:      ScalersRawToDigi
00005 // 
00012 //
00013 // Original Author:  William Badgett
00014 //         Created:  Wed Nov 14 07:47:59 CDT 2006
00015 //
00016 
00017 #include <memory>
00018 
00019 #include "FWCore/Framework/interface/EDProducer.h"
00020 #include "FWCore/Framework/interface/Event.h"
00021 #include "FWCore/Framework/interface/MakerMacros.h"
00022 
00023 // FEDRawData 
00024 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00025 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00026 
00027 // Scalers classes
00028 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
00029 #include "DataFormats/Scalers/interface/L1TriggerScalers.h"
00030 #include "DataFormats/Scalers/interface/LumiScalers.h"
00031 #include "DataFormats/Scalers/interface/ScalersRaw.h"
00032 
00033 class ScalersRawToDigi : public edm::EDProducer 
00034 {
00035   public:
00036     explicit ScalersRawToDigi(const edm::ParameterSet&);
00037     ~ScalersRawToDigi();
00038 
00039     virtual void produce(edm::Event&, const edm::EventSetup&);
00040 };
00041 
00042 // Constructor
00043 ScalersRawToDigi::ScalersRawToDigi(const edm::ParameterSet& iConfig)
00044 {
00045   produces<L1AcceptBunchCrossingCollection>();
00046   produces<L1TriggerScalersCollection>();
00047   produces<LumiScalersCollection>();
00048 }
00049 
00050 // Destructor
00051 ScalersRawToDigi::~ScalersRawToDigi() {}
00052 
00053 // Method called to produce the data 
00054 void ScalersRawToDigi::produce(edm::Event& iEvent, 
00055                                const edm::EventSetup& iSetup)
00056 {
00057   using namespace edm;
00058 
00059   // Get a handle to the FED data collection
00060   edm::Handle<FEDRawDataCollection> rawdata;
00061   iEvent.getByLabel("source" , rawdata);
00062 
00063   std::auto_ptr<LumiScalersCollection> pLumi(new LumiScalersCollection());
00064 
00065   std::auto_ptr<L1TriggerScalersCollection> 
00066     pTrigger(new L1TriggerScalersCollection());
00067 
00068   std::auto_ptr<L1AcceptBunchCrossingCollection> 
00069     pBunch(new L1AcceptBunchCrossingCollection());
00070 
00072   const FEDRawData & fedData = rawdata->FEDData(ScalersRaw::SCALERS_FED_ID);
00073   unsigned short int length =  fedData.size();
00074   if ( length > 0 ) 
00075   {
00076     L1TriggerScalers triggerScalers(fedData.data());
00077     pTrigger->push_back(triggerScalers);
00078     iEvent.put(pTrigger); 
00079 
00080     LumiScalers      lumiScalers(fedData.data());
00081     pLumi->push_back(lumiScalers);
00082     iEvent.put(pLumi); 
00083 
00084     int nWords = length / 8;
00085     int nBytesExtra = length - sizeof(struct ScalersEventRecordRaw_v1);
00086     if (( nBytesExtra >= 8 ) && (( nBytesExtra % 8 ) == 0 ))
00087     {
00088       unsigned long long * data = 
00089         (unsigned long long *)fedData.data();
00090 
00091       int nWordsExtra = nBytesExtra / 8;
00092       for ( int i=0; i<nWordsExtra; i++)
00093       {
00094         int index = nWords - 5 + i;
00095         L1AcceptBunchCrossing bc(i,data[index]);
00096         pBunch->push_back(bc);
00097       }
00098       iEvent.put(pBunch); 
00099     }
00100   }
00101 }
00102 
00103 // Define this as a plug-in
00104 DEFINE_FWK_MODULE(ScalersRawToDigi);

Generated on Tue Jun 9 17:34:47 2009 for CMSSW by  doxygen 1.5.4