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
SiStripGainCosmicCalculator Class Reference

#include <SiStripGainCosmicCalculator.h>

Inheritance diagram for SiStripGainCosmicCalculator:
ConditionDBWriter< SiStripApvGain > edm::EDAnalyzer

Public Member Functions

 SiStripGainCosmicCalculator (const edm::ParameterSet &)
 
 ~SiStripGainCosmicCalculator ()
 
- Public Member Functions inherited from ConditionDBWriter< SiStripApvGain >
 ConditionDBWriter (const edm::ParameterSet &iConfig)
 
virtual ~ConditionDBWriter ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

void algoAnalyze (const edm::Event &, const edm::EventSetup &)
 
void algoBeginJob (const edm::EventSetup &)
 
void algoEndJob ()
 
SiStripApvGaingetNewObject ()
 
std::pair< double, double > getPeakOfLandau (TH1F *inputHisto)
 
double moduleThickness (const uint32_t detid, const edm::EventSetup *iSetup)
 
double moduleWidth (const uint32_t detid, const edm::EventSetup *iSetup)
 

Private Attributes

std::vector< uint32_t > detModulesToBeExcluded
 
const edm::EventSetupeventSetupCopy_
 
double ExpectedChargeDeposition
 
TObjArray * HlistAPVPairs
 
TObjArray * HlistOtherHistos
 
double MaxChi2OverNDF
 
unsigned int MinNrEntries
 
TString outputFileName
 
bool outputHistogramsInRootFile
 
bool printdebug_
 
std::vector< uint32_t > SelectedDetIds
 
std::map< uint32_t, double > thickness_map
 
uint32_t total_nr_of_events
 
std::string TrackLabel
 
std::string TrackProducer
 

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)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from ConditionDBWriter< SiStripApvGain >
void setDoStore (const bool doStore)
 When set to false the payload will not be written to the db. More...
 
void storeOnDbNow ()
 
cond::Time_t timeOfLastIOV ()
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 24 of file SiStripGainCosmicCalculator.h.

Constructor & Destructor Documentation

SiStripGainCosmicCalculator::SiStripGainCosmicCalculator ( const edm::ParameterSet iConfig)
explicit

Definition at line 40 of file SiStripGainCosmicCalculator.cc.

References detModulesToBeExcluded, ExpectedChargeDeposition, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), MaxChi2OverNDF, MinNrEntries, outputFileName, outputHistogramsInRootFile, printdebug_, and TrackLabel.

41  edm::LogInfo("SiStripGainCosmicCalculator::SiStripGainCosmicCalculator");
43  edm::LogInfo("SiStripApvGainCalculator::SiStripApvGainCalculator")<<"ExpectedChargeDeposition="<<ExpectedChargeDeposition;
44 
45  TrackProducer = iConfig.getParameter<std::string>("TrackProducer");
46  TrackLabel = iConfig.getParameter<std::string>("TrackLabel");
47 
48  detModulesToBeExcluded.clear(); detModulesToBeExcluded = iConfig.getParameter< std::vector<unsigned> >("detModulesToBeExcluded");
49  MinNrEntries = iConfig.getUntrackedParameter<unsigned>("minNrEntries", 20);
50  MaxChi2OverNDF = iConfig.getUntrackedParameter<double>("maxChi2OverNDF", 5.);
51 
52  outputHistogramsInRootFile = iConfig.getParameter<bool>("OutputHistogramsInRootFile");
53  outputFileName = iConfig.getParameter<std::string>("OutputFileName");
54 
55  edm::LogInfo("SiStripApvGainCalculator")<<"Clusters from "<<detModulesToBeExcluded.size()<<" modules will be ignored in the calibration:";
56  edm::LogInfo("SiStripApvGainCalculator")<<"The calibration for these DetIds will be set to a default value";
57  for( std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin(); imod != detModulesToBeExcluded.end(); imod++){
58  edm::LogInfo("SiStripApvGainCalculator")<<"exclude detid = "<< *imod;
59  }
60 
61  printdebug_ = iConfig.getUntrackedParameter<bool>("printDebug", false);
62 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< uint32_t > detModulesToBeExcluded
SiStripGainCosmicCalculator::~SiStripGainCosmicCalculator ( )

Definition at line 65 of file SiStripGainCosmicCalculator.cc.

65  {
66  edm::LogInfo("SiStripGainCosmicCalculator::~SiStripGainCosmicCalculator");
67 }

Member Function Documentation

void SiStripGainCosmicCalculator::algoAnalyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Reimplemented from ConditionDBWriter< SiStripApvGain >.

Definition at line 142 of file SiStripGainCosmicCalculator.cc.

