Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:46:42

0001 #ifndef CSCObjects_CSCGainsStudyHistograms_H
0002 #define CSCObjects_CSCGainsStudyHistograms_H
0003 
0004 /** \class CSCGainsStudyHistograms
0005  * Collection of histograms for plotting gain correction weights
0006  *  in each chambers used in MTCC.
0007  *
0008  * Author: D. Fortin  - UC Riverside
0009  */
0010 
0011 #include "TH1F.h"
0012 #include "TH2F.h"
0013 #include "TFile.h"
0014 #include "TString.h"
0015 #include <string>
0016 #include <iostream>
0017 
0018 class HCSCGains {
0019 public:
0020   /// Constructor from collection name
0021   HCSCGains(std::string name_) {
0022     TString N = name_.c_str();
0023     name = N;
0024     hGains = new TH1F("hGain_" + name, name, 200, 0.75, 3.75);
0025     hGaindiff = new TH1F("hGaindiff_" + name, name, 101, -0.101, 0.101);
0026     hGainvsch = new TH2F("hGainvsch_" + name, name, 600, 100.5, 700.5, 81, 0.7975, 1.225);
0027   }
0028 
0029   /// Constructor from collection name and TFile.
0030   HCSCGains(TString name_, TFile *file) {
0031     name = name_;
0032     hGains = (TH1F *)file->Get("hGains_" + name);
0033     hGaindiff = (TH1F *)file->Get("hGaindiff_" + name);
0034     hGainvsch = (TH2F *)file->Get("hGainvsch_" + name);
0035   }
0036 
0037   /// Destructor
0038   virtual ~HCSCGains() {
0039     delete hGains;
0040     delete hGaindiff;
0041     delete hGainvsch;
0042   }
0043 
0044   // Operations
0045 
0046   /// Fill all the histos
0047   void Fill(float weight, float weight0, int channel, int layer) {
0048     float weightx = weight + float(layer - 1) / 2.;
0049     hGains->Fill(weightx);
0050 
0051     int id = channel + 100 * layer;
0052     hGainvsch->Fill(id, weight);
0053     if (weight0 > 0. && weight > 0.)
0054       hGaindiff->Fill(weight - weight0);
0055   }
0056 
0057   /// Write all the histos to currently opened file
0058   void Write() {
0059     hGains->GetXaxis()->SetTitle("weight + (layer - 1)/2");
0060     hGains->Write();
0061 
0062     hGaindiff->GetXaxis()->SetTitle("weight diff between 2 adjacent strips");
0063     hGaindiff->Write();
0064 
0065     hGainvsch->GetXaxis()->SetTitle("strip_id + 100 x layer");
0066     hGainvsch->GetYaxis()->SetTitle("weight");
0067     hGainvsch->Write();
0068   }
0069 
0070   TH1F *hGains;
0071   TH1F *hGaindiff;
0072   TH2F *hGainvsch;
0073 
0074   TString name;
0075 };
0076 #endif