CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
SiStripBaselineAnalyzer Class Reference

#include <Validation/SiStripAnalyzer/src/SiStripBaselineAnalyzer.cc>

Inheritance diagram for SiStripBaselineAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 SiStripBaselineAnalyzer (const edm::ParameterSet &)
 
 ~SiStripBaselineAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob () override
 
virtual void endJob () override
 

Private Attributes

uint16_t actualModule_
 
TCanvas * Canvas_
 
edm::Service< TFileServicefs_
 
TH1F * h1APVCM_
 
TH1F * h1BadAPVperEvent_
 
TH1F * h1Baseline_
 
TH1F * h1Clusters_
 
TH1F * h1Pedestals_
 
TH1F * h1ProcessedRawDigis_
 
uint16_t nModuletoDisplay_
 
std::vector< int > pedestals
 
edm::ESHandle< SiStripPedestalspedestalsHandle
 
uint32_t peds_cache_id
 
bool plotAPVCM_
 
bool plotBaseline_
 
bool plotBaselinePoints_
 
bool plotClusters_
 
bool plotPedestals_
 
bool plotRawDigi_
 
edm::InputTag srcAPVCM_
 
edm::InputTag srcBaseline_
 
edm::InputTag srcBaselinePoints_
 
edm::InputTag srcProcessedRawDigi_
 
std::auto_ptr
< SiStripPedestalsSubtractor
subtractorPed_
 
std::vector< TH1F > vBaselineHisto_
 
std::vector< TH1F > vBaselinePointsHisto_
 
std::vector< TH1F > vClusterHisto_
 
std::vector< TH1F > vProcessedRawDigiHisto_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 74 of file SiStripBaselineAnalyzer.cc.

Constructor & Destructor Documentation

SiStripBaselineAnalyzer::SiStripBaselineAnalyzer ( const edm::ParameterSet conf)
explicit

Definition at line 123 of file SiStripBaselineAnalyzer.cc.

References SiStripRawProcessingFactory::create_SubtractorPed(), fs_, edm::ParameterSet::getParameter(), h1APVCM_, h1BadAPVperEvent_, h1Pedestals_, TFileService::make(), nModuletoDisplay_, plotAPVCM_, plotBaseline_, plotBaselinePoints_, plotClusters_, plotPedestals_, plotRawDigi_, srcAPVCM_, srcBaseline_, srcBaselinePoints_, srcProcessedRawDigi_, and subtractorPed_.

123  {
124 
125  srcBaseline_ = conf.getParameter<edm::InputTag>( "srcBaseline" );
126  srcBaselinePoints_ = conf.getParameter<edm::InputTag>( "srcBaselinePoints" );
127  srcProcessedRawDigi_ = conf.getParameter<edm::InputTag>( "srcProcessedRawDigi" );
128  srcAPVCM_ = conf.getParameter<edm::InputTag>( "srcAPVCM" );
130  nModuletoDisplay_ = conf.getParameter<uint32_t>( "nModuletoDisplay" );
131  plotClusters_ = conf.getParameter<bool>( "plotClusters" );
132  plotBaseline_ = conf.getParameter<bool>( "plotBaseline" );
133  plotBaselinePoints_ = conf.getParameter<bool>( "plotBaselinePoints" );
134  plotRawDigi_ = conf.getParameter<bool>( "plotRawDigi" );
135  plotAPVCM_ = conf.getParameter<bool>( "plotAPVCM" );
136  plotPedestals_ = conf.getParameter<bool>( "plotPedestals" );
137 
138  h1BadAPVperEvent_ = fs_->make<TH1F>("BadAPV/Event","BadAPV/Event", 2001, -0.5, 2000.5);
139  h1BadAPVperEvent_->SetXTitle("# Modules with Bad APVs");
140  h1BadAPVperEvent_->SetYTitle("Entries");
141  h1BadAPVperEvent_->SetLineWidth(2);
142  h1BadAPVperEvent_->SetLineStyle(2);
143 
144  h1APVCM_ = fs_->make<TH1F>("APV CM","APV CM", 2048, -1023.5, 1023.5);
145  h1APVCM_->SetXTitle("APV CM [adc]");
146  h1APVCM_->SetYTitle("Entries");
147  h1APVCM_->SetLineWidth(2);
148  h1APVCM_->SetLineStyle(2);
149 
150  h1Pedestals_ = fs_->make<TH1F>("Pedestals","Pedestals", 2048, -1023.5, 1023.5);
151  h1Pedestals_->SetXTitle("Pedestals [adc]");
152  h1Pedestals_->SetYTitle("Entries");
153  h1Pedestals_->SetLineWidth(2);
154  h1Pedestals_->SetLineStyle(2);
155 
156 
157 
158 }
T getParameter(std::string const &) const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::Service< TFileService > fs_
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed_
static std::auto_ptr< SiStripPedestalsSubtractor > create_SubtractorPed(const edm::ParameterSet &)
SiStripBaselineAnalyzer::~SiStripBaselineAnalyzer ( )