References funct::cos(), HcalObjRepresent::Fill(), edm::Event::getByLabel(), HlistAPVPairs, HlistOtherHistos, TrackingRecHit::localPosition(), moduleThickness(), moduleWidth(), total_nr_of_events, TrackLabel, testEve_cfg::tracks, and PV3DBase< T, PVType, FrameType >::x().

142  {
143  using namespace edm;
145 
146  //TO BE RESTORED
147  // anglefinder_->init(event,iSetup);
148 
149 
150  // get seeds
151 // edm::Handle<TrajectorySeedCollection> seedcoll;
152 // event.getByType(seedcoll);
153  // get tracks
154  Handle<reco::TrackCollection> trackCollection; iEvent.getByLabel(TrackProducer, TrackLabel, trackCollection);
155  const reco::TrackCollection *tracks=trackCollection.product();
156 
157 // // get magnetic field
158 // edm::ESHandle<MagneticField> esmagfield;
159 // es.get<IdealMagneticFieldRecord>().get(esmagfield);
160 // magfield=&(*esmagfield);
161  // loop over tracks
162  for(reco::TrackCollection::const_iterator itr = tracks->begin(); itr != tracks->end(); itr++){ // looping over tracks
163 
164  //TO BE RESTORED
165  // std::vector<std::pair<const TrackingRecHit *,float> >hitangle =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
166  std::vector<std::pair<const TrackingRecHit *,float> >hitangle;// =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
167 
168  for(std::vector<std::pair<const TrackingRecHit *,float> >::const_iterator hitangle_iter=hitangle.begin();hitangle_iter!=hitangle.end();hitangle_iter++){
169  const TrackingRecHit * trechit = hitangle_iter->first;
170  float local_angle=hitangle_iter->second;
171  LocalPoint local_position= trechit->localPosition();
172  const SiStripRecHit2D* sistripsimplehit=dynamic_cast<const SiStripRecHit2D*>(trechit);
173  const SiStripMatchedRecHit2D* sistripmatchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>(trechit);
174 // std::cout<<" hit/matched "<<std::ios::hex<<sistripsimplehit<<" "<<sistripmatchedhit<<std::endl;
175  ((TH1F*) HlistOtherHistos->FindObject("LocalAngle"))->Fill(local_angle);
176  ((TH1F*) HlistOtherHistos->FindObject("LocalAngleAbsoluteCosine"))->Fill(fabs(cos(local_angle)));
177  if(sistripsimplehit){
178  ((TH1F*) HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(1.);
179  const SiStripRecHit2D::ClusterRef & cluster=sistripsimplehit->cluster();
180  const std::vector<uint8_t>& ampls = cluster->amplitudes();
181 // const std::vector<uint16_t>& ampls = cluster->amplitudes();
182  uint32_t thedetid = cluster->geographicalId();
183  double module_width = moduleWidth(thedetid, &iSetup);
184  ((TH1F*) HlistOtherHistos->FindObject("LocalPosition_cm"))->Fill(local_position.x());
185  ((TH1F*) HlistOtherHistos->FindObject("LocalPosition_normalized"))->Fill(local_position.x()/module_width);
186  double module_thickness = moduleThickness(thedetid, &iSetup);
187  int ifirststrip= cluster->firstStrip();
188  int theapvpairid = int(float(ifirststrip)/256.);
189  TH1F* histopointer = (TH1F*) HlistAPVPairs->FindObject(Form("ChargeAPVPair_%i_%i",thedetid,theapvpairid));
190  if( histopointer ){
191  short cCharge = 0;
192  for(unsigned int iampl = 0; iampl<ampls.size(); iampl++){
193  cCharge += ampls[iampl];
194  }
195  double cluster_charge_over_path = ((double)cCharge) * fabs(cos(local_angle)) / ( 10. * module_thickness);
196  histopointer->Fill(cluster_charge_over_path);
197  }
198  }else{
199  if(sistripmatchedhit) ((TH1F*) HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(2.);
200  }
201  }
202  }
203 }
double moduleThickness(const uint32_t detid, const edm::EventSetup *iSetup)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
double moduleWidth(const uint32_t detid, const edm::EventSetup *iSetup)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple tracks
Definition: testEve_cfg.py:39
T x() const
Definition: PV3DBase.h:61
virtual LocalPoint localPosition() const =0
void SiStripGainCosmicCalculator::algoBeginJob ( const edm::EventSetup iSetup)
privatevirtual

Reimplemented from ConditionDBWriter< SiStripApvGain >.

Definition at line 72 of file SiStripGainCosmicCalculator.cc.

References gather_cfg::cout, cond::rpcobgas::detid, detModulesToBeExcluded, eventSetupCopy_, cmsRelvalreport::exit, edm::EventSetup::get(), SiStripSubStructure::getTIBDetectors(), SiStripSubStructure::getTOBDetectors(), HlistAPVPairs, HlistOtherHistos, StripTopology::nstrips(), AlCaHLTBitMon_ParallelJobs::p, SelectedDetIds, thickness_map, and total_nr_of_events.

73 {
74  eventSetupCopy_ = &iSetup;
75  std::cout<<"SiStripGainCosmicCalculator::algoBeginJob called"<<std::endl;
77  HlistAPVPairs = new TObjArray(); HlistOtherHistos = new TObjArray();
78  //
79  HlistOtherHistos->Add(new TH1F( Form("APVPairCorrections"), Form("APVPairCorrections"), 50,-1.,4.));
80  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB1mono"),Form("APVPairCorrectionsTIB1mono"),50,-1.,4.));
81  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB1stereo"),Form("APVPairCorrectionsTIB1stereo"),50,-1.,4.));
82  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB2"),Form("APVPairCorrectionsTIB2"),50,-1.,4.));
83  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB1"),Form("APVPairCorrectionsTOB1"),50,-1.,4.));
84  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB2"),Form("APVPairCorrectionsTOB2"),50,-1.,4.));
85  HlistOtherHistos->Add(new TH1F(Form("LocalAngle"),Form("LocalAngle"),70,-0.1,3.4));
86  HlistOtherHistos->Add(new TH1F(Form("LocalAngleAbsoluteCosine"),Form("LocalAngleAbsoluteCosine"),48,-0.1,1.1));
87  HlistOtherHistos->Add(new TH1F(Form("LocalPosition_cm"),Form("LocalPosition_cm"),100,-5.,5.));
88  HlistOtherHistos->Add(new TH1F(Form("LocalPosition_normalized"),Form("LocalPosition_normalized"),100,-1.1,1.1));
89  TH1F* local_histo = new TH1F(Form("SiStripRecHitType"),Form("SiStripRecHitType"),2,0.5,2.5); HlistOtherHistos->Add(local_histo);
90  local_histo->GetXaxis()->SetBinLabel(1,"simple"); local_histo->GetXaxis()->SetBinLabel(2,"matched");
91 
92  // get cabling and find out list of active detectors
93  edm::ESHandle<SiStripDetCabling> siStripDetCabling; iSetup.get<SiStripDetCablingRcd>().get(siStripDetCabling);
94  std::vector<uint32_t> activeDets; activeDets.clear();
95  SelectedDetIds.clear();
96  siStripDetCabling->addActiveDetectorsRawIds(activeDets);
97 // SelectedDetIds = activeDets; // all active detector modules
98  // use SiStripSubStructure for selecting certain regions
99  SiStripSubStructure substructure;
100  substructure.getTIBDetectors(activeDets, SelectedDetIds, 0, 0, 0, 0); // this adds rawDetIds to SelectedDetIds
101  substructure.getTOBDetectors(activeDets, SelectedDetIds, 0, 0, 0); // this adds rawDetIds to SelectedDetIds
102  // get tracker geometry and find nr. of apv pairs for each active detector
103  edm::ESHandle<TrackerGeometry> tkGeom; iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
104  for(TrackerGeometry::DetContainer::const_iterator it = tkGeom->dets().begin(); it != tkGeom->dets().end(); it++){ // loop over detector modules
105  if( dynamic_cast<StripGeomDetUnit*>((*it))!=0){
106  uint32_t detid=((*it)->geographicalId()).rawId();
107  // get thickness for all detector modules, not just for active, this is strange
108  double module_thickness = (*it)->surface().bounds().thickness(); // get thickness of detector from GeomDet (DetContainer == vector<GeomDet*>)
109  thickness_map.insert(std::make_pair(detid,module_thickness));
110  //
111  bool is_active_detector = false;
112  for(std::vector<uint32_t>::iterator iactive = SelectedDetIds.begin(); iactive != SelectedDetIds.end(); iactive++){
113  if( *iactive == detid ){
114  is_active_detector = true;
115  break; // leave for loop if found matching detid
116  }
117  }
118  //
119  bool exclude_this_detid = false;
120  for( std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin(); imod != detModulesToBeExcluded.end(); imod++ ){
121  if(*imod == detid) exclude_this_detid = true; // found in exclusion list
122  break;
123  }
124  //
125  if(is_active_detector && (!exclude_this_detid)){ // check whether is active detector and that should not be excluded
126  const StripTopology& p = dynamic_cast<StripGeomDetUnit*>((*it))->specificTopology();
127  unsigned short NAPVPairs = p.nstrips()/256;
128  if( NAPVPairs<2 || NAPVPairs>3 ) {
129  edm::LogError("SiStripGainCosmicCalculator")<<"Problem with Number of strips in detector: "<<p.nstrips()<<" Exiting program";
130  exit(1);
131  }
132  for(int iapp = 0; iapp<NAPVPairs; iapp++){
133  TString hid = Form("ChargeAPVPair_%i_%i",detid,iapp);
134  HlistAPVPairs->Add(new TH1F(hid,hid,45,0.,1350.)); // multiply by 3 to take into account division by width
135  }
136  }
137  }
138  }
139 }
virtual int nstrips() const =0
std::vector< uint32_t > SelectedDetIds
const edm::EventSetup * eventSetupCopy_
void getTOBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tobDetRawIds, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t rod=0) const
std::map< uint32_t, double > thickness_map
const T & get() const
Definition: EventSetup.h:55
std::vector< uint32_t > detModulesToBeExcluded
void getTIBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tibDetRawIds, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t int_ext=0, uint32_t string=0) const
tuple cout
Definition: gather_cfg.py:121
void SiStripGainCosmicCalculator::algoEndJob ( )
privatevirtual

