CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalCalibFEDSelector.cc
Go to the documentation of this file.
1 // user include files
4 
7 
10 
16 
18 public:
21 
22 
23 private:
24  virtual void beginJob() ;
25  virtual void produce(edm::Event&, const edm::EventSetup&);
26  virtual void endJob() ;
27 
28  // ----------member data ---------------------------
30  std::vector<int> extraFEDs_ ;
31 
32 };
33 
34 
35 
36 
38 {
39  rawInLabel_ = iConfig.getParameter<edm::InputTag>("rawInputLabel");
40  extraFEDs_ = iConfig.getParameter< std::vector<int> >("extraFEDsToKeep") ;
41  produces<FEDRawDataCollection>();
42 }
43 
45 {
46 }
47 
48 void
50 {
51 
52  std::auto_ptr<FEDRawDataCollection> producedData(new FEDRawDataCollection);
53 
55  iEvent.getByLabel(rawInLabel_,rawIn);
56 
57  std::vector<int> selFEDs;
58 
59  //--- Get the list of FEDs to be kept ---//
60  int calibType = -1 ;
63  const FEDRawData& fedData = rawIn->FEDData(i) ;
64  if ( fedData.size() < 24 ) continue ; // FED is empty
65  int value = ((const HcalDCCHeader*)(fedData.data()))->getCalibType() ;
66  if ( calibType < 0 ) {
67  calibType = value ;
68  } else {
69  if ( calibType != value )
70  edm::LogWarning("HcalCalibFEDSelector") << "Conflicting calibration types found: "
71  << calibType << " vs. " << value
72  << ". Staying with " << calibType ;
73  }
74  }
75 
76  HcalFEDList calibFeds(calibType) ;
77  selFEDs = calibFeds.getListOfFEDs() ;
78  for (unsigned int i=0; i<extraFEDs_.size(); i++) {
79  bool duplicate = false ;
80  for (unsigned int j=0; j<selFEDs.size(); j++) {
81  if (extraFEDs_.at(i) == selFEDs.at(j)) {
82  duplicate = true ;
83  break ;
84  }
85  }
86  if ( !duplicate ) selFEDs.push_back( extraFEDs_.at(i) ) ;
87  }
88 
89  // Copying:
90  const FEDRawDataCollection *rdc=rawIn.product();
91 
92  for ( int j=0; j< FEDNumbering::lastFEDId(); ++j )
93  {
94  bool rightFED=false;
95  for (uint32_t k=0; k<selFEDs.size(); k++)
96  {
97  if (j==selFEDs[k])
98  {
99  rightFED=true;
100  }
101  }
102  if (!rightFED) continue;
103  const FEDRawData & fedData = rdc->FEDData(j);
104  size_t size=fedData.size();
105 
106  if ( size > 0 )
107  {
108  // this fed has data -- lets copy it
109  FEDRawData & fedDataProd = producedData->FEDData(j);
110  if ( fedDataProd.size() != 0 ) {
111  continue;
112  }
113  fedDataProd.resize(size);
114  unsigned char *dataProd=fedDataProd.data();
115  const unsigned char *data=fedData.data();
116  // memcpy is at-least-as-fast as assignment and can be much faster
117  memcpy(dataProd, data, size);
118  }
119  }
120 
121  iEvent.put(producedData);
122 }
123 
124 
125 // ------------ method called once each job just before starting event loop ------------
127 {
128 }
129 
130 // ------------ method called once each job just after ending the event loop ------------
132 }
133 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
tuple calibType
Definition: diJetCalib.py:20
HcalCalibFEDSelector(const edm::ParameterSet &)
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
std::vector< int > extraFEDs_
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:84
int j
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
int k[5][pyjets_maxn]
static int lastFEDId()
Definition: FEDNumbering.cc:19
T const * product() const
Definition: Handle.h:74
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:29
std::vector< int > getListOfFEDs()
Definition: HcalFEDList.h:17
tuple size
Write out results.