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

Public Member Functions

 SiStripBaselineAnalyzer (const edm::ParameterSet &)
 
 ~SiStripBaselineAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

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

Private Attributes

uint16_t actualModule_
 
TCanvas * Canvas_
 
TH1F * h1BadAPVperEvent_
 
TH1F * h1Baseline_
 
TH1F * h1Clusters_
 
TH1F * h1ProcessedRawDigis_
 
uint16_t nModuletoDisplay_
 
TFile * oFile_
 
std::string outputFile_
 
edm::InputTag srcBaseline_
 
edm::InputTag srcProcessedRawDigi_
 
std::auto_ptr
< SiStripPedestalsSubtractor
subtractorPed_
 
std::vector< TH1F > vBaselineHisto_
 
std::vector< TH1F > vClusterHisto_
 
std::vector< TH1F > vProcessedRawDigiHisto_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

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

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

Definition at line 68 of file SiStripBaselineAnalyzer.cc.

Constructor & Destructor Documentation

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

Definition at line 103 of file SiStripBaselineAnalyzer.cc.

References SiStripRawProcessingFactory::create_SubtractorPed(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), h1BadAPVperEvent_, nModuletoDisplay_, outputFile_, srcBaseline_, srcProcessedRawDigi_, subtractorPed_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

103  {
104 
105  outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "StripHistos.root");
106  srcBaseline_ = conf.getParameter<edm::InputTag>( "srcBaseline" );
107  srcProcessedRawDigi_ = conf.getParameter<edm::InputTag>( "srcProcessedRawDigi" );
109  nModuletoDisplay_ = conf.getParameter<uint32_t>( "nModuletoDisplay_" );
110 
111  vProcessedRawDigiHisto_.clear();
112  vProcessedRawDigiHisto_.reserve(10000);
113 
114  vBaselineHisto_.clear();
115  vBaselineHisto_.reserve(10000);
116 
117  vClusterHisto_.clear();
118  vClusterHisto_.reserve(10000);
119 
120 
121  h1BadAPVperEvent_ = new TH1F("BadAPV/Event","BadAPV/Event", 2001, -0.5, 2000.5);
122  h1BadAPVperEvent_->SetXTitle("# Modules with Bad APVs");
123  h1BadAPVperEvent_->SetYTitle("Entries");
124  h1BadAPVperEvent_->SetLineWidth(2);
125  h1BadAPVperEvent_->SetLineStyle(2);
126 
127 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< TH1F > vBaselineHisto_
std::vector< TH1F > vClusterHisto_
std::vector< TH1F > vProcessedRawDigiHisto_
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed_
static std::auto_ptr< SiStripPedestalsSubtractor > create_SubtractorPed(const edm::ParameterSet &)
SiStripBaselineAnalyzer::~SiStripBaselineAnalyzer ( )

Definition at line 130 of file SiStripBaselineAnalyzer.cc.

