CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions | Variables
HFLightCalRand.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/HFLightCalRand.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 Fit3Peak (Double_t *x, Double_t *par)
 
void HistSpec (TH1F *hist, Double_t &mean, Double_t &rms, Double_t range=4)
 

Variables

Int_t EventN =0
 
Int_t NEvents
 
Int_t runNumb =0
 

Function Documentation

Double_t Fit3Peak ( Double_t *  x,
Double_t *  par 
)

Definition at line 142 of file HFLightCalRand.cc.

References reco::ParticleMasses::k0, NEvents, diffTwoXMLs::r1, diffTwoXMLs::r2, and mathSSE::sqrt().

Referenced by HFLightCalRand::endJob().

142  {
143 // Spectra fit function: Pedestal Gaussian + asymmetric 1PE + 2PE +3PE peaks
144 
145  Double_t sum,xx,A0,C0,r0,sigma0,mean1,sigma1,A1,C1,r1,mean2,sigma2,A2,C2,r2,mean3,sigma3,A3,C3,r3;
146 
147  const Double_t k0=2.0,k1=1.6, k2=2.0;
148 
149  xx=x[0];
150  sigma0 = par[2];
151  A0 = 2*NEvents/(2+2*par[0]+par[0]*par[0]+par[0]*par[0]*par[0]/3);
152  r0 = ((xx-par[1])/sigma0);
153  C0 = 1/(sigma0* TMath::Exp(-k0*k0/2)/k0 +
154  sigma0*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k0/1.41421)));
155  //sum = 1/(sqrt(2*3.14159)*par[2])*A0*TMath::Exp(-0.5*r0*r0);
156  if(r0 < k0) sum = C0*A0*TMath::Exp(-0.5*r0*r0);
157  else sum = C0*A0*TMath::Exp(0.5*k0*k0-k0*r0);
158 
159  mean1 = par[1]+par[3];
160  //sigma1 = par[4];
161  sigma1 = 1.547+0.125*par[3]+0.004042*par[3]*par[3];
162  sigma1 = (sigma1+(9.1347e-3+3.845e-2*par[3])*par[4]*2.0)*par[2];
163  A1 = A0*par[0];
164  C1 = 1/(sigma1* TMath::Exp(-k1*k1/2)/k1 +
165  sigma1*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k1/1.41421)));
166  r1 = ((xx-mean1)/sigma1);
167  if(r1 < k1) sum += C1*A1*TMath::Exp(-0.5*r1*r1);
168  else sum += C1*A1*TMath::Exp(0.5*k1*k1-k1*r1);
169 
170  mean2 = 2*par[3]+par[1];
171  sigma2 = sqrt(2*sigma1*sigma1 - par[2]*par[2]);
172  //A2 = A0*par[5]*par[0]*par[0]/2;
173  A2 = A0*par[0]*par[0]/2;
174  C2 = 1/(sigma2* TMath::Exp(-k2*k2/2)/k2 +
175  sigma2*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k2/1.41421)));
176  r2 = ((xx-mean2)/sigma2);
177  if(r2 < k2) sum += C2*A2*TMath::Exp(-0.5*r2*r2);
178  else sum += C2*A2*TMath::Exp(0.5*k2*k2-k2*r2);
179 
180  mean3 = 3*par[3]+par[1];
181  sigma3 = sqrt(3*sigma1*sigma1 - 2*par[2]*par[2]);
182  A3 = A0*par[0]*par[0]*par[0]/6;
183  C3 = 1/(sigma3*sqrt(2*3.14159));
184  r3 = ((xx-mean3)/sigma3);
185  sum += C3*A3*TMath::Exp(-0.5*r3*r3);
186 
187  return sum;
188 }
T sqrt(T t)
Definition: SSEVec.h:46
Int_t NEvents
Definition: DDAxes.h:10
void HistSpec ( TH1F *  hist,
Double_t &  mean,
Double_t &  rms,
Double_t  range = 4 
)

Definition at line 126 of file HFLightCalRand.cc.

Referenced by HFLightCalRand::endJob().

126  {
127  Double_t xmin,xmax;
128  mean=hist->GetMean();
129  rms=hist->GetRMS();
130  xmin=hist->GetXaxis()->GetXmin();
131  xmax=hist->GetXaxis()->GetXmax();
132  hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
133  mean=hist->GetMean();
134  rms=hist->GetRMS();
135  hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
136  mean=hist->GetMean();
137  rms=hist->GetRMS();
138  hist->SetAxisRange(xmin,xmax);
139  return;
140 }

Variable Documentation

Int_t EventN =0

Definition at line 36 of file HFLightCalRand.cc.

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

Int_t NEvents

Definition at line 36 of file HFLightCalRand.cc.

Referenced by HFLightCalRand::endJob(), and Fit3Peak().

Int_t runNumb =0

Definition at line 36 of file HFLightCalRand.cc.

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