Reimplemented from ConditionDBWriter< SiStripApvGain >.

Definition at line 69 of file SiStripGainCosmicCalculator.cc.

69  {
70 }
SiStripApvGain * SiStripGainCosmicCalculator::getNewObject ( )
privatevirtual

Implements ConditionDBWriter< SiStripApvGain >.

Definition at line 278 of file SiStripGainCosmicCalculator.cc.

References FedChannelConnection::ccuAddr(), FedChannelConnection::ccuChan(), gather_cfg::cout, eventSetupCopy_, FedChannelConnection::fecCrate(), FedChannelConnection::fecRing(), FedChannelConnection::fecSlot(), HcalObjRepresent::Fill(), edm::EventSetup::get(), getPeakOfLandau(), SiStripDetId::glued(), HlistAPVPairs, HlistOtherHistos, TOBDetId::layer(), TIBDetId::layer(), moduleThickness(), moduleWidth(), getGTfromDQMFile::obj, estimatePileup_makeJSON::outputfile, outputFileName, outputHistogramsInRootFile, DetId::rawId(), SiStripDetId::stereo(), DetId::subdetId(), StripSubdetector::TIB, StripSubdetector::TOB, and total_nr_of_events.

278  {
279  std::cout<<"SiStripGainCosmicCalculator::getNewObject called"<<std::endl;
280 
281  std::cout<<"total_nr_of_events="<<total_nr_of_events<<std::endl;
282  // book some more histograms
283  TH1F *ChargeOfEachAPVPair = new TH1F("ChargeOfEachAPVPair","ChargeOfEachAPVPair",1,0,1); ChargeOfEachAPVPair->SetBit(TH1::kCanRebin);
284  TH1F *EntriesApvPairs = new TH1F("EntriesApvPairs","EntriesApvPairs",1,0,1); EntriesApvPairs->SetBit(TH1::kCanRebin);
285  TH1F * NrOfEntries = new TH1F("NrOfEntries","NrOfEntries",351,-0.5,350.5);// NrOfEntries->SetBit(TH1::kCanRebin);
286  TH1F * ModuleThickness = new TH1F("ModuleThickness","ModuleThickness",2,0.5,2.5); HlistOtherHistos->Add(ModuleThickness);
287  ModuleThickness->GetXaxis()->SetBinLabel(1,"320mu"); ModuleThickness->GetXaxis()->SetBinLabel(2,"500mu"); ModuleThickness->SetYTitle("Nr APVPairs");
288  TH1F * ModuleWidth = new TH1F("ModuleWidth","ModuleWidth",5,0.5,5.5); HlistOtherHistos->Add(ModuleWidth);
289  ModuleWidth->GetXaxis()->SetBinLabel(1,"6.144cm"); ModuleWidth->GetXaxis()->SetBinLabel(2,"7.14cm");
290  ModuleWidth->GetXaxis()->SetBinLabel(3,"9.3696cm"); ModuleWidth->GetXaxis()->SetBinLabel(4,"10.49cm");
291  ModuleWidth->GetXaxis()->SetBinLabel(5,"12.03cm");
292  ModuleWidth->SetYTitle("Nr APVPairs");
293  // loop over single histograms and extract peak value of charge
294  HlistAPVPairs->Sort(); // sort alfabetically
295  TIter hiterator(HlistAPVPairs);
296  double MeanCharge = 0.;
297  double NrOfApvPairs = 0.;
298  TH1F *MyHisto = (TH1F*)hiterator();
299  while( MyHisto ){
300  TString histo_title = MyHisto->GetTitle();
301  if(histo_title.Contains("ChargeAPVPair_")){
302  std::pair<double,double> two_values = getPeakOfLandau(MyHisto);
303  double local_nrofadcs = two_values.first;
304  double local_sigma = two_values.second;
305  ChargeOfEachAPVPair->Fill(histo_title, local_nrofadcs);
306  int ichbin = ChargeOfEachAPVPair->GetXaxis()->FindBin(histo_title.Data());
307  ChargeOfEachAPVPair->SetBinError(ichbin,local_sigma);
308  EntriesApvPairs->Fill(histo_title, MyHisto->GetEntries());
309  NrOfEntries->Fill(MyHisto->GetEntries());
310  if(local_nrofadcs > 0){ // if nr of adcs is negative, the fitting routine could not extract meaningfull numbers
311  MeanCharge += local_nrofadcs;
312  NrOfApvPairs += 1.; // count nr of apv pairs since do not know whether nr of bins of histogram is the same
313  }
314  }
315  MyHisto = (TH1F*)hiterator();
316  }
317  ChargeOfEachAPVPair->LabelsDeflate("X"); EntriesApvPairs->LabelsDeflate("X"); // trim nr. of bins to match active labels
318  HlistOtherHistos->Add(ChargeOfEachAPVPair);
319  HlistOtherHistos->Add(EntriesApvPairs);
320  HlistOtherHistos->Add(NrOfEntries);
321  MeanCharge = MeanCharge / NrOfApvPairs;
322  // calculate correction
323  TH1F* CorrectionOfEachAPVPair = (TH1F*) ChargeOfEachAPVPair->Clone("CorrectionOfEachAPVPair");
324  TH1F *ChargeOfEachAPVPairControlView = new TH1F("ChargeOfEachAPVPairControlView","ChargeOfEachAPVPairControlView",1,0,1); ChargeOfEachAPVPairControlView->SetBit(TH1::kCanRebin);
325 TH1F *CorrectionOfEachAPVPairControlView = new TH1F("CorrectionOfEachAPVPairControlView","CorrectionOfEachAPVPairControlView",1,0,1); CorrectionOfEachAPVPairControlView->SetBit(TH1::kCanRebin);
326  std::ofstream APVPairTextOutput("apvpair_corrections.txt");
327  APVPairTextOutput<<"# MeanCharge = "<<MeanCharge<<std::endl;
328  APVPairTextOutput<<"# Nr. of APVPairs = "<<NrOfApvPairs<<std::endl;
329  for(int ibin=1; ibin <= ChargeOfEachAPVPair->GetNbinsX(); ibin++){
330  TString local_bin_label = ChargeOfEachAPVPair->GetXaxis()->GetBinLabel(ibin);
331  double local_charge_over_path = ChargeOfEachAPVPair->GetBinContent(ibin);
332  if(local_bin_label.Contains("ChargeAPVPair_") && local_charge_over_path > 0.0000001){ // calculate correction only for meaningful numbers
333  uint32_t extracted_detid; std::istringstream read_label((local_bin_label(14,9)).Data()); read_label >> extracted_detid;
334  unsigned short extracted_apvpairid; std::istringstream read_apvpair((local_bin_label(24,1)).Data()); read_apvpair >> extracted_apvpairid;
335  double local_error_of_charge = ChargeOfEachAPVPair->GetBinError(ibin);
336  double local_correction = -0.5;
337  double local_error_correction = 0.;
338  local_correction = MeanCharge / local_charge_over_path; // later use ExpectedChargeDeposition instead of MeanCharge
339  local_error_correction = local_correction * local_error_of_charge / local_charge_over_path;
340  if(local_error_correction>1.8){ // understand why error too large sometimes
341  std::cout<<"too large error "<<local_error_correction<<" for histogram "<<local_bin_label<<std::endl;
342  }
343  double nr_of_entries = EntriesApvPairs->GetBinContent(ibin);
344  APVPairTextOutput<<local_bin_label<<" "<<local_correction<<" "<<local_charge_over_path<<" "<<nr_of_entries<<std::endl;
345  CorrectionOfEachAPVPair->SetBinContent(ibin, local_correction);
346  CorrectionOfEachAPVPair->SetBinError(ibin, local_error_correction);
347  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrections"))->Fill(local_correction);
348  DetId thedetId = DetId(extracted_detid);
349  unsigned int generalized_layer = 0;
350  // calculate generalized_layer: 31,32 = TIB1, 33 = TIB2, 33 = TIB3, 51 = TOB1, 52 = TOB2, 60 = TEC
351  if(thedetId.subdetId()==StripSubdetector::TIB){
352  TIBDetId ptib = TIBDetId(thedetId.rawId());
353  generalized_layer = 10*thedetId.subdetId() + ptib.layer() + ptib.stereo();
354  if(ptib.layer()==2){
355  generalized_layer++;
356  if (ptib.glued()) edm::LogError("ClusterMTCCFilter")<<"WRONGGGG"<<std::endl;
357  }
358  }else{
359  generalized_layer = 10*thedetId.subdetId();
360  if(thedetId.subdetId()==StripSubdetector::TOB){
361  TOBDetId ptob = TOBDetId(thedetId.rawId());
362  generalized_layer += ptob.layer();
363  }
364  }
365  if(generalized_layer==31){
366  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB1mono"))->Fill(local_correction);
367  }
368  if(generalized_layer==32){
369  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB1stereo"))->Fill(local_correction);
370  }
371  if(generalized_layer==33){
372  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB2"))->Fill(local_correction);
373  }
374  if(generalized_layer==51){
375  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTOB1"))->Fill(local_correction);
376  }
377  if(generalized_layer==52){
378  ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTOB2"))->Fill(local_correction);
379  }
380  // control view
381  edm::ESHandle<SiStripDetCabling> siStripDetCabling; eventSetupCopy_->get<SiStripDetCablingRcd>().get(siStripDetCabling);
382  const FedChannelConnection& fedchannelconnection = siStripDetCabling->getConnection( extracted_detid, extracted_apvpairid );
383  std::ostringstream local_key;
384  // in S. Mersi's analysis the APVPair id seems to be used instead of the lldChannel, hence use the same here
385  local_key<<"fecCrate"<<fedchannelconnection.fecCrate()<<"_fecSlot"<<fedchannelconnection.fecSlot()<<"_fecRing"<<fedchannelconnection.fecRing()<<"_ccuAddr"<<fedchannelconnection.ccuAddr()<<"_ccuChan"<<fedchannelconnection.ccuChan()<<"_apvPair"<<extracted_apvpairid;
386  TString control_key = local_key.str();
387  ChargeOfEachAPVPairControlView->Fill(control_key,local_charge_over_path);
388  int ibin1 = ChargeOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
389  ChargeOfEachAPVPairControlView->SetBinError(ibin1,local_error_of_charge);
390  CorrectionOfEachAPVPairControlView->Fill(control_key, local_correction);
391  int ibin2 = CorrectionOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
392  CorrectionOfEachAPVPairControlView->SetBinError(ibin2, local_error_correction);
393  // thickness of each module
394  double module_thickness = moduleThickness(extracted_detid, eventSetupCopy_);
395  if( fabs(module_thickness - 0.032)<0.001 ) ModuleThickness->Fill(1);
396  if( fabs(module_thickness - 0.05)<0.001 ) ModuleThickness->Fill(2);
397  // width of each module
398  double module_width = moduleWidth(extracted_detid, eventSetupCopy_);
399  if(fabs(module_width-6.144)<0.01) ModuleWidth->Fill(1);
400  if(fabs(module_width-7.14)<0.01) ModuleWidth->Fill(2);
401  if(fabs(module_width-9.3696)<0.01) ModuleWidth->Fill(3);
402  if(fabs(module_width-10.49)<0.01) ModuleWidth->Fill(4);
403  if(fabs(module_width-12.03)<0.01) ModuleWidth->Fill(5);
404  }
405  }
406  HlistOtherHistos->Add(CorrectionOfEachAPVPair);
407  ChargeOfEachAPVPairControlView->LabelsDeflate("X");
408  CorrectionOfEachAPVPairControlView->LabelsDeflate("X");
409  HlistOtherHistos->Add(ChargeOfEachAPVPairControlView);
410  HlistOtherHistos->Add(CorrectionOfEachAPVPairControlView);
411  // output histograms to file
412 
413 
415  TFile *outputfile = new TFile(outputFileName,"RECREATE");
416  HlistAPVPairs->Write();
417  HlistOtherHistos->Write();
418  outputfile->Close();
419  }
420 
422 
423 // for(std::map<uint32_t,OptoScanAnalysis*>::const_iterator it = analyses.begin(); it != analyses.end(); it++){
424 // //Generate Gain for det detid
425 // std::vector<float> theSiStripVector;
426 // for(unsigned short j=0; j<it->second; j++){
427 // float gain;
428 
429 // // if(sigmaGain_/meanGain_ < 0.00001) gain = meanGain_;
430 // // else{
431 // gain = CLHEP::RandGauss::shoot(meanGain_, sigmaGain_);
432 // if(gain<=minimumPosValue_) gain=minimumPosValue_;
433 // // }
434 
435 // if (printdebug_)
436 // edm::LogInfo("SiStripGainCalculator") << "detid " << it->first << " \t"
437 // << " apv " << j << " \t"
438 // << gain << " \t"
439 // << std::endl;
440 // theSiStripVector.push_back(gain);
441 // }
442 // SiStripApvGain::Range range(theSiStripVector.begin(),theSiStripVector.end());
443 // if ( ! obj->put(it->first,range) )
444 // edm::LogError("SiStripGainCalculator")<<"[SiStripGainCalculator::beginJob] detid already exists"<<std::endl;
445 // }
446 
447  return obj;
448 }
const uint16_t & fecSlot() const
const uint16_t & fecCrate() const
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
uint32_t stereo() const
Definition: SiStripDetId.h:162
double moduleThickness(const uint32_t detid, const edm::EventSetup *iSetup)
const edm::EventSetup * eventSetupCopy_
double moduleWidth(const uint32_t detid, const edm::EventSetup *iSetup)
std::pair< double, double > getPeakOfLandau(TH1F *inputHisto)
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
const uint16_t & fecRing() const
Class containning control, module, detector and connection information, at the level of a FED channel...
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
const uint16_t & ccuChan() const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
const uint16_t & ccuAddr() const
Definition: DetId.h:20
uint32_t glued() const
Definition: SiStripDetId.h:154
const T & get() const
Definition: EventSetup.h:55
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
tuple cout
Definition: gather_cfg.py:121
std::pair< double, double > SiStripGainCosmicCalculator::getPeakOfLandau ( TH1F *  inputHisto)
private

