CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/CalibTracker/SiStripChannelGain/plugins/SiStripGainCosmicCalculator.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 // Package:    SiStripChannelGain
00003 // Class:      SiStripGainCosmicCalculator
00004 // Original Author:  G. Bruno, D. Kcira
00005 //         Created:  Mon May 20 10:04:31 CET 2007
00006 // $Id: SiStripGainCosmicCalculator.cc,v 1.13 2013/01/11 05:51:19 wmtan Exp $
00007 #include "CalibTracker/SiStripChannelGain/plugins/SiStripGainCosmicCalculator.h"
00008 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
00009 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00014 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00015 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00016 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00017 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
00018 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00019 #include "CLHEP/Random/RandFlat.h"
00020 #include "CLHEP/Random/RandGauss.h"
00021 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00022 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00023 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetType.h"
00024 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00025 #include "DataFormats/TrackReco/interface/Track.h"
00026 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00027 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00028 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00029 #include "DataFormats/SiStripDetId/interface/SiStripSubStructure.h"
00030 #include "DataFormats/DetId/interface/DetId.h"
00031 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00032 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00033 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00034 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00035 //#include "DQM/SiStripCommon/interface/SiStripGenerateKey.h"
00036 
00037 //---------------------------------------------------------------------------------------------------------
00038 SiStripGainCosmicCalculator::SiStripGainCosmicCalculator(const edm::ParameterSet& iConfig) : ConditionDBWriter<SiStripApvGain>(iConfig){
00039   edm::LogInfo("SiStripGainCosmicCalculator::SiStripGainCosmicCalculator");
00040   ExpectedChargeDeposition = 200.;
00041   edm::LogInfo("SiStripApvGainCalculator::SiStripApvGainCalculator")<<"ExpectedChargeDeposition="<<ExpectedChargeDeposition;
00042 
00043   TrackProducer =  iConfig.getParameter<std::string>("TrackProducer");
00044   TrackLabel    =  iConfig.getParameter<std::string>("TrackLabel");
00045 
00046   detModulesToBeExcluded.clear(); detModulesToBeExcluded = iConfig.getParameter< std::vector<unsigned> >("detModulesToBeExcluded");
00047   MinNrEntries = iConfig.getUntrackedParameter<unsigned>("minNrEntries", 20);
00048   MaxChi2OverNDF = iConfig.getUntrackedParameter<double>("maxChi2OverNDF", 5.);
00049 
00050   outputHistogramsInRootFile = iConfig.getParameter<bool>("OutputHistogramsInRootFile");
00051   outputFileName = iConfig.getParameter<std::string>("OutputFileName");
00052 
00053   edm::LogInfo("SiStripApvGainCalculator")<<"Clusters from "<<detModulesToBeExcluded.size()<<" modules will be ignored in the calibration:";
00054   edm::LogInfo("SiStripApvGainCalculator")<<"The calibration for these DetIds will be set to a default value";
00055   for( std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin(); imod != detModulesToBeExcluded.end(); imod++){
00056     edm::LogInfo("SiStripApvGainCalculator")<<"exclude detid = "<< *imod;
00057   }
00058 
00059   printdebug_ = iConfig.getUntrackedParameter<bool>("printDebug", false);
00060   tTopo = nullptr;
00061 }
00062 
00063 
00064 SiStripGainCosmicCalculator::~SiStripGainCosmicCalculator(){
00065   edm::LogInfo("SiStripGainCosmicCalculator::~SiStripGainCosmicCalculator");
00066 }
00067 
00068 void SiStripGainCosmicCalculator::algoEndJob(){
00069 }
00070 
00071 void SiStripGainCosmicCalculator::algoBeginJob(const edm::EventSetup& iSetup)
00072 {
00073    //Retrieve tracker topology from geometry
00074    edm::ESHandle<TrackerTopology> tTopoHandle;
00075    iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
00076    tTopo = tTopoHandle.product();
00077  
00078    eventSetupCopy_ = &iSetup;
00079    std::cout<<"SiStripGainCosmicCalculator::algoBeginJob called"<<std::endl;
00080    total_nr_of_events = 0;
00081    HlistAPVPairs = new TObjArray(); HlistOtherHistos = new TObjArray();
00082    //
00083    HlistOtherHistos->Add(new TH1F( Form("APVPairCorrections"), Form("APVPairCorrections"), 50,-1.,4.));
00084    HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB1mono"),Form("APVPairCorrectionsTIB1mono"),50,-1.,4.));
00085    HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB1stereo"),Form("APVPairCorrectionsTIB1stereo"),50,-1.,4.));
00086    HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB2"),Form("APVPairCorrectionsTIB2"),50,-1.,4.));
00087    HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB1"),Form("APVPairCorrectionsTOB1"),50,-1.,4.));
00088    HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB2"),Form("APVPairCorrectionsTOB2"),50,-1.,4.));
00089    HlistOtherHistos->Add(new TH1F(Form("LocalAngle"),Form("LocalAngle"),70,-0.1,3.4));
00090    HlistOtherHistos->Add(new TH1F(Form("LocalAngleAbsoluteCosine"),Form("LocalAngleAbsoluteCosine"),48,-0.1,1.1));
00091    HlistOtherHistos->Add(new TH1F(Form("LocalPosition_cm"),Form("LocalPosition_cm"),100,-5.,5.));
00092    HlistOtherHistos->Add(new TH1F(Form("LocalPosition_normalized"),Form("LocalPosition_normalized"),100,-1.1,1.1));
00093    TH1F* local_histo = new TH1F(Form("SiStripRecHitType"),Form("SiStripRecHitType"),2,0.5,2.5); HlistOtherHistos->Add(local_histo);
00094    local_histo->GetXaxis()->SetBinLabel(1,"simple"); local_histo->GetXaxis()->SetBinLabel(2,"matched");
00095 
00096    // get cabling and find out list of active detectors
00097    edm::ESHandle<SiStripDetCabling> siStripDetCabling; iSetup.get<SiStripDetCablingRcd>().get(siStripDetCabling);
00098    std::vector<uint32_t> activeDets; activeDets.clear();
00099    SelectedDetIds.clear();
00100    siStripDetCabling->addActiveDetectorsRawIds(activeDets);
00101 //    SelectedDetIds = activeDets; // all active detector modules
00102    // use SiStripSubStructure for selecting certain regions
00103    SiStripSubStructure substructure;
00104    substructure.getTIBDetectors(activeDets, SelectedDetIds, 0, 0, 0, 0); // this adds rawDetIds to SelectedDetIds
00105    substructure.getTOBDetectors(activeDets, SelectedDetIds, 0, 0, 0);    // this adds rawDetIds to SelectedDetIds
00106    // get tracker geometry and find nr. of apv pairs for each active detector 
00107    edm::ESHandle<TrackerGeometry> tkGeom; iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );     
00108    for(TrackerGeometry::DetContainer::const_iterator it = tkGeom->dets().begin(); it != tkGeom->dets().end(); it++){ // loop over detector modules
00109      if( dynamic_cast<StripGeomDetUnit*>((*it))!=0){
00110        uint32_t detid=((*it)->geographicalId()).rawId();
00111        // get thickness for all detector modules, not just for active, this is strange 
00112        double module_thickness = (*it)->surface().bounds().thickness(); // get thickness of detector from GeomDet (DetContainer == vector<GeomDet*>)
00113        thickness_map.insert(std::make_pair(detid,module_thickness));
00114        //
00115        bool is_active_detector = false;
00116        for(std::vector<uint32_t>::iterator iactive = SelectedDetIds.begin(); iactive != SelectedDetIds.end(); iactive++){
00117          if( *iactive == detid ){
00118            is_active_detector = true;
00119            break; // leave for loop if found matching detid
00120          }
00121        }
00122        //
00123        bool exclude_this_detid = false;
00124        for( std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin(); imod != detModulesToBeExcluded.end(); imod++ ){
00125            if(*imod == detid) exclude_this_detid = true; // found in exclusion list
00126            break;
00127        }
00128        //
00129        if(is_active_detector && (!exclude_this_detid)){ // check whether is active detector and that should not be excluded
00130          const StripTopology& p = dynamic_cast<StripGeomDetUnit*>((*it))->specificTopology();
00131          unsigned short NAPVPairs = p.nstrips()/256;
00132          if( NAPVPairs<2 || NAPVPairs>3 ) {
00133            edm::LogError("SiStripGainCosmicCalculator")<<"Problem with Number of strips in detector: "<<p.nstrips()<<" Exiting program";
00134            exit(1);
00135          }
00136          for(int iapp = 0; iapp<NAPVPairs; iapp++){
00137            TString hid = Form("ChargeAPVPair_%i_%i",detid,iapp);
00138            HlistAPVPairs->Add(new TH1F(hid,hid,45,0.,1350.)); // multiply by 3 to take into account division by width
00139          }
00140        }
00141      }
00142    }
00143 }
00144 
00145 //---------------------------------------------------------------------------------------------------------
00146 void SiStripGainCosmicCalculator::algoAnalyze(const edm::Event & iEvent, const edm::EventSetup& iSetup){
00147   using namespace edm;
00148   total_nr_of_events++;
00149 
00150   //TO BE RESTORED
00151   //  anglefinder_->init(event,iSetup);
00152 
00153 
00154   // get seeds
00155 //  edm::Handle<TrajectorySeedCollection> seedcoll;
00156 //  event.getByType(seedcoll);
00157   // get tracks
00158   Handle<reco::TrackCollection> trackCollection; iEvent.getByLabel(TrackProducer, TrackLabel, trackCollection);
00159   const reco::TrackCollection *tracks=trackCollection.product();
00160 
00161 //  // get magnetic field
00162 //  edm::ESHandle<MagneticField> esmagfield;
00163 //  es.get<IdealMagneticFieldRecord>().get(esmagfield);
00164 //  magfield=&(*esmagfield);
00165   // loop over tracks
00166   for(reco::TrackCollection::const_iterator itr = tracks->begin(); itr != tracks->end(); itr++){ // looping over tracks
00167 
00168     //TO BE RESTORED
00169     //    std::vector<std::pair<const TrackingRecHit *,float> >hitangle =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
00170     std::vector<std::pair<const TrackingRecHit *,float> >hitangle;// =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
00171 
00172     for(std::vector<std::pair<const TrackingRecHit *,float> >::const_iterator hitangle_iter=hitangle.begin();hitangle_iter!=hitangle.end();hitangle_iter++){
00173       const TrackingRecHit * trechit = hitangle_iter->first;
00174       float local_angle=hitangle_iter->second;
00175       LocalPoint local_position= trechit->localPosition();
00176       const SiStripRecHit2D* sistripsimplehit=dynamic_cast<const SiStripRecHit2D*>(trechit);
00177       const SiStripMatchedRecHit2D* sistripmatchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>(trechit);
00178 //      std::cout<<" hit/matched "<<std::ios::hex<<sistripsimplehit<<" "<<sistripmatchedhit<<std::endl;
00179       ((TH1F*) HlistOtherHistos->FindObject("LocalAngle"))->Fill(local_angle);
00180       ((TH1F*) HlistOtherHistos->FindObject("LocalAngleAbsoluteCosine"))->Fill(fabs(cos(local_angle)));
00181       if(sistripsimplehit){
00182         ((TH1F*) HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(1.);
00183         const SiStripRecHit2D::ClusterRef & cluster=sistripsimplehit->cluster();
00184         const std::vector<uint8_t>& ampls = cluster->amplitudes();
00185 //        const std::vector<uint16_t>& ampls = cluster->amplitudes();
00186         uint32_t thedetid  = cluster->geographicalId();
00187         double module_width = moduleWidth(thedetid, &iSetup);
00188         ((TH1F*) HlistOtherHistos->FindObject("LocalPosition_cm"))->Fill(local_position.x());
00189         ((TH1F*) HlistOtherHistos->FindObject("LocalPosition_normalized"))->Fill(local_position.x()/module_width);
00190         double module_thickness = moduleThickness(thedetid, &iSetup);
00191         int ifirststrip= cluster->firstStrip();
00192         int theapvpairid = int(float(ifirststrip)/256.);
00193         TH1F* histopointer = (TH1F*) HlistAPVPairs->FindObject(Form("ChargeAPVPair_%i_%i",thedetid,theapvpairid));
00194         if( histopointer ){
00195           short cCharge = 0;
00196           for(unsigned int iampl = 0; iampl<ampls.size(); iampl++){
00197           cCharge += ampls[iampl];
00198           }
00199           double cluster_charge_over_path = ((double)cCharge) * fabs(cos(local_angle)) / ( 10. * module_thickness);
00200           histopointer->Fill(cluster_charge_over_path);
00201         }
00202       }else{
00203            if(sistripmatchedhit) ((TH1F*) HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(2.);
00204       }
00205     }
00206   }
00207 }
00208 
00209 
00210 //---------------------------------------------------------------------------------------------------------
00211 std::pair<double,double> SiStripGainCosmicCalculator::getPeakOfLandau( TH1F * inputHisto){ // automated fitting with finding of the appropriate nr. of ADCs
00212   // set some default dummy value and return if no entries
00213   double adcs = -0.5; double error = 0.; double nr_of_entries = inputHisto->GetEntries();
00214   if(nr_of_entries < MinNrEntries){
00215     return std::make_pair(adcs,error);
00216   }
00217 //
00218 //  // fit with initial setting of  parameter values
00219 //  double rms_of_histogram = inputHisto->GetRMS();
00220 //  TF1 *landaufit = new TF1("landaufit","landau",0.,450.);
00221 //  landaufit->SetParameters(nr_of_entries,mean_of_histogram,rms_of_histogram);
00222 //  inputHisto->Fit("landaufit","0Q+");
00223 //  delete landaufit;
00224 //
00225   // perform fit with standard landau
00226   inputHisto->Fit("landau","0Q");
00227   TF1 * fitfunction = (TF1*) inputHisto->GetListOfFunctions()->First();
00228   adcs = fitfunction->GetParameter("MPV");
00229   error = fitfunction->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
00230   double chi2 = fitfunction->GetChisquare();
00231   double ndf = fitfunction->GetNDF();
00232   double chi2overndf = chi2 / ndf;
00233   // in case things went wrong, try to refit in smaller range
00234   if(adcs< 2. || (error/adcs)>1.8 ){
00235      inputHisto->Fit("landau","0Q",0,0.,400.);
00236      TF1 * fitfunction2 = (TF1*) inputHisto->GetListOfFunctions()->First();
00237      std::cout<<"refitting landau for histogram "<<inputHisto->GetTitle()<<std::endl;
00238      std::cout<<"initial error/adcs ="<<error<<" / "<<adcs<<std::endl;
00239      std::cout<<"new     error/adcs ="<<fitfunction2->GetParError(1)<<" / "<<fitfunction2->GetParameter("MPV")<<std::endl;
00240      adcs = fitfunction2->GetParameter("MPV");
00241      error = fitfunction2->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
00242      chi2 = fitfunction2->GetChisquare();
00243      ndf = fitfunction2->GetNDF();
00244      chi2overndf = chi2 / ndf;
00245    }
00246    // if still wrong, give up
00247    if(adcs<2. || chi2overndf>MaxChi2OverNDF){
00248      adcs = -0.5; error = 0.;
00249    }
00250   return std::make_pair(adcs,error);
00251 }
00252 
00253 //---------------------------------------------------------------------------------------------------------
00254 double SiStripGainCosmicCalculator::moduleWidth(const uint32_t detid, const edm::EventSetup* iSetup) // get width of the module detid
00255 { //dk: copied from A. Giammanco and hacked,  module_width values : 10.49 12.03 6.144 7.14 9.3696
00256   edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );     
00257   double module_width=0.;
00258   const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
00259   if (dynamic_cast<const StripGeomDetUnit*>(it)==0 && dynamic_cast<const PixelGeomDetUnit*>(it)==0) {
00260     std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
00261   }else{
00262     module_width = it->surface().bounds().width();
00263   }
00264   return module_width;
00265 }
00266 
00267 //---------------------------------------------------------------------------------------------------------
00268 double SiStripGainCosmicCalculator::moduleThickness(const uint32_t detid, const edm::EventSetup* iSetup) // get thickness of the module detid
00269 { //dk: copied from A. Giammanco and hacked
00270   edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );
00271   double module_thickness=0.;
00272   const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
00273   if (dynamic_cast<const StripGeomDetUnit*>(it)==0 && dynamic_cast<const PixelGeomDetUnit*>(it)==0) {
00274     std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
00275   }else{
00276     module_thickness = it->surface().bounds().thickness();
00277   }
00278   return module_thickness;
00279 }
00280 
00281 //---------------------------------------------------------------------------------------------------------
00282 SiStripApvGain * SiStripGainCosmicCalculator::getNewObject() {
00283   std::cout<<"SiStripGainCosmicCalculator::getNewObject called"<<std::endl;
00284 
00285   std::cout<<"total_nr_of_events="<<total_nr_of_events<<std::endl;
00286   // book some more histograms
00287   TH1F *ChargeOfEachAPVPair = new TH1F("ChargeOfEachAPVPair","ChargeOfEachAPVPair",1,0,1); ChargeOfEachAPVPair->SetBit(TH1::kCanRebin);
00288   TH1F *EntriesApvPairs = new TH1F("EntriesApvPairs","EntriesApvPairs",1,0,1); EntriesApvPairs->SetBit(TH1::kCanRebin);
00289   TH1F * NrOfEntries = new TH1F("NrOfEntries","NrOfEntries",351,-0.5,350.5);// NrOfEntries->SetBit(TH1::kCanRebin);
00290   TH1F * ModuleThickness = new TH1F("ModuleThickness","ModuleThickness",2,0.5,2.5); HlistOtherHistos->Add(ModuleThickness);
00291   ModuleThickness->GetXaxis()->SetBinLabel(1,"320mu"); ModuleThickness->GetXaxis()->SetBinLabel(2,"500mu"); ModuleThickness->SetYTitle("Nr APVPairs");
00292   TH1F * ModuleWidth = new TH1F("ModuleWidth","ModuleWidth",5,0.5,5.5); HlistOtherHistos->Add(ModuleWidth);
00293   ModuleWidth->GetXaxis()->SetBinLabel(1,"6.144cm"); ModuleWidth->GetXaxis()->SetBinLabel(2,"7.14cm");
00294   ModuleWidth->GetXaxis()->SetBinLabel(3,"9.3696cm"); ModuleWidth->GetXaxis()->SetBinLabel(4,"10.49cm");
00295   ModuleWidth->GetXaxis()->SetBinLabel(5,"12.03cm");
00296   ModuleWidth->SetYTitle("Nr APVPairs");
00297   // loop over single histograms and extract peak value of charge
00298   HlistAPVPairs->Sort(); // sort alfabetically
00299   TIter hiterator(HlistAPVPairs);
00300   double MeanCharge = 0.;
00301   double NrOfApvPairs = 0.;
00302   TH1F *MyHisto = (TH1F*)hiterator();
00303   while( MyHisto ){
00304     TString histo_title = MyHisto->GetTitle();
00305     if(histo_title.Contains("ChargeAPVPair_")){
00306       std::pair<double,double> two_values = getPeakOfLandau(MyHisto);
00307       double local_nrofadcs = two_values.first;
00308       double local_sigma = two_values.second;
00309       ChargeOfEachAPVPair->Fill(histo_title, local_nrofadcs);
00310       int ichbin  = ChargeOfEachAPVPair->GetXaxis()->FindBin(histo_title.Data());
00311       ChargeOfEachAPVPair->SetBinError(ichbin,local_sigma);
00312       EntriesApvPairs->Fill(histo_title, MyHisto->GetEntries());
00313       NrOfEntries->Fill(MyHisto->GetEntries());
00314       if(local_nrofadcs > 0){ // if nr of adcs is negative, the fitting routine could not extract meaningfull numbers
00315        MeanCharge += local_nrofadcs;
00316        NrOfApvPairs += 1.; // count nr of apv pairs since do not know whether nr of bins of histogram is the same
00317       }
00318     }
00319     MyHisto = (TH1F*)hiterator();
00320   }
00321   ChargeOfEachAPVPair->LabelsDeflate("X"); EntriesApvPairs->LabelsDeflate("X"); // trim nr. of bins to match active labels
00322   HlistOtherHistos->Add(ChargeOfEachAPVPair);
00323   HlistOtherHistos->Add(EntriesApvPairs);
00324   HlistOtherHistos->Add(NrOfEntries); 
00325   MeanCharge = MeanCharge / NrOfApvPairs;
00326   // calculate correction
00327   TH1F* CorrectionOfEachAPVPair = (TH1F*) ChargeOfEachAPVPair->Clone("CorrectionOfEachAPVPair");
00328   TH1F *ChargeOfEachAPVPairControlView = new TH1F("ChargeOfEachAPVPairControlView","ChargeOfEachAPVPairControlView",1,0,1); ChargeOfEachAPVPairControlView->SetBit(TH1::kCanRebin);
00329 TH1F *CorrectionOfEachAPVPairControlView = new TH1F("CorrectionOfEachAPVPairControlView","CorrectionOfEachAPVPairControlView",1,0,1); CorrectionOfEachAPVPairControlView->SetBit(TH1::kCanRebin);
00330   std::ofstream APVPairTextOutput("apvpair_corrections.txt");
00331   APVPairTextOutput<<"# MeanCharge = "<<MeanCharge<<std::endl;
00332   APVPairTextOutput<<"# Nr. of APVPairs = "<<NrOfApvPairs<<std::endl;
00333   for(int ibin=1; ibin <= ChargeOfEachAPVPair->GetNbinsX(); ibin++){
00334      TString local_bin_label = ChargeOfEachAPVPair->GetXaxis()->GetBinLabel(ibin);
00335      double local_charge_over_path = ChargeOfEachAPVPair->GetBinContent(ibin);
00336      if(local_bin_label.Contains("ChargeAPVPair_") && local_charge_over_path > 0.0000001){ // calculate correction only for meaningful numbers
00337        uint32_t extracted_detid; std::istringstream read_label((local_bin_label(14,9)).Data()); read_label >> extracted_detid; 
00338        unsigned short extracted_apvpairid; std::istringstream read_apvpair((local_bin_label(24,1)).Data()); read_apvpair >> extracted_apvpairid; 
00339        double local_error_of_charge = ChargeOfEachAPVPair->GetBinError(ibin);
00340        double local_correction = -0.5;
00341        double local_error_correction = 0.;
00342        local_correction = MeanCharge / local_charge_over_path; // later use ExpectedChargeDeposition instead of MeanCharge
00343        local_error_correction = local_correction * local_error_of_charge / local_charge_over_path;
00344        if(local_error_correction>1.8){ // understand why error too large sometimes
00345          std::cout<<"too large error "<<local_error_correction<<" for histogram "<<local_bin_label<<std::endl;
00346        }
00347        double nr_of_entries = EntriesApvPairs->GetBinContent(ibin);
00348        APVPairTextOutput<<local_bin_label<<" "<<local_correction<<" "<<local_charge_over_path<<" "<<nr_of_entries<<std::endl;
00349        CorrectionOfEachAPVPair->SetBinContent(ibin, local_correction);
00350        CorrectionOfEachAPVPair->SetBinError(ibin, local_error_correction);
00351        ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrections"))->Fill(local_correction);
00352        DetId thedetId = DetId(extracted_detid);
00353        unsigned int generalized_layer = 0;
00354        // calculate generalized_layer:  31,32 = TIB1, 33 = TIB2, 33 = TIB3, 51 = TOB1, 52 = TOB2, 60 = TEC
00355        if(thedetId.subdetId()==StripSubdetector::TIB){
00356           
00357           generalized_layer = 10*thedetId.subdetId() + tTopo->tibLayer(thedetId.rawId()) + tTopo->tibStereo(thedetId.rawId());
00358           if(tTopo->tibLayer(thedetId.rawId())==2){
00359             generalized_layer++;
00360             if (tTopo->tibGlued(thedetId.rawId())) edm::LogError("ClusterMTCCFilter")<<"WRONGGGG"<<std::endl;
00361           }
00362         }else{
00363           generalized_layer = 10*thedetId.subdetId();
00364           if(thedetId.subdetId()==StripSubdetector::TOB){
00365             
00366             generalized_layer += tTopo->tobLayer(thedetId.rawId());
00367           }
00368         }
00369        if(generalized_layer==31){
00370          ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB1mono"))->Fill(local_correction);
00371        }
00372        if(generalized_layer==32){
00373          ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB1stereo"))->Fill(local_correction);
00374        }
00375        if(generalized_layer==33){
00376         ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTIB2"))->Fill(local_correction);
00377        }
00378        if(generalized_layer==51){
00379         ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTOB1"))->Fill(local_correction);
00380        }
00381        if(generalized_layer==52){
00382         ((TH1F*) HlistOtherHistos->FindObject("APVPairCorrectionsTOB2"))->Fill(local_correction);
00383        }
00384        // control view
00385        edm::ESHandle<SiStripDetCabling> siStripDetCabling; eventSetupCopy_->get<SiStripDetCablingRcd>().get(siStripDetCabling);
00386        const FedChannelConnection& fedchannelconnection = siStripDetCabling->getConnection( extracted_detid, extracted_apvpairid );
00387        std::ostringstream local_key;
00388        // in S. Mersi's analysis the APVPair id seems to be used instead of the lldChannel, hence use the same here
00389        local_key<<"fecCrate"<<fedchannelconnection.fecCrate()<<"_fecSlot"<<fedchannelconnection.fecSlot()<<"_fecRing"<<fedchannelconnection.fecRing()<<"_ccuAddr"<<fedchannelconnection.ccuAddr()<<"_ccuChan"<<fedchannelconnection.ccuChan()<<"_apvPair"<<extracted_apvpairid;
00390        TString control_key = local_key.str();
00391        ChargeOfEachAPVPairControlView->Fill(control_key,local_charge_over_path);
00392        int ibin1  = ChargeOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
00393        ChargeOfEachAPVPairControlView->SetBinError(ibin1,local_error_of_charge);
00394        CorrectionOfEachAPVPairControlView->Fill(control_key, local_correction);
00395        int ibin2  = CorrectionOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
00396        CorrectionOfEachAPVPairControlView->SetBinError(ibin2, local_error_correction);
00397        // thickness of each module
00398        double module_thickness = moduleThickness(extracted_detid, eventSetupCopy_);
00399        if( fabs(module_thickness - 0.032)<0.001 ) ModuleThickness->Fill(1);
00400        if( fabs(module_thickness - 0.05)<0.001 )  ModuleThickness->Fill(2);
00401        // width of each module     
00402        double module_width = moduleWidth(extracted_detid, eventSetupCopy_);
00403        if(fabs(module_width-6.144)<0.01) ModuleWidth->Fill(1);
00404        if(fabs(module_width-7.14)<0.01) ModuleWidth->Fill(2);
00405        if(fabs(module_width-9.3696)<0.01) ModuleWidth->Fill(3);
00406        if(fabs(module_width-10.49)<0.01) ModuleWidth->Fill(4);
00407        if(fabs(module_width-12.03)<0.01) ModuleWidth->Fill(5);
00408      }
00409   }
00410   HlistOtherHistos->Add(CorrectionOfEachAPVPair);
00411   ChargeOfEachAPVPairControlView->LabelsDeflate("X");
00412   CorrectionOfEachAPVPairControlView->LabelsDeflate("X");
00413   HlistOtherHistos->Add(ChargeOfEachAPVPairControlView);
00414   HlistOtherHistos->Add(CorrectionOfEachAPVPairControlView);
00415   // output histograms to file
00416 
00417 
00418   if(outputHistogramsInRootFile){
00419     TFile *outputfile = new TFile(outputFileName,"RECREATE");
00420     HlistAPVPairs->Write();
00421     HlistOtherHistos->Write();
00422     outputfile->Close();
00423   }
00424 
00425   SiStripApvGain * obj = new SiStripApvGain();
00426 
00427 //   for(std::map<uint32_t,OptoScanAnalysis*>::const_iterator it = analyses.begin(); it != analyses.end(); it++){
00428 //     //Generate Gain for det detid
00429 //     std::vector<float> theSiStripVector;
00430 //     for(unsigned short j=0; j<it->second; j++){
00431 //       float gain;
00432 
00433 //       //      if(sigmaGain_/meanGain_ < 0.00001) gain = meanGain_;
00434 //       //      else{
00435 //       gain = CLHEP::RandGauss::shoot(meanGain_, sigmaGain_);
00436 //       if(gain<=minimumPosValue_) gain=minimumPosValue_;
00437 //       //      }
00438 
00439 //       if (printdebug_)
00440 //      edm::LogInfo("SiStripGainCalculator") << "detid " << it->first << " \t"
00441 //                                            << " apv " << j << " \t"
00442 //                                            << gain    << " \t" 
00443 //                                            << std::endl;         
00444 //       theSiStripVector.push_back(gain);
00445 //     }
00446 //     SiStripApvGain::Range range(theSiStripVector.begin(),theSiStripVector.end());
00447 //     if ( ! obj->put(it->first,range) )
00448 //       edm::LogError("SiStripGainCalculator")<<"[SiStripGainCalculator::beginJob] detid already exists"<<std::endl;
00449 //   }
00450   
00451   return obj;
00452 }
00453