CMS 3D CMS Logo

testChannel.cc
Go to the documentation of this file.
1 
10 
11 #include "TFile.h"
12 #include "TString.h"
13 
15 
18  : m_digiProducerToken(
19  consumes<EBDigiCollection>(edm::InputTag(paramSet.getParameter<std::string>("digiProducer")))),
20  m_headerProducerToken(
21  consumes<EcalRawDataCollection>(edm::InputTag(paramSet.getParameter<std::string>("headerProducer")))),
22  m_xmlFile(paramSet.getParameter<std::string>("xmlFile")),
23  m_DACmin(paramSet.getParameter<int>("DACmin")),
24  m_DACmax(paramSet.getParameter<int>("DACmax")),
25  m_RMSmax(paramSet.getParameter<double>("RMSmax")),
26  m_bestPed(paramSet.getParameter<int>("bestPed")),
27  m_xtal(paramSet.getParameter<int>("xtal")),
28  m_pedVSDAC("pedVSDAC", "pedVSDAC", 100, 150, 250, m_DACmax - m_DACmin, m_DACmin, m_DACmax),
29  m_singlePedVSDAC_1("singlePedVSDAC_1",
30  "pedVSDAC (g1) for xtal " + TString(m_xtal),
31  100,
32  150,
33  250,
34  m_DACmax - m_DACmin,
35  m_DACmin,
36  m_DACmax),
37  m_singlePedVSDAC_2("singlePedVSDAC_2",
38  "pedVSDAC (g2) for xtal " + TString(m_xtal),
39  100,
40  150,
41  250,
42  m_DACmax - m_DACmin,
43  m_DACmin,
44  m_DACmax),
45  m_singlePedVSDAC_3("singlePedVSDAC_3",
46  "pedVSDAC (g3) for xtal " + TString(m_xtal),
47  100,
48  150,
49  250,
50  m_DACmax - m_DACmin,
51  m_DACmin,
52  m_DACmax) {
53  edm::LogInfo("testChannel") << " reading "
54  << " m_DACmin: " << m_DACmin << " m_DACmax: " << m_DACmax << " m_RMSmax: " << m_RMSmax
55  << " m_bestPed: " << m_bestPed;
56 }
57 
60 
62 void testChannel::beginJob() { LogDebug("testChannel") << "entering beginJob ..."; }
63 
66  LogDebug("testChannel") << "entering analyze ...";
67 
68  // get the headers
69  // (one header for each supermodule)
70  const edm::Handle<EcalRawDataCollection> &DCCHeaders = event.getHandle(m_headerProducerToken);
71  if (!DCCHeaders.isValid()) {
72  edm::LogError("testChannel") << "Error! can't get the product for EcalRawDataCollection";
73  }
74 
75  std::map<int, int> DACvalues;
76 
77  // loop over the headers
78  for (EcalRawDataCollection::const_iterator headerItr = DCCHeaders->begin(); headerItr != DCCHeaders->end();
79  ++headerItr) {
80  EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings();
81  DACvalues[getHeaderSMId(headerItr->id())] = settings.ped_offset;
82  // std::cout << "DCCid: " << headerItr->id () << "" ;
83  // std::cout << "Ped offset DAC: " << settings.ped_offset << "" ;
84  }
85 
86  // get the digis
87  // (one digi for each crystal)
88  const edm::Handle<EBDigiCollection> &pDigis = event.getHandle(m_digiProducerToken);
89  if (!pDigis.isValid()) {
90  edm::LogError("testChannel") << "Error! can't get the product for EBDigiCollection";
91  }
92 
93  // loop over the digis
94  for (EBDigiCollection::const_iterator itdigi = pDigis->begin(); itdigi != pDigis->end(); ++itdigi) {
95  EBDataFrame df(*itdigi);
96  int gainId = df.sample(0).gainId();
97  int crystalId = EBDetId(itdigi->id()).ic();
98  int smId = EBDetId(itdigi->id()).ism();
99 
100  edm::LogInfo("testChannel") << "channel " << event.id() << "\tcry: " << crystalId << "\tG: " << gainId
101  << "\tDAC: " << DACvalues[smId];
102 
103  // loop over the samples
104  for (int iSample = 0; iSample < EBDataFrame::MAXSAMPLES; ++iSample) {
105  edm::LogInfo("testChannel") << "\t`-->" << df.sample(iSample).adc();
106  m_pedVSDAC.Fill(df.sample(iSample).adc(), DACvalues[smId]);
107  if (crystalId == m_xtal) {
108  if (gainId == 1)
109  m_singlePedVSDAC_1.Fill(df.sample(iSample).adc(), DACvalues[smId]);
110  if (gainId == 2)
111  m_singlePedVSDAC_2.Fill(df.sample(iSample).adc(), DACvalues[smId]);
112  if (gainId == 3)
113  m_singlePedVSDAC_3.Fill(df.sample(iSample).adc(), DACvalues[smId]);
114  }
115  } // loop over the samples
116  } // loop over the digis
117 }
118 
121  char ccout[80];
122  sprintf(ccout, "out_%d.root", m_xtal);
123  TFile out(ccout, "RECREATE");
124  out.cd();
125  m_pedVSDAC.Write();
126  m_singlePedVSDAC_1.Write();
127  m_singlePedVSDAC_2.Write();
128  m_singlePedVSDAC_3.Write();
129  TProfile *profilo1 = m_singlePedVSDAC_1.ProfileX();
130  TProfile *profilo2 = m_singlePedVSDAC_2.ProfileX();
131  TProfile *profilo3 = m_singlePedVSDAC_3.ProfileX();
132  profilo1->Write("singleProfile_1");
133  profilo2->Write("singleProfile_2");
134  profilo3->Write("singleProfile_3");
135  out.Close();
136 }
137 
138 /*
139 void testChannel::writeDb (EcalCondDBInterface* econn,
140  MonRunIOV* moniov)
141 {}
142 */
143 
144 int testChannel::getHeaderSMId(const int headerId) {
145  // PG FIXME temporary solution
146  // PG FIXME check it is consistent with the TB!
147  return 1;
148 }
149 
151 
153 
void endJob(void) override
EndJob.
Definition: testChannel.cc:120
std::vector< T >::const_iterator const_iterator
void analyze(edm::Event const &event, edm::EventSetup const &eventSetup) override
! Analyze
Definition: testChannel.cc:65
const int m_bestPed
Definition: testChannel.h:69
TH2F m_pedVSDAC
Definition: testChannel.h:73
Log< level::Error, false > LogError
int getHeaderSMId(const int headerId)
Definition: testChannel.cc:144
const double m_RMSmax
Definition: testChannel.h:68
void unsubscribe(void)
Definition: testChannel.cc:154
~testChannel() override
Destructor.
Definition: testChannel.cc:59
void beginJob() override
BeginJob.
Definition: testChannel.cc:62
TH2F m_singlePedVSDAC_3
Definition: testChannel.h:76
TH2F m_singlePedVSDAC_1
Definition: testChannel.h:74
void subscribe(void)
Subscribe/Unsubscribe to Monitoring Elements.
Definition: testChannel.cc:150
const_iterator begin() const
const int m_DACmin
name of the xml file to be saved
Definition: testChannel.h:66
const_iterator end() const
const_iterator end() const
Log< level::Info, false > LogInfo
const int m_DACmax
Definition: testChannel.h:67
const edm::EDGetTokenT< EBDigiCollection > m_digiProducerToken
Definition: testChannel.h:61
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
const_iterator begin() const
The iterator returned can not safely be used across threads.
const int m_xtal
Definition: testChannel.h:71
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
testChannel(const edm::ParameterSet &ps)
Constructor.
Definition: testChannel.cc:17
static constexpr int MAXSAMPLES
Definition: EcalDataFrame.h:48
int ism(int ieta, int iphi)
Definition: EcalPyUtils.cc:49
TH2F m_singlePedVSDAC_2
Definition: testChannel.h:75
void subscribeNew(void)
Definition: testChannel.cc:152
const edm::EDGetTokenT< EcalRawDataCollection > m_headerProducerToken
Token to access digis.
Definition: testChannel.h:62
Definition: event.py:1
#define LogDebug(id)