Definition at line 207 of file SiStripGainCosmicCalculator.cc.

References gather_cfg::cout, error, MaxChi2OverNDF, and MinNrEntries.

Referenced by getNewObject().

207  { // automated fitting with finding of the appropriate nr. of ADCs
208  // set some default dummy value and return if no entries
209  double adcs = -0.5; double error = 0.; double nr_of_entries = inputHisto->GetEntries();
210  if(nr_of_entries < MinNrEntries){
211  return std::make_pair(adcs,error);
212  }
213 //
214 // // fit with initial setting of parameter values
215 // double rms_of_histogram = inputHisto->GetRMS();
216 // TF1 *landaufit = new TF1("landaufit","landau",0.,450.);
217 // landaufit->SetParameters(nr_of_entries,mean_of_histogram,rms_of_histogram);
218 // inputHisto->Fit("landaufit","0Q+");
219 // delete landaufit;
220 //
221  // perform fit with standard landau
222  inputHisto->Fit("landau","0Q");
223  TF1 * fitfunction = (TF1*) inputHisto->GetListOfFunctions()->First();
224  adcs = fitfunction->GetParameter("MPV");
225  error = fitfunction->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
226  double chi2 = fitfunction->GetChisquare();
227  double ndf = fitfunction->GetNDF();
228  double chi2overndf = chi2 / ndf;
229  // in case things went wrong, try to refit in smaller range
230  if(adcs< 2. || (error/adcs)>1.8 ){
231  inputHisto->Fit("landau","0Q",0,0.,400.);
232  TF1 * fitfunction2 = (TF1*) inputHisto->GetListOfFunctions()->First();
233  std::cout<<"refitting landau for histogram "<<inputHisto->GetTitle()<<std::endl;
234  std::cout<<"initial error/adcs ="<<error<<" / "<<adcs<<std::endl;
235  std::cout<<"new error/adcs ="<<fitfunction2->GetParError(1)<<" / "<<fitfunction2->GetParameter("MPV")<<std::endl;
236  adcs = fitfunction2->GetParameter("MPV");
237  error = fitfunction2->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
238  chi2 = fitfunction2->GetChisquare();
239  ndf = fitfunction2->GetNDF();
240  chi2overndf = chi2 / ndf;
241  }
242  // if still wrong, give up
243  if(adcs<2. || chi2overndf>MaxChi2OverNDF){
244  adcs = -0.5; error = 0.;
245  }
246  return std::make_pair(adcs,error);
247 }
tuple cout
Definition: gather_cfg.py:121
double SiStripGainCosmicCalculator::moduleThickness ( const uint32_t  detid,
const edm::EventSetup iSetup 
)
private

