CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ScalersRawToDigi.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EventFilter/ScalersRawToDigi
4 // Class: ScalersRawToDigi
5 //
12 //
13 // Original Author: William Badgett
14 // Created: Wed Nov 14 07:47:59 CDT 2006
15 //
16 
17 #include <memory>
18 
25 
26 // FEDRawData
29 
30 // Scalers classes
39 
41 {
42  public:
43  explicit ScalersRawToDigi(const edm::ParameterSet&);
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46 
47  virtual void produce(edm::Event&, const edm::EventSetup&) override;
48 
49  private:
52 
53 };
54 
55 // Constructor
57  inputTag_((char const *)"rawDataCollector")
58 {
59  produces<L1AcceptBunchCrossingCollection>();
60  produces<L1TriggerScalersCollection>();
61  produces<Level1TriggerScalersCollection>();
62  produces<LumiScalersCollection>();
63  produces<BeamSpotOnlineCollection>();
64  produces<DcsStatusCollection>();
65  if ( iConfig.exists("scalersInputTag") )
66  {
67  inputTag_ = iConfig.getParameter<edm::InputTag>("scalersInputTag");
68  }
69  fedToken_=consumes<FEDRawDataCollection>(inputTag_);
70 
71 }
72 
73 // Destructor
75 
78  desc.add<edm::InputTag>("scalersInputTag",edm::InputTag("rawDataCollector"));
79  descriptions.add("scalersRawToDigi",desc);
80 }
81 
82 // Method called to produce the data
84  const edm::EventSetup& iSetup)
85 {
86  using namespace edm;
87 
88  // Get a handle to the FED data collection
90  iEvent.getByToken(fedToken_, rawdata);
91 
92  std::auto_ptr<LumiScalersCollection> pLumi(new LumiScalersCollection());
93 
94  std::auto_ptr<L1TriggerScalersCollection>
95  pOldTrigger(new L1TriggerScalersCollection());
96 
97  std::auto_ptr<Level1TriggerScalersCollection>
98  pTrigger(new Level1TriggerScalersCollection());
99 
100  std::auto_ptr<L1AcceptBunchCrossingCollection>
101  pBunch(new L1AcceptBunchCrossingCollection());
102 
103  std::auto_ptr<BeamSpotOnlineCollection> pBeamSpotOnline(new BeamSpotOnlineCollection());
104  std::auto_ptr<DcsStatusCollection> pDcsStatus(new DcsStatusCollection());
105 
107  const FEDRawData & fedData = rawdata->FEDData(ScalersRaw::SCALERS_FED_ID);
108  unsigned short int length = fedData.size();
109  if ( length > 0 )
110  {
111  int nWords = length / 8;
112  int nBytesExtra = 0;
113 
114  const ScalersEventRecordRaw_v6 * raw
115  = (struct ScalersEventRecordRaw_v6 *)fedData.data();
116  if ( ( raw->version == 1 ) || ( raw->version == 2 ) )
117  {
118  L1TriggerScalers oldTriggerScalers(fedData.data());
119  pOldTrigger->push_back(oldTriggerScalers);
120  nBytesExtra = length - sizeof(struct ScalersEventRecordRaw_v1);
121  }
122  else if ( raw->version >= 3 )
123  {
124  Level1TriggerScalers triggerScalers(fedData.data());
125  pTrigger->push_back(triggerScalers);
126  if ( raw->version >= 6 )
127  {
128  nBytesExtra = ScalersRaw::N_BX_v6 * sizeof(unsigned long long);
129  }
130  else
131  {
132  nBytesExtra = ScalersRaw::N_BX_v2 * sizeof(unsigned long long);
133  }
134  }
135 
136  LumiScalers lumiScalers(fedData.data());
137  pLumi->push_back(lumiScalers);
138 
139  if (( nBytesExtra >= 8 ) && (( nBytesExtra % 8 ) == 0 ))
140  {
141  unsigned long long * data =
142  (unsigned long long *)fedData.data();
143 
144  int nWordsExtra = nBytesExtra / 8;
145  for ( int i=0; i<nWordsExtra; i++)
146  {
147  int index = nWords - (nWordsExtra + 1) + i;
148  L1AcceptBunchCrossing bc(i,data[index]);
149  pBunch->push_back(bc);
150  }
151  }
152 
153  if ( raw->version >= 4 )
154  {
155  BeamSpotOnline beamSpotOnline(fedData.data());
156  pBeamSpotOnline->push_back(beamSpotOnline);
157 
158  DcsStatus dcsStatus(fedData.data());
159  pDcsStatus->push_back(dcsStatus);
160  }
161  }
162  iEvent.put(pOldTrigger);
163  iEvent.put(pTrigger);
164  iEvent.put(pLumi);
165  iEvent.put(pBunch);
166  iEvent.put(pBeamSpotOnline);
167  iEvent.put(pDcsStatus);
168 }
169 
170 // Define this as a plug-in
std::vector< BeamSpotOnline > BeamSpotOnlineCollection
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::InputTag inputTag_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
ScalersRawToDigi(const edm::ParameterSet &)
std::vector< L1AcceptBunchCrossing > L1AcceptBunchCrossingCollection
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< DcsStatus > DcsStatusCollection
Definition: DcsStatus.h:116
virtual void produce(edm::Event &, const edm::EventSetup &) override
string const
Definition: compareJSON.py:14
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< Level1TriggerScalers > Level1TriggerScalersCollection
edm::EDGetTokenT< FEDRawDataCollection > fedToken_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
std::vector< LumiScalers > LumiScalersCollection
Definition: LumiScalers.h:160
std::vector< L1TriggerScalers > L1TriggerScalersCollection
dictionary rawdata
Definition: lumiPlot.py:393