File indexing completed on 2024-04-06 11:59:29
0001 #ifndef ZITERATIVEALGORITHMWITHFIT_H
0002 #define ZITERATIVEALGORITHMWITHFIT_H
0003
0004
0005
0006
0007
0008
0009
0010 #include <TROOT.h>
0011 #include <TClonesArray.h>
0012 #include <vector>
0013 #include <string>
0014 #include <TH1.h>
0015
0016 #include "Calibration/Tools/interface/CalibElectron.h"
0017
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019
0020
0021
0022
0023
0024
0025 #define nMaxIterations 50
0026 #define nMaxChannels 250
0027
0028 class ZIterativeAlgorithmWithFit {
0029 public:
0030 struct ZIterativeAlgorithmWithFitPlots {
0031 TH1* weightedRescaleFactor[nMaxIterations][nMaxChannels];
0032 TH1* unweightedRescaleFactor[nMaxIterations][nMaxChannels];
0033 TH1* weight[nMaxIterations][nMaxChannels];
0034 };
0035
0036
0037 ZIterativeAlgorithmWithFit();
0038
0039
0040 ZIterativeAlgorithmWithFit(const edm::ParameterSet& ps);
0041
0042
0043
0044 ZIterativeAlgorithmWithFit& operator=(const ZIterativeAlgorithmWithFit& r) { return *this; }
0045
0046
0047 virtual ~ZIterativeAlgorithmWithFit();
0048
0049 bool resetIteration();
0050
0051 bool iterate();
0052
0053 bool addEvent(calib::CalibElectron*, calib::CalibElectron*, float);
0054
0055 const ZIterativeAlgorithmWithFitPlots* getHistos() const { return thePlots_; }
0056
0057 int getNumberOfIterations() const { return numberOfIterations_; }
0058
0059 int getNumberOfChannels() const { return channels_; }
0060
0061 const std::vector<float>& getOptimizedCoefficients() const { return optimizedCoefficients_; }
0062
0063 const std::vector<float>& getOptimizedCoefficientsError() const { return optimizedCoefficientsError_; }
0064
0065 const std::vector<float>& getOptimizedChiSquare() const { return optimizedChiSquare_; }
0066
0067 const std::vector<int>& getOptimizedIterations() const { return optimizedIterations_; }
0068
0069 const std::vector<float>& getWeightSum() const { return weight_sum_; }
0070
0071 const std::vector<float>& getEpsilonSum() const { return calib_fac_; }
0072
0073
0074
0075 static inline float invMassCalc(float Energy1, float Eta1, float Phi1, float Energy2, float Eta2, float Phi2) {
0076 return (sqrt(2 * Energy1 * Energy2 * (1 - cosTheta12(Eta1, Phi1, Eta2, Phi2))));
0077 }
0078
0079 static inline float cosTheta12(float Eta1, float Phi1, float Eta2, float Phi2) {
0080 return ((cos(Phi1 - Phi2) + sinh(Eta1) * sinh(Eta2)) / (cosh(Eta1) * cosh(Eta2)));
0081 }
0082
0083
0084
0085
0086
0087
0088
0089 static void gausfit(
0090 TH1F* histoou, double* par, double* errpar, float nsigmalow, float nsigmaup, double* mychi2, int* iterations);
0091
0092 private:
0093 void addWeightsCorrections(unsigned int event_id);
0094
0095 void getStatWeights(const std::string& file);
0096
0097 float getEventWeight(unsigned int event_id);
0098
0099 void recalculateWeightsEnergies();
0100
0101 void recalculateMasses();
0102
0103 void recalculateWeightsEnergies(calib::CalibElectron* electron);
0104
0105 void getWeight(unsigned int evid, std::pair<calib::CalibElectron*, calib::CalibElectron*>, float);
0106
0107 void getWeight(unsigned int evid, calib::CalibElectron* ele, float);
0108
0109 void bookHistograms();
0110
0111 ZIterativeAlgorithmWithFitPlots* thePlots_;
0112
0113 int nCrystalCut_;
0114
0115 unsigned int channels_;
0116 unsigned int totalEvents_;
0117 unsigned int numberOfIterations_;
0118
0119 unsigned int currentEvent_;
0120 unsigned int currentIteration_;
0121
0122 std::vector<std::pair<calib::CalibElectron*, calib::CalibElectron*> > electrons_;
0123
0124 std::vector<float> optimizedCoefficients_;
0125 std::vector<float> optimizedCoefficientsError_;
0126 std::vector<float> calib_fac_;
0127 std::vector<float> weight_sum_;
0128 std::vector<float> massReco_;
0129 std::vector<float> optimizedChiSquare_;
0130 std::vector<int> optimizedIterations_;
0131
0132 std::string massMethod;
0133
0134 bool UseStatWeights_;
0135 std::string WeightFileName_;
0136
0137 std::vector<float> StatWeights_;
0138 std::vector<float> Event_Weight_;
0139
0140 TString calibType_;
0141
0142 static const double M_Z_;
0143 };
0144
0145 #endif