Definition at line 264 of file SiStripGainCosmicCalculator.cc.

References BoundSurface::bounds(), gather_cfg::cout, edm::EventSetup::get(), GeomDet::surface(), and Bounds::thickness().

Referenced by algoAnalyze(), and getNewObject().

265 { //dk: copied from A. Giammanco and hacked
266  edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );
267  double module_thickness=0.;
268  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
269  if (dynamic_cast<const StripGeomDetUnit*>(it)==0 && dynamic_cast<const PixelGeomDetUnit*>(it)==0) {
270  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
271  }else{
272  module_thickness = it->surface().bounds().thickness();
273  }
274  return module_thickness;
275 }
virtual float thickness() const =0
Definition: DetId.h:20
const Bounds & bounds() const
Definition: BoundSurface.h:89
const T & get() const
Definition: EventSetup.h:55
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
tuple cout
Definition: gather_cfg.py:121
double SiStripGainCosmicCalculator::moduleWidth ( const uint32_t  detid,
const edm::EventSetup iSetup 
)
private

Definition at line 250 of file SiStripGainCosmicCalculator.cc.

References BoundSurface::bounds(), gather_cfg::cout, edm::EventSetup::get(), GeomDet::surface(), and Bounds::width().

Referenced by algoAnalyze(), and getNewObject().

