00001 #ifndef BackgroundHandler_cc
00002 #define BackgroundHandler_cc
00003
00004 #include "MuonAnalysis/MomentumScaleCalibration/interface/BackgroundHandler.h"
00005 #include <algorithm>
00006 #include <TF1.h>
00007 #include <iostream>
00008 #include <boost/foreach.hpp>
00009
00010 typedef reco::Particle::LorentzVector lorentzVector;
00011
00012 BackgroundHandler::BackgroundHandler( const std::vector<int> & identifiers,
00013 const std::vector<double> & leftWindowBorders,
00014 const std::vector<double> & rightWindowBorders,
00015 const double * ResMass,
00016 const double * massWindowHalfWidth )
00017 {
00018
00019
00020
00021
00022
00023 regToResHW_[0] = 0;
00024 regToResHW_[1] = 3;
00025 regToResHW_[2] = 5;
00026
00027
00028 resToReg_[0] = 0;
00029 resToReg_[1] = 1;
00030 resToReg_[2] = 1;
00031 resToReg_[3] = 1;
00032 resToReg_[4] = 2;
00033 resToReg_[5] = 2;
00034
00035
00036 consistencyCheck(identifiers, leftWindowBorders, rightWindowBorders);
00037
00038
00039 for( unsigned int i=0; i<6; ++i ) {
00040 double mass = ResMass[i];
00041 double lowerLimit = mass - massWindowHalfWidth[i];
00042 double upperLimit = mass + massWindowHalfWidth[i];
00043 resonanceWindow_.push_back( MassWindow( mass, lowerLimit, upperLimit,
00044 std::vector<unsigned int>(1,i),
00045 backgroundFunctionService(identifiers[resToReg_[i]],
00046 lowerLimit,
00047 upperLimit) ) );
00048 }
00049
00050
00051
00052
00053 double resMassForRegion[3];
00054 resMassForRegion[0] = ResMass[0];
00055 resMassForRegion[1] = (ResMass[1]+ResMass[2]+ResMass[3])/3;
00056 resMassForRegion[2] = (ResMass[4]+ResMass[5])/2;
00057
00058
00059 std::vector<std::vector<unsigned int> > indexes;
00060 indexes.push_back(std::vector<unsigned int>(1,0));
00061 indexes.push_back(std::vector<unsigned int>());
00062 for( int i=1; i<=3; ++i ) { indexes[1].push_back(i); }
00063 indexes.push_back(std::vector<unsigned int>());
00064 for( int i=4; i<=5; ++i ) { indexes[2].push_back(i); }
00065
00066 unsigned int i=0;
00067 typedef std::vector<unsigned int> indexVec;
00068 BOOST_FOREACH(const indexVec & index, indexes) {
00069
00070
00071
00072
00073 backgroundWindow_.push_back( MassWindow( resMassForRegion[i], leftWindowBorders[i], rightWindowBorders[i], index,
00074 backgroundFunctionService(identifiers[i], leftWindowBorders[i], rightWindowBorders[i] ) ) );
00075 ++i;
00076 }
00077
00078 initializeParNums();
00079 }
00080
00081 BackgroundHandler::~BackgroundHandler()
00082 {
00083 }
00084
00085 void BackgroundHandler::initializeParNums()
00086 {
00087
00088 parNumsRegions_[0] = 0;
00089 for( unsigned int i=1; i<backgroundWindow_.size() ; ++i ) {
00090 parNumsRegions_[i] = parNumsRegions_[i-1] + backgroundWindow_[i-1].backgroundFunction()->parNum();
00091 }
00092 parNumsResonances_[0] = parNumsRegions_[2]+backgroundWindow_[2].backgroundFunction()->parNum();
00093 for( unsigned int i=1; i<resonanceWindow_.size() ; ++i ) {
00094 parNumsResonances_[i] = parNumsResonances_[i-1] + resonanceWindow_[i-1].backgroundFunction()->parNum();
00095 }
00096 }
00097
00098 void BackgroundHandler::setParameters(double* Start, double* Step, double* Mini, double* Maxi, int* ind, TString* parname, const std::vector<double> & parBgr, const std::vector<int> & parBgrOrder, const int muonType)
00099 {
00100 std::vector<double>::const_iterator parBgrIt = parBgr.begin();
00101 std::vector<int>::const_iterator parBgrOrderIt = parBgrOrder.begin();
00102
00103 for( int iReg = 0; iReg < 3; ++iReg ) {
00104 int shift = parNumsRegions_[iReg];
00105 backgroundWindow_[iReg].backgroundFunction()->setParameters( &(Start[shift]), &(Step[shift]), &(Mini[shift]),
00106 &(Maxi[shift]), &(ind[shift]), &(parname[shift]),
00107 parBgrIt+shift, parBgrOrderIt+shift, muonType );
00108 }
00109 for( int iRes = 0; iRes < 6; ++iRes ) {
00110
00111 int shift = parNumsResonances_[iRes];
00112 resonanceWindow_[iRes].backgroundFunction()->setParameters( &(Start[shift]), &(Step[shift]), &(Mini[shift]),
00113 &(Maxi[shift]), &(ind[shift]), &(parname[shift]),
00114 parBgrIt+shift, parBgrOrderIt+shift, muonType );
00115 }
00116 }
00117
00118 bool BackgroundHandler::unlockParameter(const std::vector<int> & resfind, const unsigned int ipar)
00119 {
00120
00121 if( ipar < unsigned(parNumsRegions_[1]) && resfind[0] > 0 ) {
00122 return true;
00123 }
00124 if( ipar >= unsigned(parNumsRegions_[1]) && ipar < unsigned(parNumsRegions_[2]) && ( resfind[1] > 0 || resfind[2] > 0 || resfind[3] > 0 ) ) {
00125 return true;
00126 }
00127
00128 if( ipar >= unsigned(parNumsRegions_[2]) && ipar < unsigned(parNumsResonances_[0]) && ( resfind[4] > 0 || resfind[5] > 0 ) ) {
00129 return true;
00130 }
00131 return false;
00132 }
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146 std::pair<double, double> BackgroundHandler::windowBorders( const bool doBackgroundFit, const int ires )
00147 {
00148 if( doBackgroundFit ) {
00149
00150 return std::make_pair(backgroundWindow_[resToReg_[ires]].lowerBound(), backgroundWindow_[resToReg_[ires]].upperBound());
00151 }
00152 else {
00153
00154
00155 return std::make_pair(resonanceWindow_[ires].lowerBound(), resonanceWindow_[ires].upperBound());
00156 }
00157 }
00158
00159 double BackgroundHandler::resMass( const bool doBackgroundFit, const int ires )
00160 {
00161 if( doBackgroundFit ) {
00162
00163 return backgroundWindow_[resToReg_[ires]].mass();
00164 }
00165 else {
00166
00167 return resonanceWindow_[ires].mass();
00168 }
00169 }
00170
00171 void BackgroundHandler::rescale( std::vector<double> & parBgr, const double * ResMass, const double * massWindowHalfWidth,
00172 const std::vector<std::pair<reco::Particle::LorentzVector,reco::Particle::LorentzVector> > & muonPairs,
00173 const double & weight )
00174 {
00175 countEventsInAllWindows(muonPairs, weight);
00176
00177
00178
00179 unsigned int iRegion = 0;
00180 BOOST_FOREACH(MassWindow & backgroundWindow, backgroundWindow_)
00181 {
00182
00183 std::vector<double>::const_iterator parBgrIt = (parBgr.begin()+parNumsRegions_[iRegion]);
00184 TF1 * backgroundFunctionForIntegral = backgroundWindow.backgroundFunction()->functionForIntegral(parBgrIt);
00185
00186 double kOld = *parBgrIt;
00187 double Nbw = backgroundWindow.events();
00188 double Ibw = backgroundFunctionForIntegral->Integral(backgroundWindow.lowerBound(),
00189 backgroundWindow.upperBound());
00190
00191
00192 BOOST_FOREACH( unsigned int index, *(backgroundWindow.indexes()) )
00193 {
00194
00195 for( int iPar = 0; iPar < resonanceWindow_[index].backgroundFunction()->parNum(); ++iPar ) {
00196 parBgr[parNumsResonances_[index]+iPar] = parBgr[parNumsRegions_[resToReg_[index]]+iPar];
00197 }
00198
00199 double Irw = backgroundFunctionForIntegral->Integral(resonanceWindow_[index].lowerBound(),
00200 resonanceWindow_[index].upperBound());
00201 double Nrw = resonanceWindow_[index].events();
00202
00203
00204 if( Nrw != 0 ) parBgr[parNumsResonances_[index]] = kOld*Nbw/Nrw*Irw;
00205 else parBgr[parNumsResonances_[index]] = 0.;
00206
00207
00208 if( parBgr[parNumsResonances_[index]] > 1. ) parBgr[parNumsResonances_[index]] = 1.;
00209
00210 double kNew = parBgr[parNumsResonances_[index]];
00211 std::cout << "For resonance = " << index << std::endl;
00212 std::cout << "backgroundWindow.lowerBound = " << backgroundWindow.lowerBound() << std::endl;
00213 std::cout << "backgroundWindow.upperBound = " << backgroundWindow.upperBound() << std::endl;
00214 std::cout << "parNumsResonances_["<<index<<"] = " << parNumsResonances_[index] << std::endl;
00215 std::cout << "Nbw = " << Nbw << ", Ibw = " << Ibw << std::endl;
00216 std::cout << "Nrw = " << Nrw << ", Irw = " << Irw << std::endl;
00217 std::cout << "k = " << kOld << ", k' = " << parBgr[parNumsResonances_[index]] << std::endl;
00218 std::cout << "background fraction in background window = Nbw*k = " << Nbw*kOld << std::endl;
00219 std::cout << "background fraction in resonance window = Nrw*k' = " << Nrw*kNew << std::endl;
00220 }
00221 ++iRegion;
00222 delete backgroundFunctionForIntegral;
00223 }
00224 }
00225
00226 std::pair<double, double> BackgroundHandler::backgroundFunction( const bool doBackgroundFit,
00227 const double * parval, const int resTotNum, const int ires,
00228 const bool * resConsidered, const double * ResMass, const double ResHalfWidth[],
00229 const int MuonType, const double & mass,
00230 const double & eta1, const double & eta2 )
00231 {
00232 if( doBackgroundFit ) {
00233
00234 int iReg = resToReg_[ires];
00235
00236 return std::make_pair( parval[parNumsRegions_[iReg]] * backgroundWindow_[iReg].backgroundFunction()->fracVsEta(&(parval[parNumsRegions_[iReg]]), eta1, eta2 ),
00237 (*(backgroundWindow_[iReg].backgroundFunction()))( &(parval[parNumsRegions_[iReg]]), mass, eta1, eta2 ) );
00238
00239
00240 }
00241
00242
00243
00244 return std::make_pair( parval[parNumsResonances_[ires]] * resonanceWindow_[ires].backgroundFunction()->fracVsEta(&(parval[parNumsResonances_[ires]]), eta1, eta2),
00245 (*(resonanceWindow_[ires].backgroundFunction()))( &(parval[parNumsResonances_[ires]]), mass, eta1, eta2 ) );
00246 }
00247
00248 void BackgroundHandler::countEventsInAllWindows(const std::vector<std::pair<reco::Particle::LorentzVector,reco::Particle::LorentzVector> > & muonPairs,
00249 const double & weight)
00250 {
00251
00252 BOOST_FOREACH(MassWindow & resonanceWindow, resonanceWindow_) {
00253 resonanceWindow.resetCounter();
00254 }
00255
00256 BOOST_FOREACH(MassWindow & backgroundWindow, backgroundWindow_) {
00257 backgroundWindow.resetCounter();
00258 }
00259
00260
00261 std::pair<lorentzVector,lorentzVector> muonPair;
00262 BOOST_FOREACH(muonPair, muonPairs) {
00263
00264 BOOST_FOREACH(MassWindow & resonanceWindow, resonanceWindow_) {
00265 resonanceWindow.count((muonPair.first + muonPair.second).mass(), weight);
00266 }
00267
00268 BOOST_FOREACH(MassWindow & backgroundWindow, backgroundWindow_) {
00269 backgroundWindow.count((muonPair.first + muonPair.second).mass(), weight);
00270 }
00271 }
00272 }
00273
00274 void BackgroundHandler::consistencyCheck(const std::vector<int> & identifiers,
00275 const std::vector<double> & leftWindowBorders,
00276 const std::vector<double> & rightWindowBorders) const throw(cms::Exception)
00277 {
00278 if( leftWindowBorders.size() != rightWindowBorders.size() ) {
00279 throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: leftWindowBorders.size() = " << leftWindowBorders.size()
00280 << " != rightWindowBorders.size() = " << rightWindowBorders.size() << std::endl;
00281 }
00282 if( leftWindowBorders.size() != 3 ) {
00283 throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: leftWindowBorders.size() = rightWindowBorders.size() = "
00284 << leftWindowBorders.size() << " != 3" << std::endl;
00285 }
00286 if( identifiers.size() != 3 ) {
00287 throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: identifiers must match the number of regions = 3" << std::endl;
00288 }
00289 }
00290
00291 #endif // BackgroundHandler_cc