131 {
132 
133 
134 
135 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 138 of file SiStripBaselineAnalyzer.cc.

References actualModule_, edm::DetSet< T >::begin(), edmNew::DetSetVector< T >::begin(), gather_cfg::cout, edmNew::DetSetVector< T >::end(), edm::EventID::event(), event(), edm::Event::getByLabel(), h1BadAPVperEvent_, h1Baseline_, h1Clusters_, h1ProcessedRawDigis_, edm::EventBase::id(), nModuletoDisplay_, edm::EventID::run(), DTTTrigCorrFirst::run, getRunRegistry::runs, srcBaseline_, srcProcessedRawDigi_, strip(), subtractorPed_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

139 {
140  using namespace edm;
141 
142 
143  subtractorPed_->init(es);
144 
146  e.getByLabel(srcBaseline_, moduleBaseline);
147 
149  e.getByLabel(srcProcessedRawDigi_,moduleRawDigi);
150 
152  edm::InputTag clusLabel("siStripClusters");
153  e.getByLabel(clusLabel, clusters);
154 
155 
156  char detIds[20];
157  char evs[20];
158  char runs[20];
159 
160  edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itDSBaseline = moduleBaseline->begin();
161  edm::DetSetVector<SiStripRawDigi>::const_iterator itRawDigis = moduleRawDigi->begin();
162 
163  uint32_t NBabAPVs = moduleRawDigi->size();
164  std::cout<< "Number of module with HIP in this event: " << NBabAPVs << std::endl;
165  h1BadAPVperEvent_->Fill(NBabAPVs);
166 
167  for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis, ++itDSBaseline) {
168  if(actualModule_ > nModuletoDisplay_) return;
169  uint32_t detId = itRawDigis->id;
170 
171 
172  if(itDSBaseline->id != detId){
173  std::cout << "Collections out of Synch. Something of fishy is going on ;-)" << std::endl;
174  return;
175  }
176 
177  actualModule_++;
178  uint32_t event = e.id().event();
179  uint32_t run = e.id().run();
180  //std::cout << "processing module N: " << actualModule_<< " detId: " << detId << " event: "<< event << std::endl;
181 
182  sprintf(detIds,"%ul", detId);
183  sprintf(evs,"%ul", event);
184  sprintf(runs,"%ul", run);
185  char* dHistoName = Form("Id:%s_run:%s_ev:%s",detIds, runs, evs);
186  h1ProcessedRawDigis_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5);
187  h1Baseline_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5);
188  h1Clusters_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5);
189 
190 
191  h1ProcessedRawDigis_->SetXTitle("strip#");
192  h1ProcessedRawDigis_->SetYTitle("ADC");
193  h1ProcessedRawDigis_->SetMaximum(1024.);
194  h1ProcessedRawDigis_->SetMinimum(-300.);
195  h1ProcessedRawDigis_->SetLineWidth(2);
196 
197 
198  h1Baseline_->SetXTitle("strip#");
199  h1Baseline_->SetYTitle("ADC");
200  h1Baseline_->SetMaximum(1024.);
201  h1Baseline_->SetMinimum(-300.);
202  h1Baseline_->SetLineWidth(2);
203  h1Baseline_->SetLineStyle(2);
204  h1Baseline_->SetLineColor(2);
205 
206  h1Clusters_->SetXTitle("strip#");
207  h1Clusters_->SetYTitle("ADC");
208  h1Clusters_->SetMaximum(1024.);
209  h1Clusters_->SetMinimum(-300.);
210  h1Clusters_->SetLineWidth(2);
211  h1Clusters_->SetLineStyle(2);
212  h1Clusters_->SetLineColor(3);
213 
215  std::vector<int16_t>::const_iterator itProcessedRawDigis;
216 
217 
218  std::vector<int16_t> ProcessedRawDigis(itRawDigis->size());
219  subtractorPed_->subtract( *itRawDigis, ProcessedRawDigis);
220 
221 
222  int strip =0;
223  for(itProcessedRawDigis = ProcessedRawDigis.begin(), itBaseline = itDSBaseline->begin();itProcessedRawDigis != ProcessedRawDigis.end(); ++itProcessedRawDigis, ++itBaseline){
224  h1ProcessedRawDigis_->Fill(strip, *itProcessedRawDigis);
225  h1Baseline_->Fill(strip, itBaseline->adc());
226  ++strip;
227  }
228 
229  edmNew::DetSetVector<SiStripCluster>::const_iterator itClusters = clusters->begin();
230  for ( ; itClusters != clusters->end(); ++itClusters ){
231  for ( edmNew::DetSet<SiStripCluster>::const_iterator clus = itClusters->begin(); clus != itClusters->end(); ++clus){
232  if(clus->geographicalId() == detId){
233  int firststrip = clus->firstStrip();
234  //std::cout << "Found cluster in detId " << detId << " " << firststrip << " " << clus->amplitudes().size() << " -----------------------------------------------" << std::endl;
235  strip=0;
236  for( std::vector<uint8_t>::const_iterator itAmpl = clus->amplitudes().begin(); itAmpl != clus->amplitudes().end(); ++itAmpl){
237  h1Clusters_->Fill(firststrip+strip, *itAmpl);
238  ++strip;
239  }
240  }
241  }
242  }
243 
244 
245 
247  vBaselineHisto_.push_back(*h1Baseline_);
248  vClusterHisto_.push_back(*h1Clusters_);
249 
250  }
251 
252 
253 }
RunNumber_t run() const
Definition: EventID.h:42
EventNumber_t event() const
Definition: EventID.h:44
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator begin() const
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
std::vector< TH1F > vBaselineHisto_
data_type const * const_iterator
Definition: DetSetNew.h:25
const_iterator end() const
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
std::vector< TH1F > vClusterHisto_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
iterator begin()
Definition: DetSet.h:58
edm::EventID id() const
Definition: EventBase.h:56
std::vector< TH1F > vProcessedRawDigiHisto_
tuple cout
Definition: gather_cfg.py:41
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed_
void SiStripBaselineAnalyzer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 257 of file SiStripBaselineAnalyzer.cc.