251 { //dk: copied from A. Giammanco and hacked, module_width values : 10.49 12.03 6.144 7.14 9.3696
252  edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );
253  double module_width=0.;
254  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
255  if (dynamic_cast<const StripGeomDetUnit*>(it)==0 && dynamic_cast<const PixelGeomDetUnit*>(it)==0) {
256  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
257  }else{
258  module_width = it->surface().bounds().width();
259  }
260  return module_width;
261 }
Definition: DetId.h:20
const Bounds & bounds() const
Definition: BoundSurface.h:89
const T & get() const
Definition: EventSetup.h:55
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
tuple cout
Definition: gather_cfg.py:121
virtual float width() const =0

Member Data Documentation

std::vector<uint32_t> SiStripGainCosmicCalculator::detModulesToBeExcluded
private

Definition at line 47 of file SiStripGainCosmicCalculator.h.

Referenced by algoBeginJob(), and SiStripGainCosmicCalculator().

const edm::EventSetup* SiStripGainCosmicCalculator::eventSetupCopy_
private

Definition at line 48 of file SiStripGainCosmicCalculator.h.

Referenced by algoBeginJob(), and getNewObject().

double SiStripGainCosmicCalculator::ExpectedChargeDeposition
private

Definition at line 44 of file SiStripGainCosmicCalculator.h.