Definition at line 161 of file SiStripBaselineAnalyzer.cc.

162 {
163 
164 
165 
166 }

Member Function Documentation

void SiStripBaselineAnalyzer::analyze ( const edm::Event e,
const edm::EventSetup es 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 169 of file SiStripBaselineAnalyzer.cc.

References actualModule_, ecalMGPA::adc(), sistrip::APV, edm::DetSet< T >::begin(), edm::DetSetVector< T >::begin(), edmNew::DetSetVector< T >::begin(), gather_cfg::cout, edm::DetSet< T >::end(), edm::DetSetVector< T >::end(), edmNew::DetSetVector< T >::end(), edm::EventID::event(), event(), fs_, edm::EventSetup::get(), edm::Event::getByLabel(), h1APVCM_, h1BadAPVperEvent_, h1Baseline_, h1Clusters_, h1Pedestals_, h1ProcessedRawDigis_, i, edm::EventBase::id(), edmNew::DetSetVector< T >::id(), TFileService::mkdir(), nModuletoDisplay_, pedestals, pedestalsHandle, peds_cache_id, plotAPVCM_, plotBaseline_, plotBaselinePoints_, plotClusters_, plotPedestals_, plotRawDigi_, edm::EventID::run(), DTTTrigCorrFirst::run, gather_cfg::runs, srcAPVCM_, srcBaseline_, srcProcessedRawDigi_, and subtractorPed_.

170 {
171  using namespace edm;
172  if(plotPedestals_&&actualModule_ ==0){
173  uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
174  if(p_cache_id != peds_cache_id) {
176  peds_cache_id = p_cache_id;
177  }
178 
179 
180  std::vector<uint32_t> detIdV;
181  pedestalsHandle->getDetIds(detIdV);
182 
183  for(uint32_t i=0; i < detIdV.size(); ++i){
184  pedestals.clear();
185  SiStripPedestals::Range pedestalsRange = pedestalsHandle->getRange(detIdV[i]);
186  pedestals.resize((pedestalsRange.second- pedestalsRange.first)*8/10);
187  pedestalsHandle->allPeds(pedestals, pedestalsRange);
188  for(uint32_t it=0; it < pedestals.size(); ++it) h1Pedestals_->Fill(pedestals[it]);
189  }
190  }
191 
192  if(plotAPVCM_){
194  edm::InputTag CMLabel("siStripZeroSuppression:APVCM");
195  e.getByLabel(srcAPVCM_,moduleCM);
196 
197  edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itCMDetSetV =moduleCM->begin();
198  for (; itCMDetSetV != moduleCM->end(); ++itCMDetSetV){
200  for(;itCM != itCMDetSetV->end(); ++itCM) h1APVCM_->Fill(itCM->adc());
201  }
202  }
203 
204  if(!plotRawDigi_) return;
205  subtractorPed_->init(es);
206 
207 
208 
210  e.getByLabel(srcProcessedRawDigi_,moduleRawDigi);
211 
213  if(plotBaseline_) e.getByLabel(srcBaseline_, moduleBaseline);
214 
215  edm::Handle<edm::DetSetVector<SiStripDigi> > moduleBaselinePoints;
216  if(plotBaselinePoints_) e.getByLabel(srcBaseline_, moduleBaselinePoints);
217 
219  if(plotClusters_){
220  edm::InputTag clusLabel("siStripClusters");
221  e.getByLabel(clusLabel, clusters);
222  }
223 
224  char detIds[20];
225  char evs[20];
226  char runs[20];
227 
228 
229  TFileDirectory sdProcessedRawDigis_= fs_->mkdir("ProcessedRawDigis");
230  TFileDirectory sdBaseline_= fs_->mkdir("Baseline");
231  TFileDirectory sdBaselinePoints_= fs_->mkdir("BaselinePoints");
232  TFileDirectory sdClusters_= fs_->mkdir("Clusters");
233 
234 
236  if(plotBaseline_) itDSBaseline = moduleBaseline->begin();
237  edm::DetSetVector<SiStripRawDigi>::const_iterator itRawDigis = moduleRawDigi->begin();
238 
239  uint32_t NBabAPVs = moduleRawDigi->size();
240  std::cout<< "Number of module with HIP in this event: " << NBabAPVs << std::endl;
241  h1BadAPVperEvent_->Fill(NBabAPVs);
242 
243  for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis) {
244  if(actualModule_ > nModuletoDisplay_) return;
245  uint32_t detId = itRawDigis->id;
246 
247  if(plotBaseline_){
248  //std::cout << "bas id: " << itDSBaseline->id << " raw id: " << detId << std::endl;
249  if(itDSBaseline->id != detId){
250  std::cout << "Collections out of Synch. Something of fishy is going on ;-)" << std::endl;
251  return;
252  }
253  }
254 
255 
256  actualModule_++;
257  uint32_t event = e.id().event();
258  uint32_t run = e.id().run();
259  //std::cout << "processing module N: " << actualModule_<< " detId: " << detId << " event: "<< event << std::endl;
260 
261 
262 
263  edm::DetSet<SiStripRawDigi>::const_iterator itRaw = itRawDigis->begin();
264  bool restAPV[6] = {0,0,0,0,0,0};
265  int strip =0, totADC=0;
266  int minAPVRes = 7, maxAPVRes = -1;
267  for(;itRaw != itRawDigis->end(); ++itRaw, ++strip){
268  float adc = itRaw->adc();
269  totADC+= adc;
270  if(strip%127 ==0){
271  //std::cout << "totADC " << totADC << std::endl;
272  int APV = strip/128;
273  if(totADC!= 0){
274  restAPV[APV] = true;
275  totADC =0;
276  if(APV>maxAPVRes) maxAPVRes = APV;
277  if(APV<minAPVRes) minAPVRes = APV;
278  }
279  }
280  }
281 
282  uint16_t bins =768;
283  float minx = -0.5, maxx=767.5;
284  if(minAPVRes !=7){
285  minx = minAPVRes * 128 -0.5;
286  maxx = maxAPVRes * 128 + 127.5;
287  bins = maxx-minx;
288  }
289 
290  sprintf(detIds,"%ul", detId);
291  sprintf(evs,"%ul", event);
292  sprintf(runs,"%ul", run);
293  char* dHistoName = Form("Id:%s_run:%s_ev:%s",detIds, runs, evs);
294  h1ProcessedRawDigis_ = sdProcessedRawDigis_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
295 
296  if(plotBaseline_){
297  h1Baseline_ = sdBaseline_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
298  h1Baseline_->SetXTitle("strip#");
299  h1Baseline_->SetYTitle("ADC");
300  h1Baseline_->SetMaximum(1024.);
301  h1Baseline_->SetMinimum(-300.);
302  h1Baseline_->SetLineWidth(2);
303  h1Baseline_->SetLineStyle(2);
304  h1Baseline_->SetLineColor(2);
305  }
306 
307  if(plotClusters_){
308  h1Clusters_ = sdClusters_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
309 
310  h1Clusters_->SetXTitle("strip#");
311  h1Clusters_->SetYTitle("ADC");
312  h1Clusters_->SetMaximum(1024.);
313  h1Clusters_->SetMinimum(-300.);
314  h1Clusters_->SetLineWidth(2);
315  h1Clusters_->SetLineStyle(2);
316  h1Clusters_->SetLineColor(3);
317  }
318 
319  h1ProcessedRawDigis_->SetXTitle("strip#");
320  h1ProcessedRawDigis_->SetYTitle("ADC");
321  h1ProcessedRawDigis_->SetMaximum(1024.);
322  h1ProcessedRawDigis_->SetMinimum(-300.);
323  h1ProcessedRawDigis_->SetLineWidth(2);
324 
325  std::vector<int16_t> ProcessedRawDigis(itRawDigis->size());
326  subtractorPed_->subtract( *itRawDigis, ProcessedRawDigis);
327 
329  if(plotBaseline_) itBaseline = itDSBaseline->begin();
330  std::vector<int16_t>::const_iterator itProcessedRawDigis;
331 
332  strip =0;
333  for(itProcessedRawDigis = ProcessedRawDigis.begin();itProcessedRawDigis != ProcessedRawDigis.end(); ++itProcessedRawDigis){
334  if(restAPV[strip/128]){
335  float adc = *itProcessedRawDigis;
336  h1ProcessedRawDigis_->Fill(strip, adc);
337  if(plotBaseline_){
338  h1Baseline_->Fill(strip, itBaseline->adc());
339  ++itBaseline;
340  }
341  }
342  ++strip;
343  }
344 
345  if(plotBaseline_) ++itDSBaseline;
346  if(plotClusters_){
348  for ( ; itClusters != clusters->end(); ++itClusters ){
349  for ( edmNew::DetSet<SiStripCluster>::const_iterator clus = itClusters->begin(); clus != itClusters->end(); ++clus){
350  if(itClusters->id() == detId){
351  int firststrip = clus->firstStrip();
352  //std::cout << "Found cluster in detId " << detId << " " << firststrip << " " << clus->amplitudes().size() << " -----------------------------------------------" << std::endl;
353  strip=0;
354  for( std::vector<uint8_t>::const_iterator itAmpl = clus->amplitudes().begin(); itAmpl != clus->amplitudes().end(); ++itAmpl){
355  h1Clusters_->Fill(firststrip+strip, *itAmpl);
356  ++strip;
357  }
358  }
359  }
360  }
361  }
362  }
363 
364 }
int adc(sample_type sample)
get the ADC sample (12 bits)
RunNumber_t run() const
Definition: EventID.h:39
iterator end()
Definition: DetSet.h:60
EventNumber_t event() const
Definition: EventID.h:41
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
data_type const * const_iterator
Definition: DetSetNew.h:30
std::pair< ContainerIterator, ContainerIterator > Range
id_type id(size_t cell) const
tuple runs
Definition: gather_cfg.py:87
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
edm::Service< TFileService > fs_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:402
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
iterator begin()
Definition: DetSet.h:59
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
const T & get() const
Definition: EventSetup.h:55
edm::EventID id() const
Definition: EventBase.h:56
tuple cout
Definition: gather_cfg.py:121
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:350
edm::ESHandle< SiStripPedestals > pedestalsHandle
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed_
const_iterator begin(bool update=false) const
void SiStripBaselineAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 368 of file SiStripBaselineAnalyzer.cc.