References actualModule_.

258 {
259 
260 
261  actualModule_ =0;
262 
263 
264 }
void SiStripBaselineAnalyzer::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 268 of file SiStripBaselineAnalyzer.cc.

References h1BadAPVperEvent_, oFile_, outputFile_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

268  {
269  oFile_ = new TFile((const char*)outputFile_.c_str(), "RECREATE");
270  oFile_->mkdir("ProcessedRawDigis");
271  oFile_->mkdir("Baseline");
272  oFile_->mkdir("Clusters");
273 
274 
275  std::vector<TH1F>::iterator itvProcessedRawDigis, itvBaseline, itvClusters;
276  itvProcessedRawDigis = vProcessedRawDigiHisto_.begin();
277  itvBaseline = vBaselineHisto_.begin();
278  itvClusters = vClusterHisto_.begin();
279 
280  oFile_->cd();
281  h1BadAPVperEvent_->Write();
282  for(; itvProcessedRawDigis != vProcessedRawDigiHisto_.end(); ++itvProcessedRawDigis, ++itvBaseline, ++itvClusters){
283  oFile_->cd("ProcessedRawDigis");
284  itvProcessedRawDigis->Write();
285  oFile_->cd("Baseline");
286  itvBaseline->Write();
287  oFile_->cd("Clusters");
288  itvClusters->Write();
289  }
290  oFile_->Write();
291  oFile_->Close();
292 
293 }
std::vector< TH1F > vBaselineHisto_
std::vector< TH1F > vClusterHisto_
std::vector< TH1F > vProcessedRawDigiHisto_

Member Data Documentation

uint16_t SiStripBaselineAnalyzer::actualModule_
private

Definition at line 99 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and beginJob().

TCanvas* SiStripBaselineAnalyzer::Canvas_
private

Definition at line 93 of file SiStripBaselineAnalyzer.cc.

TH1F* SiStripBaselineAnalyzer::h1BadAPVperEvent_
private

Definition at line 85 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

TH1F* SiStripBaselineAnalyzer::h1Baseline_
private

Definition at line 88 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

TH1F* SiStripBaselineAnalyzer::h1Clusters_
private

Definition at line 89 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

TH1F* SiStripBaselineAnalyzer::h1ProcessedRawDigis_
private

Definition at line 87 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

uint16_t SiStripBaselineAnalyzer::nModuletoDisplay_
private

Definition at line 98 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

TFile* SiStripBaselineAnalyzer::oFile_
private

Definition at line 91 of file SiStripBaselineAnalyzer.cc.

Referenced by endJob().

std::string SiStripBaselineAnalyzer::outputFile_
private

Definition at line 81 of file SiStripBaselineAnalyzer.cc.

Referenced by endJob(), and SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcBaseline_
private

Definition at line 82 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

edm::InputTag SiStripBaselineAnalyzer::srcProcessedRawDigi_
private

Definition at line 83 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

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

Definition at line 79 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

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

Definition at line 95 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

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

Definition at line 96 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

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

Definition at line 94 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().