CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions | Variables
HFLightCal.cc File Reference
#include <memory>
#include <string>
#include <iostream>
#include "TH1F.h"
#include "TH2F.h"
#include "TFile.h"
#include "math.h"
#include "TMath.h"
#include "TF1.h"
#include "CalibCalorimetry/HcalStandardModules/interface/HFLightCal.h"
#include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "CalibFormats/HcalObjects/interface/HcalDbService.h"
#include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
#include "CondFormats/HcalObjects/interface/HcalQIEShape.h"
#include "CondFormats/HcalObjects/interface/HcalQIECoder.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ESHandle.h"

Go to the source code of this file.

Functions

Double_t FitFun (Double_t *x, Double_t *par)
 
void HistSpecs (TH1F *hist, Double_t &mean, Double_t &rms, Double_t range=4)
 

Variables

Int_t eventN =0
 
Int_t itsmax [26][36][2]
 
Int_t itspinmax [8][3]
 
Int_t Nev
 
Int_t runN =0
 

Function Documentation

Double_t FitFun ( Double_t *  x,
Double_t *  par 
)

Definition at line 185 of file HFLightCal.cc.

References reco::ParticleMasses::k0, Nev, funct::pow(), diffTwoXMLs::r1, diffTwoXMLs::r2, and mathSSE::sqrt().

Referenced by HFLightCal::endJob().

185  {
186 // Spectra fit function: Pedestal Gaussian + asymmetric 1PE + 2PE +3PE peaks
187 
188  Double_t sum,xx,A0,C0,r0,sigma0,mean1,sigma1,A1,C1,r1,mean2,sigma2,A2,C2,r2,mean3,sigma3,A3,C3,r3;
189 
190  const Double_t k0=2.0, k1=1.6, k2=2.0;
191 
192  xx=x[0];
193  sigma0 = par[2];
194  A0 = 2*Nev/(2+2*par[0]+par[0]*par[0]+pow(par[0],3)/3+pow(par[0],4)/12+
195  pow(par[0],5)/60+pow(par[0],6)/360);
196  r0 = ((xx-par[1])/sigma0);
197  C0 = 1/(sigma0* TMath::Exp(-k0*k0/2)/k0 +
198  sigma0*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k0/1.41421)));
199  //sum = 1/(sqrt(2*3.14159)*par[2])*A0*TMath::Exp(-0.5*r0*r0);
200  if(r0 < k0) sum = C0*A0*TMath::Exp(-0.5*r0*r0);
201  else sum = C0*A0*TMath::Exp(0.5*k0*k0-k0*r0);
202 
203  mean1 = par[1]+par[3];
204  sigma1 = par[4];
205  //sigma1 = 1.547+0.125*par[3]+0.004042*par[3]*par[3];
206  //sigma1 = (sigma1+(9.1347e-3+3.845e-2*par[3])*par[4]*2.0)*par[2];
207  A1 = A0*par[0];
208  C1 = 1/(sigma1* TMath::Exp(-k1*k1/2)/k1 +
209  sigma1*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k1/1.41421)));
210  r1 = ((xx-mean1)/sigma1);
211  if(r1 < k1) sum += C1*A1*TMath::Exp(-0.5*r1*r1);
212  else sum += C1*A1*TMath::Exp(0.5*k1*k1-k1*r1);
213 
214  mean2 = 2*par[3]+par[1];
215  sigma2 = sqrt(2*sigma1*sigma1 - pow(par[2],2));
216  //A2 = A0*par[5]*par[0]*par[0]/2;
217  A2 = A0*par[0]*par[0]/2;
218  C2 = 1/(sigma2* TMath::Exp(-k2*k2/2)/k2 +
219  sigma2*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k2/1.41421)));
220  r2 = ((xx-mean2)/sigma2);
221  if(r2 < k2) sum += C2*A2*TMath::Exp(-0.5*r2*r2);
222  else sum += C2*A2*TMath::Exp(0.5*k2*k2-k2*r2);
223 
224  mean3 = 3*par[3]+par[1];
225  sigma3 = sqrt(3*sigma1*sigma1 - 2*pow(par[2],2));
226  A3 = A0*par[0]*par[0]*par[0]/6;
227  C3 = 1/(sigma3*sqrt(2*3.14159));
228  r3 = ((xx-mean3)/sigma3);
229  sum += C3*A3*TMath::Exp(-0.5*r3*r3);
230 
231  return sum;
232 }
T sqrt(T t)
Definition: SSEVec.h:46
Int_t Nev
Definition: HFLightCal.cc:35
x
Definition: VDTMath.h:216
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void HistSpecs ( TH1F *  hist,
Double_t &  mean,
Double_t &  rms,
Double_t  range = 4 
)

Definition at line 169 of file HFLightCal.cc.

Referenced by HFLightCal::endJob().

169  {
170  Double_t xmin,xmax;
171  mean=hist->GetMean();
172  rms=hist->GetRMS();
173  xmin=hist->GetXaxis()->GetXmin();
174  xmax=hist->GetXaxis()->GetXmax();
175  hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
176  mean=hist->GetMean();
177  rms=hist->GetRMS();
178  hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
179  mean=hist->GetMean();
180  rms=hist->GetRMS();
181  hist->SetAxisRange(xmin,xmax);
182  return;
183 }

Variable Documentation

Int_t eventN =0

Definition at line 35 of file HFLightCal.cc.

Referenced by HFLightCal::analyze(), and HFLightCal::endJob().

Int_t itsmax[26][36][2]

Definition at line 36 of file HFLightCal.cc.

Referenced by HFLightCal::analyze(), and HFLightCal::beginJob().

Int_t itspinmax[8][3]

Definition at line 37 of file HFLightCal.cc.

Referenced by HFLightCal::beginJob().

Int_t Nev

Definition at line 35 of file HFLightCal.cc.

Referenced by HFLightCal::endJob(), and FitFun().

Int_t runN =0

Definition at line 35 of file HFLightCal.cc.

Referenced by HFLightCal::analyze(), and HFLightCal::endJob().