CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SubdetFEDSelector.cc
Go to the documentation of this file.
1 
3 
5 {
6  getEcal_=iConfig.getParameter<bool>("getECAL");
7  getStrip_=iConfig.getParameter<bool>("getSiStrip");
8  getPixel_=iConfig.getParameter<bool>("getSiPixel");
9  getHcal_=iConfig.getParameter<bool>("getHCAL");
10  getMuon_=iConfig.getParameter<bool>("getMuon");
11  getTrigger_=iConfig.getParameter<bool>("getTrigger");
12 
13  rawInLabel_=iConfig.getParameter<edm::InputTag>("rawInputLabel");
14 
15  produces<FEDRawDataCollection>();
16 
17 }
18 
20 {
21 }
22 
23 void
25 {
26 
27  std::auto_ptr<FEDRawDataCollection> producedData(new FEDRawDataCollection);
28 
30  iEvent.getByLabel(rawInLabel_,rawIn);
31 
32  std::vector<int> selFEDs;
33 
34  if (getEcal_)
35  {
37  {
38  selFEDs.push_back(i);
39  }
41  {
42  selFEDs.push_back(i);
43  }
44  }
45 
46  if (getMuon_)
47  {
49  {
50  selFEDs.push_back(i);
51  }
53  {
54  selFEDs.push_back(i);
55  }
57  {
58  selFEDs.push_back(i);
59  }
61  {
62  selFEDs.push_back(i);
63  }
65  {
66  selFEDs.push_back(i);
67  }
69  {
70  selFEDs.push_back(i);
71  }
73  {
74  selFEDs.push_back(i);
75  }
77  {
78  selFEDs.push_back(i);
79  }
80  }
81 
82 
83  if (getHcal_)
84  {
86  {
87  selFEDs.push_back(i);
88  }
89  }
90 
91 
92  if (getStrip_)
93  {
95  {
96  selFEDs.push_back(i);
97  }
98  }
99 
100 
101  if (getPixel_)
102  {
104  {
105  selFEDs.push_back(i);
106  }
107  }
108 
109  if (getTrigger_)
110  {
112  {
113  selFEDs.push_back(i);
114  }
116  {
117  selFEDs.push_back(i);
118  }
120  {
121  selFEDs.push_back(i);
122  }
124  {
125  selFEDs.push_back(i);
126  }
128  {
129  selFEDs.push_back(i);
130  }
131 
133  {
134  selFEDs.push_back(i);
135  }
136 
138  {
139  selFEDs.push_back(i);
140  }
141 
143  {
144  selFEDs.push_back(i);
145  }
146 
148  {
149  selFEDs.push_back(i);
150  }
151 
153  {
154  selFEDs.push_back(i);
155  }
157  {
158  selFEDs.push_back(i);
159  }
160 
162  {
163  selFEDs.push_back(i);
164  }
166  {
167  selFEDs.push_back(i);
168  }
170  {
171  selFEDs.push_back(i);
172  }
173 
174  }
175 
177  {
178  selFEDs.push_back(i);
179  }
180 
181 
182  // Copying:
183  const FEDRawDataCollection *rdc=rawIn.product();
184 
185  // if ( ( rawData[i].provenance()->processName() != e.processHistory().rbegin()->processName() ) )
186  // continue ; // skip all raw collections not produced by the current process
187 
188  for ( int j=0; j< FEDNumbering::MAXFEDID; ++j )
189  {
190  bool rightFED=false;
191  for (uint32_t k=0; k<selFEDs.size(); k++)
192  {
193  if (j==selFEDs[k])
194  {
195  rightFED=true;
196  }
197  }
198  if (!rightFED) continue;
199  const FEDRawData & fedData = rdc->FEDData(j);
200  size_t size=fedData.size();
201 
202  if ( size > 0 )
203  {
204  // this fed has data -- lets copy it
205  FEDRawData & fedDataProd = producedData->FEDData(j);
206  if ( fedDataProd.size() != 0 ) {
207 // std::cout << " More than one FEDRawDataCollection with data in FED ";
208 // std::cout << j << " Skipping the 2nd\n";
209  continue;
210  }
211  fedDataProd.resize(size);
212  unsigned char *dataProd=fedDataProd.data();
213  const unsigned char *data=fedData.data();
214  for ( unsigned int k=0; k<size; ++k ) {
215  dataProd[k]=data[k];
216  }
217  }
218  }
219 
220  iEvent.put(producedData);
221 
222 }
223 
224 
225 // ------------ method called once each job just before starting event loop ------------
227 }
228 
229 // ------------ method called once each job just after ending the event loop ------------
231 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual void endJob()
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:49
edm::InputTag rawInLabel_
int iEvent
Definition: GenABIO.cc:243
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void resize(size_t newsize)
Definition: FEDRawData.cc:33
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
int j
Definition: DBlmapReader.cc:9
SubdetFEDSelector(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
int k[5][pyjets_maxn]
virtual void beginJob()
T const * product() const
Definition: Handle.h:74
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:29
virtual void produce(edm::Event &, const edm::EventSetup &)
tuple size
Write out results.