CMS 3D CMS Logo

SiStripApvGainReader.cc
Go to the documentation of this file.
5 
7 
8 #include <iostream>
9 #include <cstdio>
10 #include <sys/time.h>
11 
12 
13 using namespace cms;
14 
16  printdebug_(iConfig.getUntrackedParameter<bool>("printDebug",true)),
17  formatedOutput_(iConfig.getUntrackedParameter<std::string>("outputFile","")),
18  gainType_ (iConfig.getUntrackedParameter<uint32_t>("gainType",1)),
19  tree_(0){
20  if (fs_.isAvailable()){
21  tree_=fs_->make<TTree>("Gains","Gains");
22 
23  tree_->Branch("Index",&id_,"Index/I");
24  tree_->Branch("DetId",&detId_,"DetId/I");
25  tree_->Branch("APVId",&apvId_,"APVId/I");
26  tree_->Branch("Gain",&gain_,"Gain/D");
27  }
28 }
29 
31 
33 
34  edm::ESHandle<SiStripGain> SiStripApvGain_;
35  iSetup.get<SiStripGainRcd>().get(SiStripApvGain_);
36  edm::LogInfo("SiStripApvGainReader") << "[SiStripApvGainReader::analyze] End Reading SiStripApvGain" << std::endl;
37  std::vector<uint32_t> detid;
38  SiStripApvGain_->getDetIds(detid);
39  edm::LogInfo("Number of detids ") << detid.size() << std::endl;
40 
41  FILE* pFile=nullptr;
42  if(formatedOutput_!="")pFile=fopen(formatedOutput_.c_str(), "w");
43  for (size_t id=0;id<detid.size();id++){
44  SiStripApvGain::Range range=SiStripApvGain_->getRange(detid[id], gainType_);
45  if(printdebug_){
46  int apv=0;
47  for(int it=0;it<range.second-range.first;it++){
48  edm::LogInfo("SiStripApvGainReader") << "detid " << detid[id] << " \t " << apv++ << " \t " << SiStripApvGain_->getApvGain(it,range) << std::endl;
49  id_++;
50 
51  if (tree_){
52  detId_ = detid[id];
53  apvId_ = apv ;
54  gain_ = SiStripApvGain_->getApvGain(it,range);
55  tree_->Fill();
56  }
57  }
58  }
59 
60  if(pFile){
61  fprintf(pFile,"%i ",detid[id]);
62  for(int it=0;it<range.second-range.first;it++){
63  fprintf(pFile,"%f ", SiStripApvGain_->getApvGain(it,range) );
64  }fprintf(pFile, "\n");
65  }
66 
67  }
68 
69  if(pFile)fclose(pFile);
70 }
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
static float getApvGain(const uint16_t &apv, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:73
bool isAvailable() const
Definition: Service.h:46
std::pair< ContainerIterator, ContainerIterator > Range
void getDetIds(std::vector< uint32_t > &DetIds_) const
ATTENTION: we assume the detIds are the same as those from the first gain.
Definition: SiStripGain.cc:108
edm::Service< TFileService > fs_
void analyze(const edm::Event &, const edm::EventSetup &) override
const T & get() const
Definition: EventSetup.h:59
SiStripApvGainReader(const edm::ParameterSet &)
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:70