CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ECALRegFEDSelector.cc
Go to the documentation of this file.
1 
6 
8 {
9  tok_seed_ = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("regSeedLabel"));
10  delta_=iConfig.getParameter<double>("delta");
11 
12  tok_raw_ = consumes<FEDRawDataCollection>(iConfig.getParameter<edm::InputTag>("rawInputLabel"));
13 
15 
16  produces<FEDRawDataCollection>();
17  produces<EcalListOfFEDS>();
18 
19  for (int p=0; p<1200; p++)
20  {
21  fedSaved[p]=false;
22  }
23 }
24 
25 
27 {
28 }
29 
30 
32 {
33  for (int p=0; p<1200; p++)
34  {
35  fedSaved[p]=false;
36  }
37 
38  std::auto_ptr<FEDRawDataCollection> producedData(new FEDRawDataCollection);
39 
40  std::auto_ptr<EcalListOfFEDS> fedList(new EcalListOfFEDS);
41 
43  iEvent.getByToken(tok_seed_,trigSeedTrks);
44 
45  std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
46  trigSeedTrks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
47 
49  iEvent.getByToken(tok_raw_,rawIn);
50 
51  // std::vector<int> EC_FED_IDs;
52 
53  for (uint32_t p=0; p<isoPixTrackRefs.size(); p++)
54  {
55  double etaObj_=isoPixTrackRefs[p]->track()->eta();
56  double phiObj_=isoPixTrackRefs[p]->track()->phi();
57 
58  EcalEtaPhiRegion ecEtaPhi(etaObj_-delta_,etaObj_+delta_,phiObj_-delta_,phiObj_+delta_);
59 
60  const std::vector<int> EC_FED_IDs=ec_mapping->GetListofFEDs(ecEtaPhi);
61 
62  const FEDRawDataCollection *rdc=rawIn.product();
63 
64  for ( int j=0; j< FEDNumbering::MAXFEDID; j++ )
65  {
66  bool rightFED=false;
67  for (uint32_t k=0; k<EC_FED_IDs.size(); k++)
68  {
69  if (j==EcalRegionCabling::fedIndex(EC_FED_IDs[k]))
70  {
71  if (!fedSaved[j])
72  {
73  fedList->AddFED(j);
74  rightFED=true;
75  fedSaved[j]=true;
76  }
77  }
78  }
80  {
81  fedSaved[j]=true;
82  rightFED=true;
83  }
84  if (!rightFED) continue;
85  const FEDRawData & fedData = rdc->FEDData(j);
86  size_t size=fedData.size();
87 
88  if ( size > 0 )
89  {
90  // this fed has data -- lets copy it
91  FEDRawData & fedDataProd = producedData->FEDData(j);
92  if ( fedDataProd.size() != 0 )
93  {
94 // std::cout << " More than one FEDRawDataCollection with data in FED ";
95 // std::cout << j << " Skipping the 2nd\n";
96  continue;
97  }
98  fedDataProd.resize(size);
99  unsigned char *dataProd=fedDataProd.data();
100  const unsigned char *data=fedData.data();
101  for ( unsigned int k=0; k<size; ++k )
102  {
103  dataProd[k]=data[k];
104  }
105  }
106 
107  }
108  }
109 
110  iEvent.put(producedData);
111  iEvent.put(fedList);
112 
113 }
114 
115 
117 }
118 
119 
121 }
T getParameter(std::string const &) const
edm::EDGetTokenT< FEDRawDataCollection > tok_raw_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_seed_
ECALRegFEDSelector(const edm::ParameterSet &)
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
const EcalElectronicsMapping * ec_mapping
static int fedIndex(const uint32_t index)
int iEvent
Definition: GenABIO.cc:230
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void resize(size_t newsize)
Definition: FEDRawData.cc:32
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
int j
Definition: DBlmapReader.cc:9
T const * product() const
Definition: Handle.h:81
virtual void produce(edm::Event &, const edm::EventSetup &)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
std::vector< int > GetListofFEDs(const EcalEtaPhiRegion &region) const
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
tuple size
Write out results.