ZIterativeAlgorithmWithFit

ZIterativeAlgorithmWithFitPlots

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145
#ifndef ZITERATIVEALGORITHMWITHFIT_H
#define ZITERATIVEALGORITHMWITHFIT_H

/* ******************************************
 * ZIterativeAlgorithmWithFit.h 
 *
 * Paolo Meridiani 06/03/2003
 ********************************************/

#include <TROOT.h>
#include <TClonesArray.h>
#include <vector>
#include <string>
#include <TH1.h>

#include "Calibration/Tools/interface/CalibElectron.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
/// Class that implements an iterative in situ calibration algorithm
/// using Z events
/** \class ZIterativeAlgorithmWithFit 
    Author: paolo.meridiani@roma1.infn.it
*/

#define nMaxIterations 50
#define nMaxChannels 250

class ZIterativeAlgorithmWithFit {
public:
  struct ZIterativeAlgorithmWithFitPlots {
    TH1* weightedRescaleFactor[nMaxIterations][nMaxChannels];
    TH1* unweightedRescaleFactor[nMaxIterations][nMaxChannels];
    TH1* weight[nMaxIterations][nMaxChannels];
  };

  /// Default constructor
  ZIterativeAlgorithmWithFit();

  /// Constructor with explicit iterations & exponent
  ZIterativeAlgorithmWithFit(const edm::ParameterSet& ps);
  //, unsigned int events);

  /// Assignment operator
  ZIterativeAlgorithmWithFit& operator=(const ZIterativeAlgorithmWithFit& r) { return *this; }

  /// Destructor
  virtual ~ZIterativeAlgorithmWithFit();

  bool resetIteration();

  bool iterate();

  bool addEvent(calib::CalibElectron*, calib::CalibElectron*, float);

  const ZIterativeAlgorithmWithFitPlots* getHistos() const { return thePlots_; }

  int getNumberOfIterations() const { return numberOfIterations_; }

  int getNumberOfChannels() const { return channels_; }

  const std::vector<float>& getOptimizedCoefficients() const { return optimizedCoefficients_; }

  const std::vector<float>& getOptimizedCoefficientsError() const { return optimizedCoefficientsError_; }

  const std::vector<float>& getOptimizedChiSquare() const { return optimizedChiSquare_; }

  const std::vector<int>& getOptimizedIterations() const { return optimizedIterations_; }

  const std::vector<float>& getWeightSum() const { return weight_sum_; }

  const std::vector<float>& getEpsilonSum() const { return calib_fac_; }

  //Helper Methods

  static inline float invMassCalc(float Energy1, float Eta1, float Phi1, float Energy2, float Eta2, float Phi2) {
    return (sqrt(2 * Energy1 * Energy2 * (1 - cosTheta12(Eta1, Phi1, Eta2, Phi2))));
  }

  static inline float cosTheta12(float Eta1, float Phi1, float Eta2, float Phi2) {
    return ((cos(Phi1 - Phi2) + sinh(Eta1) * sinh(Eta2)) / (cosh(Eta1) * cosh(Eta2)));
  }

  /*
  static TF1* gausfit(TH1F * histoou,double* par,double* errpar) {
    return gausfit(histoou,par,errpar,1.,2.);
  }
  */

  static void gausfit(
      TH1F* histoou, double* par, double* errpar, float nsigmalow, float nsigmaup, double* mychi2, int* iterations);

private:
  void addWeightsCorrections(unsigned int event_id);

  void getStatWeights(const std::string& file);

  float getEventWeight(unsigned int event_id);

  void recalculateWeightsEnergies();

  void recalculateMasses();

  void recalculateWeightsEnergies(calib::CalibElectron* electron);

  void getWeight(unsigned int evid, std::pair<calib::CalibElectron*, calib::CalibElectron*>, float);

  void getWeight(unsigned int evid, calib::CalibElectron* ele, float);

  void bookHistograms();

  ZIterativeAlgorithmWithFitPlots* thePlots_;

  int nCrystalCut_;

  unsigned int channels_;
  unsigned int totalEvents_;
  unsigned int numberOfIterations_;

  unsigned int currentEvent_;
  unsigned int currentIteration_;

  std::vector<std::pair<calib::CalibElectron*, calib::CalibElectron*> > electrons_;

  std::vector<float> optimizedCoefficients_;
  std::vector<float> optimizedCoefficientsError_;
  std::vector<float> calib_fac_;
  std::vector<float> weight_sum_;
  std::vector<float> massReco_;
  std::vector<float> optimizedChiSquare_;
  std::vector<int> optimizedIterations_;

  std::string massMethod;

  bool UseStatWeights_;
  std::string WeightFileName_;

  std::vector<float> StatWeights_;
  std::vector<float> Event_Weight_;

  TString calibType_;

  static const double M_Z_;
};

#endif  // ZIterativeAlgorithmWithFit_H