CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/EventFilter/ScalersRawToDigi/src/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 #include "FWCore/Utilities/interface/InputTag.h"
00023 
00024 // FEDRawData 
00025 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00026 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00027 
00028 // Scalers classes
00029 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
00030 #include "DataFormats/Scalers/interface/L1TriggerScalers.h"
00031 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
00032 #include "DataFormats/Scalers/interface/Level1TriggerRates.h"
00033 #include "DataFormats/Scalers/interface/LumiScalers.h"
00034 #include "DataFormats/Scalers/interface/BeamSpotOnline.h"
00035 #include "DataFormats/Scalers/interface/DcsStatus.h"
00036 #include "DataFormats/Scalers/interface/ScalersRaw.h"
00037 
00038 class ScalersRawToDigi : public edm::EDProducer 
00039 {
00040   public:
00041     explicit ScalersRawToDigi(const edm::ParameterSet&);
00042     ~ScalersRawToDigi();
00043 
00044     virtual void produce(edm::Event&, const edm::EventSetup&);
00045 
00046   private:
00047     edm::InputTag inputTag_;
00048 };
00049 
00050 // Constructor
00051 ScalersRawToDigi::ScalersRawToDigi(const edm::ParameterSet& iConfig):
00052   inputTag_((char const *)"source")
00053 {
00054   produces<L1AcceptBunchCrossingCollection>();
00055   produces<L1TriggerScalersCollection>();
00056   produces<Level1TriggerScalersCollection>();
00057   produces<LumiScalersCollection>();
00058   produces<BeamSpotOnlineCollection>();
00059   produces<DcsStatusCollection>();
00060   if ( iConfig.exists("scalersInputTag") )
00061   {
00062     inputTag_ = iConfig.getParameter<edm::InputTag>("scalersInputTag");
00063   }
00064 }
00065 
00066 // Destructor
00067 ScalersRawToDigi::~ScalersRawToDigi() {}
00068 
00069 // Method called to produce the data 
00070 void ScalersRawToDigi::produce(edm::Event& iEvent, 
00071                                const edm::EventSetup& iSetup)
00072 {
00073   using namespace edm;
00074 
00075   // Get a handle to the FED data collection
00076   edm::Handle<FEDRawDataCollection> rawdata;
00077   iEvent.getByLabel(inputTag_, rawdata);
00078 
00079   std::auto_ptr<LumiScalersCollection> pLumi(new LumiScalersCollection());
00080 
00081   std::auto_ptr<L1TriggerScalersCollection> 
00082     pOldTrigger(new L1TriggerScalersCollection());
00083 
00084   std::auto_ptr<Level1TriggerScalersCollection> 
00085     pTrigger(new Level1TriggerScalersCollection());
00086 
00087   std::auto_ptr<L1AcceptBunchCrossingCollection> 
00088     pBunch(new L1AcceptBunchCrossingCollection());
00089 
00090   std::auto_ptr<BeamSpotOnlineCollection> pBeamSpotOnline(new BeamSpotOnlineCollection());
00091   std::auto_ptr<DcsStatusCollection> pDcsStatus(new DcsStatusCollection());
00092 
00094   const FEDRawData & fedData = rawdata->FEDData(ScalersRaw::SCALERS_FED_ID);
00095   unsigned short int length =  fedData.size();
00096   if ( length > 0 ) 
00097   {
00098     int nWords = length / 8;
00099     int nBytesExtra = 0;
00100 
00101     const ScalersEventRecordRaw_v5 * raw 
00102              = (struct ScalersEventRecordRaw_v5 *)fedData.data();
00103     if ( ( raw->version == 1 ) || ( raw->version == 2 ) )
00104     {
00105       L1TriggerScalers oldTriggerScalers(fedData.data());
00106       pOldTrigger->push_back(oldTriggerScalers);
00107       nBytesExtra = length - sizeof(struct ScalersEventRecordRaw_v1);
00108     }
00109     else if ( raw->version >= 3 )
00110     {
00111       Level1TriggerScalers triggerScalers(fedData.data());
00112       pTrigger->push_back(triggerScalers);
00113       nBytesExtra = ScalersRaw::N_BX_v2 * sizeof(unsigned long long);
00114     }
00115 
00116     LumiScalers      lumiScalers(fedData.data());
00117     pLumi->push_back(lumiScalers);
00118 
00119     if (( nBytesExtra >= 8 ) && (( nBytesExtra % 8 ) == 0 ))
00120     {
00121       unsigned long long * data = 
00122         (unsigned long long *)fedData.data();
00123 
00124       int nWordsExtra = nBytesExtra / 8;
00125       for ( int i=0; i<nWordsExtra; i++)
00126       {
00127         int index = nWords - 5 + i;
00128         L1AcceptBunchCrossing bc(i,data[index]);
00129         pBunch->push_back(bc);
00130       }
00131     }
00132 
00133     if ( raw->version >= 4 )
00134     {
00135       BeamSpotOnline beamSpotOnline(fedData.data());
00136       pBeamSpotOnline->push_back(beamSpotOnline);
00137 
00138       DcsStatus dcsStatus(fedData.data());
00139       pDcsStatus->push_back(dcsStatus);
00140     }
00141   }
00142   iEvent.put(pOldTrigger); 
00143   iEvent.put(pTrigger); 
00144   iEvent.put(pLumi); 
00145   iEvent.put(pBunch);
00146   iEvent.put(pBeamSpotOnline);
00147   iEvent.put(pDcsStatus);
00148 }
00149 
00150 // Define this as a plug-in
00151 DEFINE_FWK_MODULE(ScalersRawToDigi);