Referenced by SiStripGainCosmicCalculator().

TObjArray* SiStripGainCosmicCalculator::HlistAPVPairs
private

Definition at line 41 of file SiStripGainCosmicCalculator.h.

Referenced by algoAnalyze(), algoBeginJob(), and getNewObject().

TObjArray* SiStripGainCosmicCalculator::HlistOtherHistos
private

Definition at line 42 of file SiStripGainCosmicCalculator.h.

Referenced by algoAnalyze(), algoBeginJob(), and getNewObject().

double SiStripGainCosmicCalculator::MaxChi2OverNDF
private

Definition at line 50 of file SiStripGainCosmicCalculator.h.

Referenced by getPeakOfLandau(), and SiStripGainCosmicCalculator().

unsigned int SiStripGainCosmicCalculator::MinNrEntries
private

Definition at line 49 of file SiStripGainCosmicCalculator.h.

Referenced by getPeakOfLandau(), and SiStripGainCosmicCalculator().

TString SiStripGainCosmicCalculator::outputFileName
private

Definition at line 52 of file SiStripGainCosmicCalculator.h.

Referenced by getNewObject(), and SiStripGainCosmicCalculator().

bool SiStripGainCosmicCalculator::outputHistogramsInRootFile
private

Definition at line 51 of file SiStripGainCosmicCalculator.h.

Referenced by getNewObject(), and SiStripGainCosmicCalculator().

bool SiStripGainCosmicCalculator::printdebug_
private

Definition at line 53 of file SiStripGainCosmicCalculator.h.

Referenced by SiStripGainCosmicCalculator().

std::vector<uint32_t> SiStripGainCosmicCalculator::SelectedDetIds
private

Definition at line 46 of file SiStripGainCosmicCalculator.h.

Referenced by algoBeginJob().

std::map<uint32_t, double> SiStripGainCosmicCalculator::thickness_map
private

Definition at line 45 of file SiStripGainCosmicCalculator.h.

Referenced by algoBeginJob().

uint32_t SiStripGainCosmicCalculator::total_nr_of_events
private

Definition at line 43 of file SiStripGainCosmicCalculator.h.

Referenced by algoAnalyze(), algoBeginJob(), and getNewObject().

std::string SiStripGainCosmicCalculator::TrackLabel
private

Definition at line 39 of file SiStripGainCosmicCalculator.h.

Referenced by algoAnalyze(), and SiStripGainCosmicCalculator().

std::string SiStripGainCosmicCalculator::TrackProducer
private

Definition at line 38 of file SiStripGainCosmicCalculator.h.