References actualModule_.

369 {
370 
371 
372 actualModule_ =0;
373 
374 
375 }
void SiStripBaselineAnalyzer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 379 of file SiStripBaselineAnalyzer.cc.

379  {
380 
381 }

Member Data Documentation

uint16_t SiStripBaselineAnalyzer::actualModule_
private

Definition at line 119 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and beginJob().

TCanvas* SiStripBaselineAnalyzer::Canvas_
private

Definition at line 112 of file SiStripBaselineAnalyzer.cc.

edm::Service<TFileService> SiStripBaselineAnalyzer::fs_
private

Definition at line 102 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

TH1F* SiStripBaselineAnalyzer::h1APVCM_
private

Definition at line 109 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

TH1F* SiStripBaselineAnalyzer::h1BadAPVperEvent_
private

Definition at line 104 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

TH1F* SiStripBaselineAnalyzer::h1Baseline_
private

Definition at line 107 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

TH1F* SiStripBaselineAnalyzer::h1Clusters_
private

Definition at line 108 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

TH1F* SiStripBaselineAnalyzer::h1Pedestals_
private

Definition at line 110 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

TH1F* SiStripBaselineAnalyzer::h1ProcessedRawDigis_
private

Definition at line 106 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

uint16_t SiStripBaselineAnalyzer::nModuletoDisplay_
private

Definition at line 118 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

std::vector<int> SiStripBaselineAnalyzer::pedestals
private

Definition at line 87 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

edm::ESHandle<SiStripPedestals> SiStripBaselineAnalyzer::pedestalsHandle
private

Definition at line 86 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

uint32_t SiStripBaselineAnalyzer::peds_cache_id
private

Definition at line 88 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

bool SiStripBaselineAnalyzer::plotAPVCM_
private

Definition at line 94 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

bool SiStripBaselineAnalyzer::plotBaseline_
private

Definition at line 91 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

bool SiStripBaselineAnalyzer::plotBaselinePoints_
private

Definition at line 92 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

bool SiStripBaselineAnalyzer::plotClusters_
private

Definition at line 90 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

bool SiStripBaselineAnalyzer::plotPedestals_
private

Definition at line 95 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

bool SiStripBaselineAnalyzer::plotRawDigi_
private

Definition at line 93 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcAPVCM_
private

Definition at line 99 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcBaseline_
private

Definition at line 97 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcBaselinePoints_
private

Definition at line 98 of file SiStripBaselineAnalyzer.cc.

Referenced by SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcProcessedRawDigi_
private

Definition at line 100 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

std::auto_ptr<SiStripPedestalsSubtractor> SiStripBaselineAnalyzer::subtractorPed_
private

Definition at line 85 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

std::vector<TH1F> SiStripBaselineAnalyzer::vBaselineHisto_
private

Definition at line 114 of file SiStripBaselineAnalyzer.cc.

std::vector<TH1F> SiStripBaselineAnalyzer::vBaselinePointsHisto_
private

Definition at line 115 of file SiStripBaselineAnalyzer.cc.

std::vector<TH1F> SiStripBaselineAnalyzer::vClusterHisto_
private

Definition at line 116 of file SiStripBaselineAnalyzer.cc.

std::vector<TH1F> SiStripBaselineAnalyzer::vProcessedRawDigiHisto_
private

Definition at line 113 of file SiStripBaselineAnalyzer.cc.