CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/MuonAnalysis/MomentumScaleCalibration/plugins/MuScleFitUtils.cc

Go to the documentation of this file.
00001 
00007 // Some notes:
00008 // - M(Z) after simulation needs to be extracted as a function of |y_Z| in order to be
00009 //   a better reference point for calibration. In fact, the variation of PDF with y_Z
00010 //   in production is sizable <---- need to check though.
00011 // - ResHalfWidth needs to be optimized - this depends on the level of background.
00012 // - Background parametrization still to be worked on, so far only a constant (type=1, and
00013 //   parameter 2 fixed to 0) works.
00014 // - weights have to be assigned to dimuon mass values in regions where different resonances
00015 //   overlap, and one has to decide which resonance mass to assign the event to - this until
00016 //   we implement in the fitter a sum of probabilities of an event to belong to different
00017 //   resonances. The weight has to depend on the mass and has relative cross sections of
00018 //   Y(1S), 2S, 3S as parameters. Some overlap is also expected in the J/psi-Psi(2S) region
00019 //   when reconstructing masses with Standalone muons.
00020 //
00021 //   MODS 7/7/08 TD:
00022 //   - changed parametrization of resolution in Pt: from sigma_pt = a*Pt + b*|eta| to
00023 //                                                       sigma_pt = (a*Pt + b*|eta|)*Pt
00024 //                                                  which is more correct (I hope)
00025 //   - changed parametrization of resolution in cotgth: from sigma_cotgth = f(eta) to f(cotgth)
00026 // --------------------------------------------------------------------------------------------
00027 
00028 #include "MuScleFitUtils.h"
00029 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00030 #include "SimDataFormats/Track/interface/SimTrack.h"
00031 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00032 #include "DataFormats/Math/interface/LorentzVector.h"
00033 #include "TString.h"
00034 #include "TFile.h"
00035 #include "TTree.h"
00036 #include "TCanvas.h"
00037 #include "TH2F.h"
00038 #include "TF1.h"
00039 #include "TF2.h"
00040 #include <iostream>
00041 #include <fstream>
00042 #include <memory> // to use the auto_ptr
00043 
00044 // Includes the definitions of all the bias and scale functions
00045 // These functions are selected in the constructor according
00046 // to the input parameters.
00047 #include "MuonAnalysis/MomentumScaleCalibration/interface/Functions.h"
00048 
00049 // To use callgrind for code profiling uncomment also the following define.
00050 //#define USE_CALLGRIND
00051 #ifdef USE_CALLGRIND
00052 #include "valgrind/callgrind.h"
00053 #endif
00054 
00055 // Lorenzian Peak function
00056 // -----------------------
00057 Double_t lorentzianPeak (Double_t *x, Double_t *par) {
00058   return (0.5*par[0]*par[1]/TMath::Pi()) /
00059   TMath::Max(1.e-10,(x[0]-par[2])*(x[0]-par[2]) + .25*par[1]*par[1]);
00060 }
00061 
00062 // Gaussian function
00063 // -----------------
00064 Double_t Gaussian (Double_t *x, Double_t *par) {
00065   return par[0]*exp(-0.5*((x[0]-par[1])/par[2])*((x[0]-par[1])/par[2]));
00066 }
00067 
00068 // Array with number of parameters in the fitting functions
00069 // (not currently in use)
00070 // --------------------------------------------------------
00071 //const int nparsResol[2] = {6, 4};
00072 //const int nparsScale[13] = {2, 2, 2, 3, 3, 3, 4, 4, 2, 3, 4, 6, 8};
00073 //const int nparsBgr[3] = {1, 2, 3};
00074 
00075 // Quantities used for h-value computation
00076 // ---------------------------------------
00077 double mzsum;
00078 double isum;
00079 double f[11][100];
00080 double g[11][100];
00081 
00082 // Lorentzian convoluted with a gaussian:
00083 // --------------------------------------
00084 TF1 * GL = new TF1 ("GL",
00085                     "0.5/3.1415926*[0]/(pow(x-[1],2)+pow(0.5*[0],2))*exp(-0.5*pow((x-[2])/[3],2))/([3]*sqrt(6.283185))",
00086                     0, 1000);
00087 
00088 TF2 * GL2= new TF2 ("GL2",
00089   "0.5/3.1415926*[0]/(pow(x-[1],2)+pow(0.5*[0],2))*exp(-0.5*pow((x-y)/[2],2))/([2]*sqrt(6.283185))",
00090   0, 200, 0, 200);
00091 
00092 // // Lorentzian convoluted with a gaussian over a linear background:
00093 // // ---------------------------------------------------------------
00094 // TF1 * GLBL = new TF1 ("GLBL",
00095 //   "0.5/3.1415926*[0]/(pow(x-[1],2)+pow(0.5*[0],2))*exp(-0.5*pow((x-[2])/[3],2))/([3]*sqrt(6.283185))+[4]+[5]*x",
00096 //   0, 1000);
00097 
00098 // // Lorentzian convoluted with a gaussian over an exponential background:
00099 // // ---------------------------------------------------------------
00100 // TF1 * GLBE = new TF1 ("GLBE",
00101 //   "0.5/3.1415926*[0]/(pow(x-[1],2)+pow(0.5*[0],2))*exp(-0.5*pow((x-[2])/[3],2))/([3]*sqrt(6.283185))+exp([4]+[5]*x)",
00102 //   0, 1000);
00103 
00104 std::vector<int> MuScleFitUtils::doResolFit;
00105 std::vector<int> MuScleFitUtils::doScaleFit;
00106 std::vector<int> MuScleFitUtils::doCrossSectionFit;
00107 std::vector<int> MuScleFitUtils::doBackgroundFit;
00108 
00109 int MuScleFitUtils::minuitLoop_ = 0;
00110 TH1D* MuScleFitUtils::likelihoodInLoop_ = 0;
00111 TH1D* MuScleFitUtils::signalProb_ = 0;
00112 TH1D* MuScleFitUtils::backgroundProb_ = 0;
00113 
00114 bool MuScleFitUtils::duringMinos_ = false;
00115 
00116 const int MuScleFitUtils::totalResNum = 6;
00117 
00118 int MuScleFitUtils::SmearType = 0;
00119 smearFunctionBase * MuScleFitUtils::smearFunction = 0;
00120 int MuScleFitUtils::BiasType  = 0;
00121 // No error, we take functions from the same group for bias and scale.
00122 scaleFunctionBase<std::vector<double> > * MuScleFitUtils::biasFunction = 0;
00123 int MuScleFitUtils::ResolFitType = 0;
00124 resolutionFunctionBase<double *> * MuScleFitUtils::resolutionFunction = 0;
00125 resolutionFunctionBase<std::vector<double> > * MuScleFitUtils::resolutionFunctionForVec = 0;
00126 int MuScleFitUtils::ScaleFitType = 0;
00127 scaleFunctionBase<double*> * MuScleFitUtils::scaleFunction = 0;
00128 scaleFunctionBase<std::vector<double> > * MuScleFitUtils::scaleFunctionForVec = 0;
00129 int MuScleFitUtils::BgrFitType   = 0;
00130 
00131 CrossSectionHandler * MuScleFitUtils::crossSectionHandler;
00132 // const int MuScleFitUtils::backgroundFunctionsRegions = 3;
00133 // backgroundFunctionBase * MuScleFitUtils::backgroundFunctionForRegion[MuScleFitUtils::backgroundFunctionsRegions];
00134 // backgroundFunctionBase * MuScleFitUtils::backgroundFunction[MuScleFitUtils::totalResNum];
00135 BackgroundHandler * MuScleFitUtils::backgroundHandler;
00136 std::vector<double> MuScleFitUtils::parBias;
00137 std::vector<double> MuScleFitUtils::parSmear;
00138 std::vector<double> MuScleFitUtils::parResol;
00139 std::vector<double> MuScleFitUtils::parResolStep;
00140 std::vector<double> MuScleFitUtils::parResolMin;
00141 std::vector<double> MuScleFitUtils::parResolMax;
00142 std::vector<double> MuScleFitUtils::parScale;
00143 std::vector<double> MuScleFitUtils::parScaleStep;
00144 std::vector<double> MuScleFitUtils::parScaleMin;
00145 std::vector<double> MuScleFitUtils::parScaleMax;
00146 std::vector<double> MuScleFitUtils::parCrossSection;
00147 std::vector<double> MuScleFitUtils::parBgr;
00148 std::vector<int> MuScleFitUtils::parResolFix;
00149 std::vector<int> MuScleFitUtils::parScaleFix;
00150 std::vector<int> MuScleFitUtils::parCrossSectionFix;
00151 std::vector<int> MuScleFitUtils::parBgrFix;
00152 std::vector<int> MuScleFitUtils::parResolOrder;
00153 std::vector<int> MuScleFitUtils::parScaleOrder;
00154 std::vector<int> MuScleFitUtils::parCrossSectionOrder;
00155 std::vector<int> MuScleFitUtils::parBgrOrder;
00156 
00157 std::vector<int> MuScleFitUtils::resfind;
00158 int MuScleFitUtils::debug = 0;
00159 
00160 bool MuScleFitUtils::ResFound = false;
00161 int MuScleFitUtils::goodmuon = 0;
00162 int MuScleFitUtils::counter_resprob = 0;
00163 
00164 std::vector<std::vector<double> > MuScleFitUtils::parvalue;
00165 
00166 int MuScleFitUtils::FitStrategy = 1; // Strategy in likelihood fit (1 or 2)
00167 bool MuScleFitUtils::speedup = false; // Whether to cut corners (no sim study, fewer histos)
00168 
00169 std::vector<std::pair<lorentzVector,lorentzVector> > MuScleFitUtils::SavedPair; // Pairs of reconstructed muons making resonances
00170 std::vector<std::pair<lorentzVector,lorentzVector> > MuScleFitUtils::ReducedSavedPair; // Pairs of reconstructed muons making resonances inside smaller windows
00171 std::vector<std::pair<lorentzVector,lorentzVector> > MuScleFitUtils::genPair; // Pairs of generated muons making resonances
00172 std::vector<std::pair<lorentzVector,lorentzVector> > MuScleFitUtils::simPair; // Pairs of simulated muons making resonances
00173 
00174 // Smearing parameters
00175 // -------------------
00176 double MuScleFitUtils::x[][10000];
00177 
00178 // Probability matrices and normalization values
00179 // ---------------------------------------------
00180 int MuScleFitUtils::nbins = 1000;
00181 double MuScleFitUtils::GLZValue[][1001][1001];
00182 double MuScleFitUtils::GLZNorm[][1001];
00183 double MuScleFitUtils::GLValue[][1001][1001];
00184 double MuScleFitUtils::GLNorm[][1001];
00185 double MuScleFitUtils::ResMaxSigma[];
00186 
00187 // Masses and widths from PDG 2006, half widths to be revised
00188 // NB in particular, halfwidths have to be made a function of muonType
00189 // -------------------------------------------------------------------
00190 const double MuScleFitUtils::mMu2 = 0.011163612;
00191 const double MuScleFitUtils::muMass = 0.105658;
00192 double MuScleFitUtils::ResHalfWidth[];
00193 double MuScleFitUtils::massWindowHalfWidth[][6];
00194 int MuScleFitUtils::MuonType;
00195 int MuScleFitUtils::MuonTypeForCheckMassWindow;
00196 
00197 double MuScleFitUtils::ResGamma[] = {2.4952, 0.000020, 0.000032, 0.000054, 0.000317, 0.0000932 };
00198 // ATTENTION:
00199 // This is left because the values are used by the BackgroundHandler to define the center of the regions windows,
00200 // but the values used in the code are read computed using the probability histograms ranges.
00201 // The histograms are read after the initialization of the BackgroundHandler (this can be improved so that
00202 // the background handler too could use the new values).
00203 // At this time the values are consistent.
00204 double MuScleFitUtils::ResMinMass[] = {-99, -99, -99, -99, -99, -99};
00205 double MuScleFitUtils::ResMass[] = {91.1876, 10.3552, 10.0233, 9.4603, 3.68609, 3.0969};
00206 // From Summer08 generator production TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/ProductionSummer2008
00207 // - Z->mumu              1.233 nb
00208 // - Upsilon3S->mumu      0.82  nb
00209 // - Upsilon2S->mumu      6.33  nb
00210 // - Upsilon1S->mumu     13.9   nb
00211 // - Prompt Psi2S->mumu   2.169 nb
00212 // - Prompt J/Psi->mumu 127.2   nb
00213 // double MuScleFitUtils::crossSection[] = {1.233, 0.82, 6.33, 13.9, 2.169, 127.2};
00214 // double MuScleFitUtils::crossSection[] = {1.233, 2.07, 6.33, 13.9, 2.169, 127.2};
00215 
00216 unsigned int MuScleFitUtils::loopCounter = 5;
00217 
00218 // According to the pythia manual, there is only a code for the Upsilon and Upsilon'. It does not distinguish
00219 // between Upsilon(2S) and Upsilon(3S)
00220 const unsigned int MuScleFitUtils::motherPdgIdArray[] = {23, 100553, 100553, 553, 100443, 443};
00221 
00222 // double MuScleFitUtils::leftWindowFactor = 1.;
00223 // double MuScleFitUtils::rightWindowFactor = 1.;
00224 
00225 // double MuScleFitUtils::internalLeftWindowFactor = 1.;
00226 // double MuScleFitUtils::internalRightWindowFactor = 1.;
00227 
00228 // int MuScleFitUtils::backgroundWindowEvents_ = 0;
00229 // int MuScleFitUtils::resonanceWindowEvents_ = 0;
00230 
00231 // double MuScleFitUtils::oldEventsOutInRatio_ = 0.;
00232 
00233 bool MuScleFitUtils::scaleFitNotDone_ = true;
00234 
00235 bool MuScleFitUtils::sherpa_ = false;
00236 
00237 bool MuScleFitUtils::useProbsFile_ = true;
00238 
00239 bool MuScleFitUtils::rapidityBinsForZ_ = true;
00240 
00241 bool MuScleFitUtils::separateRanges_ = true;
00242 double MuScleFitUtils::minMuonPt_ = 0.;
00243 double MuScleFitUtils::maxMuonPt_ = 100000000.;
00244 double MuScleFitUtils::minMuonEtaFirstRange_ = -6.;
00245 double MuScleFitUtils::maxMuonEtaFirstRange_ = 6.;
00246 double MuScleFitUtils::minMuonEtaSecondRange_ = -100.;
00247 double MuScleFitUtils::maxMuonEtaSecondRange_ = 100.;
00248 double MuScleFitUtils::deltaPhiMinCut_ = -100.;
00249 double MuScleFitUtils::deltaPhiMaxCut_ = 100.;
00250 
00251 bool MuScleFitUtils::debugMassResol_;
00252 MuScleFitUtils::massResolComponentsStruct MuScleFitUtils::massResolComponents;
00253 
00254 bool MuScleFitUtils::normalizeLikelihoodByEventNumber_ = true;
00255 TMinuit * MuScleFitUtils::rminPtr_ = 0;
00256 double MuScleFitUtils::oldNormalization_ = 0.;
00257 unsigned int MuScleFitUtils::normalizationChanged_ = 0;
00258 
00259 bool MuScleFitUtils::startWithSimplex_;
00260 bool MuScleFitUtils::computeMinosErrors_;
00261 bool MuScleFitUtils::minimumShapePlots_;
00262 
00263 int MuScleFitUtils::iev_ = 0;
00265 
00266 // Find the best simulated resonance from a vector of simulated muons (SimTracks)
00267 // and return its decay muons
00268 // ------------------------------------------------------------------------------
00269 std::pair<SimTrack,SimTrack> MuScleFitUtils::findBestSimuRes (const std::vector<SimTrack>& simMuons) {
00270 
00271   std::pair<SimTrack, SimTrack> simMuFromBestRes;
00272   double maxprob = -0.1;
00273 
00274   // Double loop on muons
00275   // --------------------
00276   for (std::vector<SimTrack>::const_iterator simMu1=simMuons.begin(); simMu1!=simMuons.end(); simMu1++) {
00277     for (std::vector<SimTrack>::const_iterator simMu2=simMu1+1; simMu2!=simMuons.end(); simMu2++) {
00278       if (((*simMu1).charge()*(*simMu2).charge())>0) {
00279         continue; // this also gets rid of simMu1==simMu2...
00280       }
00281       // Choose the best resonance using its mass. Check Z, Y(3S,2S,1S), Psi(2S), J/Psi in order
00282       // ---------------------------------------------------------------------------------------
00283       double mcomb = ((*simMu1).momentum()+(*simMu2).momentum()).mass();
00284       double Y = ((*simMu1).momentum()+(*simMu2).momentum()).Rapidity();
00285       for (int ires=0; ires<6; ires++) {
00286         if (resfind[ires]>0) {
00287           double prob = massProb( mcomb, Y, ires, 0. );
00288           if (prob>maxprob) {
00289             simMuFromBestRes.first = (*simMu1);
00290             simMuFromBestRes.second = (*simMu2);
00291             maxprob = prob;
00292           }
00293         }
00294       }
00295     }
00296   }
00297 
00298   // Return most likely combination of muons making a resonance
00299   // ----------------------------------------------------------
00300   return simMuFromBestRes;
00301 }
00302 
00303 // Find the best reconstructed resonance from a collection of reconstructed muons
00304 // (MuonCollection) and return its decay muons
00305 // ------------------------------------------------------------------------------
00306 std::pair<lorentzVector,lorentzVector> MuScleFitUtils::findBestRecoRes( const std::vector<reco::LeafCandidate>& muons ){
00307   // NB this routine returns the resonance, but it also sets the ResFound flag, which
00308   // is used in MuScleFit to decide whether to use the event or not.
00309   // --------------------------------------------------------------------------------
00310   if (debug>0) std::cout << "In findBestRecoRes" << std::endl;
00311   ResFound = false;
00312   std::pair<lorentzVector, lorentzVector> recMuFromBestRes;
00313 
00314   // Choose the best resonance using its mass probability
00315   // ----------------------------------------------------
00316   double maxprob = -0.1;
00317   double minDeltaMass = 999999;
00318   std::pair<reco::LeafCandidate,reco::LeafCandidate> bestMassMuons;
00319   for (std::vector<reco::LeafCandidate>::const_iterator Muon1=muons.begin(); Muon1!=muons.end(); ++Muon1) {
00320     //rc2010
00321     if (debug>0) std::cout << "muon_1_charge:"<<(*Muon1).charge() << std::endl;
00322     for (std::vector<reco::LeafCandidate>::const_iterator Muon2=Muon1+1; Muon2!=muons.end(); ++Muon2) {
00323    //rc2010
00324       if (debug>0) std::cout << "after_2" << std::endl;
00325       if (((*Muon1).charge()*(*Muon2).charge())>0) {
00326         continue; // This also gets rid of Muon1==Muon2...
00327       }
00328       // To allow the selection of ranges at negative and positive eta independently we define two
00329       // ranges of eta: (minMuonEtaFirstRange_, maxMuonEtaFirstRange_) and (minMuonEtaSecondRange_, maxMuonEtaSecondRange_).
00330       // If the interval selected is simmetric, one only needs to specify the first range. The second has
00331       // default values that accept all muons (minMuonEtaSecondRange_ = -100., maxMuonEtaSecondRange_ = 100.).
00332       double pt1 = (*Muon1).p4().Pt();
00333       double pt2 = (*Muon2).p4().Pt();
00334       double eta1 = (*Muon1).p4().Eta();
00335       double eta2 = (*Muon2).p4().Eta();
00336       if( pt1 >= minMuonPt_ && pt1 < maxMuonPt_ &&
00337           pt2 >= minMuonPt_ && pt2 < maxMuonPt_ &&
00338           ( (eta1 >= minMuonEtaFirstRange_ && eta1 < maxMuonEtaFirstRange_ &&
00339              eta2 >= minMuonEtaFirstRange_ && eta2 < maxMuonEtaFirstRange_) ||
00340             (eta1 >= minMuonEtaSecondRange_ && eta1 < maxMuonEtaSecondRange_ &&
00341              eta2 >= minMuonEtaSecondRange_ && eta2 < maxMuonEtaSecondRange_) ) ) {
00342         double mcomb = ((*Muon1).p4()+(*Muon2).p4()).mass();
00343         double Y = ((*Muon1).p4()+(*Muon2).p4()).Rapidity();
00344         if (debug>1) {
00345           std::cout<<"muon1 "<<(*Muon1).p4().Px()<<", "<<(*Muon1).p4().Py()<<", "<<(*Muon1).p4().Pz()<<", "<<(*Muon1).p4().E()<<std::endl;
00346           std::cout<<"muon2 "<<(*Muon2).p4().Px()<<", "<<(*Muon2).p4().Py()<<", "<<(*Muon2).p4().Pz()<<", "<<(*Muon2).p4().E()<<std::endl;
00347           std::cout<<"mcomb "<<mcomb<<std::endl;}
00348         double massResol = 0.;
00349         if( useProbsFile_ ) {
00350           massResol = massResolution ((*Muon1).p4(), (*Muon2).p4(), parResol);
00351         }
00352         double prob = 0;
00353         for( int ires=0; ires<6; ires++ ) {
00354           if( resfind[ires]>0 ) {
00355             if( useProbsFile_ ) {
00356               prob = massProb( mcomb, Y, ires, massResol );
00357             }
00358             if( prob>maxprob ) {
00359               if( (*Muon1).charge()<0 ) { // store first the mu minus and then the mu plus
00360                 recMuFromBestRes.first = (*Muon1).p4();
00361                 recMuFromBestRes.second = (*Muon2).p4();
00362               } else {
00363                 recMuFromBestRes.first = (*Muon2).p4();
00364                 recMuFromBestRes.second = (*Muon1).p4();
00365               }
00366               ResFound = true; // NNBB we accept "resonances" even outside mass bounds
00367               maxprob = prob;
00368             }
00369             // if( ResMass[ires] == 0 ) {
00370             //   std::cout << "Error: ResMass["<<ires<<"] = " << ResMass[ires] << std::endl;
00371             //   exit(1);
00372             // }
00373             double deltaMass = fabs(mcomb-ResMass[ires])/ResMass[ires];
00374             if( deltaMass<minDeltaMass ){
00375               bestMassMuons = std::make_pair((*Muon1),(*Muon2));
00376               minDeltaMass = deltaMass;
00377             }
00378           }
00379         }
00380       }
00381     }
00382   }
00383   //If outside mass window (maxprob==0) then take the two muons with best invariant mass
00384   //(anyway they will not be used in the likelihood calculation, only to fill plots)
00385   if(!maxprob){
00386     if(bestMassMuons.first.charge()<0){
00387       recMuFromBestRes.first = bestMassMuons.first.p4();
00388       recMuFromBestRes.second = bestMassMuons.second.p4();
00389     }
00390     else{
00391       recMuFromBestRes.second = bestMassMuons.first.p4();
00392       recMuFromBestRes.first = bestMassMuons.second.p4();
00393     }
00394   }
00395   return recMuFromBestRes;
00396 }
00397 
00398 // Resolution smearing function called to worsen muon Pt resolution at start
00399 // -------------------------------------------------------------------------
00400 lorentzVector MuScleFitUtils::applySmearing (const lorentzVector& muon)
00401 {
00402   double pt = muon.Pt();
00403   double eta = muon.Eta();
00404   double phi = muon.Phi();
00405   double E = muon.E();
00406 
00407   double y[7];
00408   for (int i=0; i<SmearType+3; i++) {
00409     y[i] = x[i][goodmuon%10000];
00410   }
00411 
00412   // Use the smear function selected in the constructor
00413   smearFunction->smear( pt, eta, phi, y, parSmear );
00414 
00415   if (debug>9) {
00416     std::cout << "Smearing Pt,eta,phi = " << pt << " " <<  eta << " "
00417          << phi << "; x = ";
00418     for (int i=0; i<SmearType+3; i++) {
00419       std::cout << y[i];
00420     }
00421     std::cout << std::endl;
00422   }
00423 
00424   double ptEtaPhiE[4] = {pt, eta, phi, E};
00425   return( fromPtEtaPhiToPxPyPz(ptEtaPhiE) );
00426 }
00427 
00428 // Biasing function called to modify muon Pt scale at the start.
00429 // -------------------------------------------------------------
00430 lorentzVector MuScleFitUtils::applyBias( const lorentzVector& muon, const int chg )
00431 {
00432   double ptEtaPhiE[4] = {muon.Pt(),muon.Eta(),muon.Phi(),muon.E()};
00433 
00434   if (MuScleFitUtils::debug>1) std::cout << "pt before bias = " << ptEtaPhiE[0] << std::endl;
00435 
00436   // Use functors (although not with the () operator)
00437   // Note that we always pass pt, eta and phi, but internally only the needed
00438   // values are used.
00439   // The functors used are takend from the same group used for the scaling
00440   // thus the name of the method used is "scale".
00441   ptEtaPhiE[0] = biasFunction->scale(ptEtaPhiE[0], ptEtaPhiE[1], ptEtaPhiE[2], chg, MuScleFitUtils::parBias);
00442 
00443   if (MuScleFitUtils::debug>1) std::cout << "pt after bias = " << ptEtaPhiE[0] << std::endl;
00444 
00445   return( fromPtEtaPhiToPxPyPz(ptEtaPhiE) );
00446 }
00447 
00448 // Version of applyScale accepting a std::vector<double> of parameters
00449 // --------------------------------------------------------------
00450 lorentzVector MuScleFitUtils::applyScale (const lorentzVector& muon,
00451                                           const std::vector<double> & parval, const int chg)
00452 {
00453   double * p = new double[(int)(parval.size())];
00454   // Replaced by auto_ptr, which handles delete at the end
00455   // std::auto_ptr<double> p(new double[(int)(parval.size())]);
00456   // Removed auto_ptr, check massResolution for an explanation.
00457   int id = 0;
00458   for (std::vector<double>::const_iterator it=parval.begin(); it!=parval.end(); ++it, ++id) {
00459     //(&*p)[id] = *it;
00460     // Also ok would be (p.get())[id] = *it;
00461     p[id] = *it;
00462   }
00463   lorentzVector tempScaleVec( applyScale (muon, p, chg) );
00464   delete[] p;
00465   return tempScaleVec;
00466 }
00467 
00468 // This is called by the likelihood to "taste" different values for additional corrections
00469 // ---------------------------------------------------------------------------------------
00470 lorentzVector MuScleFitUtils::applyScale (const lorentzVector& muon,
00471                                           double* parval, const int chg)
00472 {
00473   double ptEtaPhiE[4] = {muon.Pt(),muon.Eta(),muon.Phi(),muon.E()};
00474   int shift = parResol.size();
00475 
00476   if (MuScleFitUtils::debug>1) std::cout << "pt before scale = " << ptEtaPhiE[0] << std::endl;
00477 
00478   // the address of parval[shift] is passed as pointer to double. Internally it is used as a normal array, thus:
00479   // array[0] = parval[shift], array[1] = parval[shift+1], ...
00480   ptEtaPhiE[0] = scaleFunction->scale(ptEtaPhiE[0], ptEtaPhiE[1], ptEtaPhiE[2], chg, &(parval[shift]));
00481 
00482   if (MuScleFitUtils::debug>1) std::cout << "pt after scale = " << ptEtaPhiE[0] << std::endl;
00483 
00484   return( fromPtEtaPhiToPxPyPz(ptEtaPhiE) );
00485 }
00486 
00487 // Useful function to convert 4-vector coordinates
00488 // -----------------------------------------------
00489 lorentzVector MuScleFitUtils::fromPtEtaPhiToPxPyPz( const double* ptEtaPhiE )
00490 {
00491   double px = ptEtaPhiE[0]*cos(ptEtaPhiE[2]);
00492   double py = ptEtaPhiE[0]*sin(ptEtaPhiE[2]);
00493   double tmp = 2*atan(exp(-ptEtaPhiE[1]));
00494   double pz = ptEtaPhiE[0]*cos(tmp)/sin(tmp);
00495   double E  = sqrt(px*px+py*py+pz*pz+muMass*muMass);
00496 
00497   // lorentzVector corrMu(px,py,pz,E);
00498   // To fix memory leaks, this is to be substituted with
00499   // std::auto_ptr<lorentzVector> corrMu(new lorentzVector(px, py, pz, E));
00500 
00501   return lorentzVector(px,py,pz,E);
00502 }
00503 
00504 // Dimuon mass
00505 // -----------
00506 inline double MuScleFitUtils::invDimuonMass( const lorentzVector& mu1,
00507                                              const lorentzVector& mu2 )
00508 {
00509   return (mu1+mu2).mass();
00510 }
00511 
00512 // Mass resolution - version accepting a std::vector<double> parval
00513 // -----------------------------------------------------------
00514 double MuScleFitUtils::massResolution( const lorentzVector& mu1,
00515                                        const lorentzVector& mu2,
00516                                        const std::vector<double> & parval )
00517 {
00518   // double * p = new double[(int)(parval.size())];
00519   // Replaced by auto_ptr, which handles delete at the end
00520   // --------- //
00521   // ATTENTION //
00522   // --------- //
00523   // auto_ptr calls delete, not delete[] and thus it must
00524   // not be used with arrays. There are alternatives see
00525   // e.g.: http://www.gotw.ca/gotw/042.htm. The best
00526   // alternative seems to be to switch to vector though.
00527   // std::auto_ptr<double> p(new double[(int)(parval.size())]);
00528 
00529   double * p = new double[(int)(parval.size())];
00530   std::vector<double>::const_iterator it = parval.begin();
00531   int id = 0;
00532   for ( ; it!=parval.end(); ++it, ++id) {
00533     // (&*p)[id] = *it;
00534     p[id] = *it;
00535   }
00536   double massRes = massResolution (mu1, mu2, p);
00537   delete[] p;
00538   return massRes;
00539 }
00540 
00558 double MuScleFitUtils::massResolution( const lorentzVector& mu1,
00559                                        const lorentzVector& mu2,
00560                                        double* parval )
00561 {
00562   double mass   = (mu1+mu2).mass();
00563   double pt1    = mu1.Pt();
00564   double phi1   = mu1.Phi();
00565   double eta1   = mu1.Eta();
00566   double theta1 = 2*atan(exp(-eta1));
00567   double pt2    = mu2.Pt();
00568   double phi2   = mu2.Phi();
00569   double eta2   = mu2.Eta();
00570   double theta2 = 2*atan(exp(-eta2));
00571 
00572   double dmdpt1  = (pt1/std::pow(sin(theta1),2)*sqrt((std::pow(pt2/sin(theta2),2)+mMu2)/(std::pow(pt1/sin(theta1),2)+mMu2))-
00573                     pt2*(cos(phi1-phi2)+cos(theta1)*cos(theta2)/(sin(theta1)*sin(theta2))))/mass;
00574   double dmdpt2  = (pt2/std::pow(sin(theta2),2)*sqrt((std::pow(pt1/sin(theta1),2)+mMu2)/(std::pow(pt2/sin(theta2),2)+mMu2))-
00575                     pt1*(cos(phi2-phi1)+cos(theta2)*cos(theta1)/(sin(theta2)*sin(theta1))))/mass;
00576   double dmdphi1 = pt1*pt2/mass*sin(phi1-phi2);
00577   double dmdphi2 = pt2*pt1/mass*sin(phi2-phi1);
00578   double dmdcotgth1 = (pt1*pt1*cos(theta1)/sin(theta1)*
00579                        sqrt((std::pow(pt2/sin(theta2),2)+mMu2)/(std::pow(pt1/sin(theta1),2)+mMu2)) -
00580                        pt1*pt2*cos(theta2)/sin(theta2))/mass;
00581   double dmdcotgth2 = (pt2*pt2*cos(theta2)/sin(theta2)*
00582                        sqrt((std::pow(pt1/sin(theta1),2)+mMu2)/(std::pow(pt2/sin(theta2),2)+mMu2)) -
00583                        pt2*pt1*cos(theta1)/sin(theta1))/mass;
00584 
00585   if( debugMassResol_ ) {
00586     massResolComponents.dmdpt1 = dmdpt1;
00587     massResolComponents.dmdpt2 = dmdpt2;
00588     massResolComponents.dmdphi1 = dmdphi1;
00589     massResolComponents.dmdphi2 = dmdphi2;
00590     massResolComponents.dmdcotgth1 = dmdcotgth1;
00591     massResolComponents.dmdcotgth2 = dmdcotgth2;
00592   }
00593 
00594   // Resolution parameters:
00595   // ----------------------
00596   double sigma_pt1 = resolutionFunction->sigmaPt( pt1,eta1,parval );
00597   double sigma_pt2 = resolutionFunction->sigmaPt( pt2,eta2,parval );
00598   double sigma_phi1 = resolutionFunction->sigmaPhi( pt1,eta1,parval );
00599   double sigma_phi2 = resolutionFunction->sigmaPhi( pt2,eta2,parval );
00600   double sigma_cotgth1 = resolutionFunction->sigmaCotgTh( pt1,eta1,parval );
00601   double sigma_cotgth2 = resolutionFunction->sigmaCotgTh( pt2,eta2,parval );
00602   double cov_pt1pt2 = resolutionFunction->covPt1Pt2( pt1, eta1, pt2, eta2, parval );
00603 
00604   // Sigma_Pt is defined as a relative sigmaPt/Pt for this reason we need to
00605   // multiply it by pt.
00606   double mass_res = sqrt(std::pow(dmdpt1*sigma_pt1*pt1,2)+std::pow(dmdpt2*sigma_pt2*pt2,2)+
00607                          std::pow(dmdphi1*sigma_phi1,2)+std::pow(dmdphi2*sigma_phi2,2)+
00608                          std::pow(dmdcotgth1*sigma_cotgth1,2)+std::pow(dmdcotgth2*sigma_cotgth2,2)+
00609                          2*dmdpt1*dmdpt2*cov_pt1pt2*sigma_pt1*sigma_pt2);
00610 
00611   if (debug>19) {
00612     std::cout << "  Pt1=" << pt1 << " phi1=" << phi1 << " cotgth1=" << cos(theta1)/sin(theta1) << " - Pt2=" << pt2
00613          << " phi2=" << phi2 << " cotgth2=" << cos(theta2)/sin(theta2) << std::endl;
00614     std::cout << " P[0]="
00615          << parval[0] << " P[1]=" << parval[1] << "P[2]=" << parval[2] << " P[3]=" << parval[3] << std::endl;
00616     std::cout << "  Dmdpt1= " << dmdpt1 << " dmdpt2= " << dmdpt2 << " sigma_pt1="
00617          << sigma_pt1 << " sigma_pt2=" << sigma_pt2 << std::endl;
00618     std::cout << "  Dmdphi1= " << dmdphi1 << " dmdphi2= " << dmdphi2 << " sigma_phi1="
00619          << sigma_phi1 << " sigma_phi2=" << sigma_phi2 << std::endl;
00620     std::cout << "  Dmdcotgth1= " << dmdcotgth1 << " dmdcotgth2= " << dmdcotgth2
00621          << " sigma_cotgth1="
00622          << sigma_cotgth1 << " sigma_cotgth2=" << sigma_cotgth2 << std::endl;
00623     std::cout << "  Mass resolution (pval) for muons of Pt = " << pt1 << " " << pt2
00624          << " : " << mass << " +- " << mass_res << std::endl;
00625   }
00626 
00627   // Debug std::cout
00628   // ----------
00629   bool didit = false;
00630   for (int ires=0; ires<6; ires++) {
00631     if (!didit && resfind[ires]>0 && fabs(mass-ResMass[ires])<ResHalfWidth[ires]) {
00632       if (mass_res>ResMaxSigma[ires] && counter_resprob<100) {
00633         counter_resprob++;
00634         LogDebug("MuScleFitUtils") << "RESOLUTION PROBLEM: ires=" << ires << std::endl;
00635         didit = true;
00636       }
00637     }
00638   }
00639 
00640   return mass_res;
00641 }
00642 
00647 double MuScleFitUtils::massResolution( const lorentzVector& mu1,
00648                                        const lorentzVector& mu2,
00649                                        const ResolutionFunction & resolFunc )
00650 {
00651   double mass   = (mu1+mu2).mass();
00652   double pt1    = mu1.Pt();
00653   double phi1   = mu1.Phi();
00654   double eta1   = mu1.Eta();
00655   double theta1 = 2*atan(exp(-eta1));
00656   double pt2    = mu2.Pt();
00657   double phi2   = mu2.Phi();
00658   double eta2   = mu2.Eta();
00659   double theta2 = 2*atan(exp(-eta2));
00660 
00661   double dmdpt1  = (pt1/std::pow(sin(theta1),2)*sqrt((std::pow(pt2/sin(theta2),2)+mMu2)/(std::pow(pt1/sin(theta1),2)+mMu2))-
00662                     pt2*(cos(phi1-phi2)+cos(theta1)*cos(theta2)/(sin(theta1)*sin(theta2))))/mass;
00663   double dmdpt2  = (pt2/std::pow(sin(theta2),2)*sqrt((std::pow(pt1/sin(theta1),2)+mMu2)/(std::pow(pt2/sin(theta2),2)+mMu2))-
00664                     pt1*(cos(phi2-phi1)+cos(theta2)*cos(theta1)/(sin(theta2)*sin(theta1))))/mass;
00665   double dmdphi1 = pt1*pt2/mass*sin(phi1-phi2);
00666   double dmdphi2 = pt2*pt1/mass*sin(phi2-phi1);
00667   double dmdcotgth1 = (pt1*pt1*cos(theta1)/sin(theta1)*
00668                        sqrt((std::pow(pt2/sin(theta2),2)+mMu2)/(std::pow(pt1/sin(theta1),2)+mMu2)) -
00669                        pt1*pt2*cos(theta2)/sin(theta2))/mass;
00670   double dmdcotgth2 = (pt2*pt2*cos(theta2)/sin(theta2)*
00671                        sqrt((std::pow(pt1/sin(theta1),2)+mMu2)/(std::pow(pt2/sin(theta2),2)+mMu2)) -
00672                        pt2*pt1*cos(theta1)/sin(theta1))/mass;
00673 
00674   // Resolution parameters:
00675   // ----------------------
00676   double sigma_pt1 = resolFunc.sigmaPt( mu1 );
00677   double sigma_pt2 = resolFunc.sigmaPt( mu2 );
00678   double sigma_phi1 = resolFunc.sigmaPhi( mu1 );
00679   double sigma_phi2 = resolFunc.sigmaPhi( mu2 );
00680   double sigma_cotgth1 = resolFunc.sigmaCotgTh( mu1 );
00681   double sigma_cotgth2 = resolFunc.sigmaCotgTh( mu2 );
00682 
00683   // Sigma_Pt is defined as a relative sigmaPt/Pt for this reason we need to
00684   // multiply it by pt.
00685   double mass_res = sqrt(std::pow(dmdpt1*sigma_pt1*pt1,2)+std::pow(dmdpt2*sigma_pt2*pt2,2)+
00686                          std::pow(dmdphi1*sigma_phi1,2)+std::pow(dmdphi2*sigma_phi2,2)+
00687                          std::pow(dmdcotgth1*sigma_cotgth1,2)+std::pow(dmdcotgth2*sigma_cotgth2,2));
00688 
00689   return mass_res;
00690 }
00691 
00692 
00693 // Mass probability - version with linear background included, accepts std::vector<double> parval
00694 // -----------------------------------------------------------------------------------------
00695 double MuScleFitUtils::massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, const std::vector<double> & parval, const bool doUseBkgrWindow, const double & eta1, const double & eta2 )
00696 {
00697 #ifdef USE_CALLGRIND
00698   CALLGRIND_START_INSTRUMENTATION;
00699 #endif
00700 
00701   double * p = new double[(int)(parval.size())];
00702   // Replaced by auto_ptr, which handles delete at the end
00703   // Removed auto_ptr, check massResolution for an explanation.
00704   // std::auto_ptr<double> p(new double[(int)(parval.size())]);
00705 
00706   std::vector<double>::const_iterator it = parval.begin();
00707   int id = 0;
00708   for ( ; it!=parval.end(); ++it, ++id) {
00709     // (&*p)[id] = *it;
00710     p[id] = *it;
00711   }
00712   // p must be passed by value as below:
00713   double massProbability = massProb( mass, resEta, rapidity, massResol, p, doUseBkgrWindow, eta1, eta2 );
00714   delete[] p;
00715 
00716 #ifdef USE_CALLGRIND
00717   CALLGRIND_STOP_INSTRUMENTATION;
00718   CALLGRIND_DUMP_STATS;
00719 #endif
00720 
00721   return massProbability;
00722 }
00723 
00729 double MuScleFitUtils::probability( const double & mass, const double & massResol,
00730                                     const double GLvalue[][1001][1001], const double GLnorm[][1001],
00731                                     const int iRes, const int iY )
00732 {
00733   if( iRes == 0 && iY > 23 ) {
00734     std::cout << "WARNING: rapidity bin selected = " << iY << " but there are only histograms for the first 24 bins" << std::endl;
00735   }
00736 
00737   double PS = 0.;
00738   bool insideProbMassWindow = true;
00739   // Interpolate the four values of GLZValue[] in the
00740   // grid square within which the (mass,sigma) values lay
00741   // ----------------------------------------------------
00742   // This must be done with respect to the width used in the computation of the probability distribution,
00743   // so that the bin 0 really matches the bin 0 of that distribution.
00744   //  double fracMass = (mass-(ResMass[iRes]-ResHalfWidth[iRes]))/(2*ResHalfWidth[iRes]);
00745   double fracMass = (mass - ResMinMass[iRes])/(2*ResHalfWidth[iRes]);
00746   if (debug>1) std::cout << std::setprecision(9)<<"mass ResMinMass[iRes] ResHalfWidth[iRes] ResHalfWidth[iRes]"
00747                     << mass << " "<<ResMinMass[iRes]<<" "<<ResHalfWidth[iRes]<<" "<<ResHalfWidth[iRes]<<std::endl;
00748   int iMassLeft  = (int)(fracMass*(double)nbins);
00749   int iMassRight = iMassLeft+1;
00750   double fracMassStep = (double)nbins*(fracMass - (double)iMassLeft/(double)nbins);
00751   if (debug>1) std::cout<<"nbins iMassLeft fracMass "<<nbins<<" "<<iMassLeft<<" "<<fracMass<<std::endl;
00752 
00753   // Simple protections for the time being: the region where we fit should not include
00754   // values outside the boundaries set by ResMass-ResHalfWidth : ResMass+ResHalfWidth
00755   // ---------------------------------------------------------------------------------
00756   if (iMassLeft<0) {
00757     edm::LogInfo("probability") << "WARNING: fracMass=" << fracMass << ", iMassLeft="
00758                            << iMassLeft << "; mass = " << mass << " and bounds are " << ResMinMass[iRes]
00759                            << ":" << ResMinMass[iRes]+2*ResHalfWidth[iRes] << " - iMassLeft set to 0" << std::endl;
00760     iMassLeft  = 0;
00761     iMassRight = 1;
00762     insideProbMassWindow = false;
00763   }
00764   if (iMassRight>nbins) {
00765     edm::LogInfo("probability") << "WARNING: fracMass=" << fracMass << ", iMassRight="
00766                            << iMassRight << "; mass = " << mass << " and bounds are " << ResMinMass[iRes]
00767                            << ":" << ResMass[iRes]+2*ResHalfWidth[iRes] << " - iMassRight set to " << nbins-1 << std::endl;
00768     iMassLeft  = nbins-1;
00769     iMassRight = nbins;
00770     insideProbMassWindow = false;
00771   }
00772   double fracSigma = (massResol/ResMaxSigma[iRes]);
00773   int iSigmaLeft = (int)(fracSigma*(double)nbins);
00774   int iSigmaRight = iSigmaLeft+1;
00775   double fracSigmaStep = (double)nbins * (fracSigma - (double)iSigmaLeft/(double)nbins);
00776   // std::cout << "massResol = " << massResol << std::endl;
00777   // std::cout << "ResMaxSigma["<<iRes<<"] = " << ResMaxSigma[iRes] << std::endl;
00778   // std::cout << "fracSigma = " << fracSigma << std::endl;
00779   // std::cout << "nbins = " << nbins << std::endl;
00780   // std::cout << "ISIGMALEFT = " << iSigmaLeft << std::endl;
00781   // std::cout << "ISIGMARIGHT = " << iSigmaRight << std::endl;
00782   // std::cout << "fracSigmaStep = " << fracSigmaStep << std::endl;
00783 
00784   // Simple protections for the time being: they should not affect convergence, since
00785   // ResMaxSigma is set to very large values, and if massResol exceeds them the fit
00786   // should not get any prize for that (for large sigma, the prob. distr. becomes flat)
00787   // ----------------------------------------------------------------------------------
00788   if (iSigmaLeft<0) {
00789     edm::LogInfo("probability") << "WARNING: fracSigma = " << fracSigma << ", iSigmaLeft="
00790                            << iSigmaLeft << ", with massResol = " << massResol << " and ResMaxSigma[iRes] = "
00791                            << ResMaxSigma[iRes] << " -  iSigmaLeft set to 0" << std::endl;
00792     iSigmaLeft  = 0;
00793     iSigmaRight = 1;
00794   }
00795   if (iSigmaRight>nbins ) {
00796     if (counter_resprob<100)
00797       edm::LogInfo("probability") << "WARNING: fracSigma = " << fracSigma << ", iSigmaRight="
00798                              << iSigmaRight << ", with massResol = " << massResol << " and ResMaxSigma[iRes] = "
00799                              << ResMaxSigma[iRes] << " -  iSigmaRight set to " << nbins-1 << std::endl;
00800     iSigmaLeft  = nbins-1;
00801     iSigmaRight = nbins;
00802   }
00803 
00804   // If f11,f12,f21,f22 are the values at the four corners, one finds by linear interpolation the
00805   // formula below for PS
00806   // --------------------------------------------------------------------------------------------
00807   if( insideProbMassWindow ) {
00808     double f11 = 0.;
00809     if (GLnorm[iY][iSigmaLeft]>0)
00810       f11 = GLvalue[iY][iMassLeft][iSigmaLeft] / GLnorm[iY][iSigmaLeft];
00811     double f12 = 0.;
00812     if (GLnorm[iY][iSigmaRight]>0)
00813       f12 = GLvalue[iY][iMassLeft][iSigmaRight] / GLnorm[iY][iSigmaRight];
00814     double f21 = 0.;
00815     if (GLnorm[iY][iSigmaLeft]>0)
00816       f21 = GLvalue[iY][iMassRight][iSigmaLeft] / GLnorm[iY][iSigmaLeft];
00817     double f22 = 0.;
00818     if (GLnorm[iY][iSigmaRight]>0)
00819       f22 = GLvalue[iY][iMassRight][iSigmaRight] / GLnorm[iY][iSigmaRight];
00820     PS = f11 + (f12-f11)*fracSigmaStep + (f21-f11)*fracMassStep +
00821       (f22-f21-f12+f11)*fracMassStep*fracSigmaStep;
00822     if (PS>0.1 || debug>1) LogDebug("MuScleFitUtils") << "iRes = " << iRes << " PS=" << PS << " f11,f12,f21,f22="
00823                                                       << f11 << " " << f12 << " " << f21 << " " << f22 << " "
00824                                                       << " fSS=" << fracSigmaStep << " fMS=" << fracMassStep << " iSL, iSR="
00825                                                       << iSigmaLeft << " " << iSigmaRight
00826                                                       << " GLvalue["<<iY<<"]["<<iMassLeft<<"] = " << GLvalue[iY][iMassLeft][iSigmaLeft]
00827                                                       << " GLnorm["<<iY<<"]["<<iSigmaLeft<<"] = " << GLnorm[iY][iSigmaLeft] << std::endl;
00828 
00829 //     if (PS>0.1) std::cout << "iRes = " << iRes << " PS=" << PS << " f11,f12,f21,f22="
00830 //                      << f11 << " " << f12 << " " << f21 << " " << f22 << " "
00831 //                      << " fSS=" << fracSigmaStep << " fMS=" << fracMassStep << " iSL, iSR="
00832 //                      << iSigmaLeft << " " << iSigmaRight << " GLV,GLN="
00833 //                      << GLvalue[iY][iMassLeft][iSigmaLeft]
00834 //                      << " " << GLnorm[iY][iSigmaLeft] << std::endl;
00835 
00836   }
00837   else {
00838     edm::LogInfo("probability") << "outside mass probability window. Setting PS["<<iRes<<"] = 0" << std::endl;
00839   }
00840 
00841 //   if( PS != PS ) {
00842 //     std::cout << "ERROR: PS = " << PS << " for iRes = " << iRes << std::endl;
00843 
00844 //     std::cout << "mass = " << mass << ", massResol = " << massResol << std::endl;
00845 //     std::cout << "fracMass = " << fracMass << ", iMassLeft = " << iMassLeft
00846 //          << ", iMassRight = " << iMassRight << ", fracMassStep = " << fracMassStep << std::endl;
00847 //     std::cout << "fracSigma = " << fracSigma << ", iSigmaLeft = " << iSigmaLeft
00848 //          << ", iSigmaRight = " << iSigmaRight << ", fracSigmaStep = " << fracSigmaStep << std::endl;
00849 //     std::cout << "ResMaxSigma["<<iRes<<"] = " << ResMaxSigma[iRes] << std::endl;
00850 
00851 //     std::cout << "GLvalue["<<iY<<"]["<<iMassLeft<<"] = " << GLvalue[iY][iMassLeft][iSigmaLeft]
00852 //          << " GLnorm["<<iY<<"]["<<iSigmaLeft<<"] = " << GLnorm[iY][iSigmaLeft] << std::endl;
00853 //   }
00854 
00855   return PS;
00856 }
00857 
00858 // Mass probability - version with linear background included
00859 // ----------------------------------------------------------
00860 double MuScleFitUtils::massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, double * parval, const bool doUseBkgrWindow, const double & eta1, const double & eta2 ) {
00861 
00862   // This routine computes the likelihood that a given measured mass "measMass" is
00863   // the result of a reference mass ResMass[] if the resolution
00864   // expected for the two muons is massResol.
00865   // This version includes two parameters (the last two in parval, by default)
00866   // to size up the background fraction and its relative normalization with respect
00867   // to the signal shape.
00868   //
00869   // We model the signal probability with a Lorentz L(M,H) of resonance mass M and natural width H
00870   // convoluted with a gaussian G(m,s) of measured mass m and expected mass resolution s,
00871   // by integrating over the intersection of the supports of L and G (which can be made to coincide with
00872   // the region where L is non-zero, given that H<<s in most cases) the product L(M,H)*G(m,s) dx as follows:
00873   //
00874   //   GL(m,s) = Int(M-10H,M+10H) [ L(x-M,H) * G(x-m,s) ] dx
00875   //
00876   // The above convolution is computed numerically by an independent root macro, Probs.C, which outputs
00877   // the values in six 1001x1001 grids, one per resonance.
00878   //
00879   // NB THe following block of explanations for background models is outdated, see detailed
00880   // explanations where the code computes PB.
00881   // +++++++++++++++++++++++
00882   //   For the background, instead, we have two choices: a linear and an exponential model.
00883   //     * For the linear model, we choose a one-parameter form whereby the line is automatically normalized
00884   //       in the support [x1,x2] where we defined our "signal region", as follows:
00885   //
00886   //         B(x;b) = 1/(x2-x1) + {x - (x2+x1)/2} * b
00887   //
00888   //       Defined as above, B(x) is a line passing for the point of coordinates (x1+x2)/2, 1/(x2-x1),
00889   //       whose slope b has as support the interval ( -2/[(x1-x2)*(x1+x2)], 2/[(x1-x2)*(x1+x2)] )
00890   //       so that B(x) is always positive-definite in [x1,x2].
00891   //
00892   //     * For the exponential model, we define B(x;b) as
00893   //
00894   //         B(x;b) = b * { exp(-b*x) / [exp(-b*x1)-exp(-b*x2)] }
00895   //
00896   //       This one-parameter definition is automatically normalized to unity in [x1,x2], with a parameter
00897   //       b which has to be positive in order for the slope to be negative.
00898   //       Please note that this model is not useful in most circumstances; a more useful form would be one
00899   //       which included a linear component.
00900   // ++++++++++++++++++++++
00901   //
00902   // Once GL(m,s) and B(x;b) are defined, we introduce a further parameter a, such that we can have the
00903   // likelihood control the relative fraction of signal and background. We first normalize GL(m,s) for
00904   // any given s by taking the integral
00905   //
00906   //   Int(x1,x2) GL(m,s) dm = K_s
00907   //
00908   // We then define the probability as
00909   //
00910   //   P(m,s,a,b) = GL(m,s)/K_s * a  +  B(x,b) * (1-a)
00911   //
00912   // with a taking on values in the interval [0,1].
00913   // Defined as above, the probability is well-behaved, in the sense that it has a value between 0 and 1,
00914   // and the four parameters m,s,a,b fully control its shape.
00915   //
00916   // It is to be noted that the formulation above requires the computation of two rather time-consuming
00917   // integrals. The one defining GL(m,s) can be stored in a TH2D and loaded by the constructor from a
00918   // file suitably prepared, and this will save loads of computing time.
00919   // ----------------------------------------------------------------------------------------------------
00920 
00921   double P = 0.;
00922   int crossSectionParShift = parResol.size() + parScale.size();
00923   // Take the relative cross sections
00924   std::vector<double> relativeCrossSections = crossSectionHandler->relativeCrossSections(&(parval[crossSectionParShift]), resfind);
00925   //   for( unsigned int i=0; i<relativeCrossSections.size(); ++i ) {
00926   //     std::cout << "relativeCrossSections["<<i<<"] = " << relativeCrossSections[i] << std::endl;
00927   //     std::cout << "parval["<<crossSectionParShift+i<<"] = " << parval[crossSectionParShift+i] << std::endl;
00928   //   }
00929 
00930   // int bgrParShift = crossSectionParShift + parCrossSection.size();
00931   int bgrParShift = crossSectionParShift + crossSectionHandler->parNum();
00932   double Bgrp1 = 0.;
00933   //   double Bgrp2 = 0.;
00934   //   double Bgrp3 = 0.;
00935 
00936   // NB defined as below, P is a non-rigorous "probability" that we observe
00937   // a dimuon mass "mass", given ResMass[], gamma, and massResol. It is what we need for the
00938   // fit which finds the best resolution parameters, though. A definition which is
00939   // more properly a probability is given below (in massProb2()), where the result
00940   // cannot be used to fit resolution parameters because the fit would always prefer
00941   // to set the res parameters to the minimum possible value (best resolution),
00942   // to have a probability close to one of observing any mass.
00943   // -------------------------------------------------------------------------------
00944 
00945   // Determine what resonance(s) we have to deal with
00946   // NB for now we assume equal xs for each resonance
00947   // so we do not assign them different weights
00948   // ------------------------------------------------
00949   double PS[6] = {0.};
00950   double PB = 0.;
00951   double PStot[6] = {0.};
00952 
00953   // Should be removed because it is not used
00954   bool resConsidered[6] = {false};
00955 
00956   bool useBackgroundWindow = (doBackgroundFit[loopCounter] || doUseBkgrWindow);
00957   // bool useBackgroundWindow = (doBackgroundFit[loopCounter]);
00958 
00959   // First check the Z, which is divided in 24 rapidity bins
00960   // NB max value of Z rapidity to be considered is 2.4 here
00961   // -------------------------------------------------------
00962 
00963   // Do this only if we want to use the rapidity bins for the Z
00964   if( MuScleFitUtils::rapidityBinsForZ_ ) {
00965     // ATTENTION: cut on Z rapidity at 2.4 since we only have histograms up to that value
00966     // std::pair<double, double> windowFactors = backgroundHandler->windowFactors( useBackgroundWindow, 0 );
00967     std::pair<double, double> windowBorders = backgroundHandler->windowBorders( useBackgroundWindow, 0 );
00968     if( resfind[0]>0
00969         // && checkMassWindow( mass, 0,
00970         //                     backgroundHandler->resMass( useBackgroundWindow, 0 ),
00971         //                     windowFactors.first, windowFactors.second )
00972         && checkMassWindow( mass, windowBorders.first, windowBorders.second )
00973         // && fabs(rapidity)<2.4
00974         ) {
00975 
00976       int iY = (int)(fabs(rapidity)*10.);
00977       if( iY > 23 ) iY = 23;
00978 
00979       if (MuScleFitUtils::debug>1) std::cout << "massProb:resFound = 0, rapidity bin =" << iY << std::endl;
00980 
00981       // In this case the last value is the rapidity bin
00982       PS[0] = probability(mass, massResol, GLZValue, GLZNorm, 0, iY);
00983 
00984       if( PS[0] != PS[0] ) {
00985         std::cout << "ERROR: PS[0] = nan, setting it to 0" << std::endl;
00986         PS[0] = 0;
00987       }
00988 
00989       // std::pair<double, double> bgrResult = backgroundHandler->backgroundFunction( doBackgroundFit[loopCounter],
00990       //                                                                                   &(parval[bgrParShift]), MuScleFitUtils::totalResNum, 0,
00991       //                                                                                   resConsidered, ResMass, ResHalfWidth, MuonType, mass, resEta );
00992 
00993       std::pair<double, double> bgrResult = backgroundHandler->backgroundFunction( doBackgroundFit[loopCounter],
00994                                                                                    &(parval[bgrParShift]), MuScleFitUtils::totalResNum, 0,
00995                                                                                    resConsidered, ResMass, ResHalfWidth, MuonType, mass, eta1, eta2 );
00996 
00997       Bgrp1 = bgrResult.first;
00998       // When fitting the background we have only one Bgrp1
00999       // When not fitting the background we have many only in a superposition region and this case is treated
01000       // separately after this loop
01001       PB = bgrResult.second;
01002       if( PB != PB ) PB = 0;
01003       PStot[0] = (1-Bgrp1)*PS[0] + Bgrp1*PB;
01004 
01005       // PStot[0] *= crossSectionHandler->crossSection(0);
01006       // PStot[0] *= parval[crossSectionParShift];
01007       PStot[0] *= relativeCrossSections[0];
01008       // std::cout << "PStot["<<0<<"] = " << "(1-"<<Bgrp1<<")*"<<PS[0]<<" + "<<Bgrp1<<"*"<<PB<<" = " << PStot[0] << std::endl;
01009     }
01010     else {
01011       if( debug > 0 ) {
01012         std::cout << "Mass = " << mass << " outside range with rapidity = " << rapidity << std::endl;
01013         std::cout << "with resMass = " << backgroundHandler->resMass( useBackgroundWindow, 0 ) << " and left border = " << windowBorders.first << " right border = " << windowBorders.second << std::endl;
01014       }
01015     }
01016   }
01017   // Next check the other resonances
01018   // -------------------------------
01019   int firstRes = 1;
01020   if( !MuScleFitUtils::rapidityBinsForZ_ ) firstRes = 0;
01021   for( int ires=firstRes; ires<6; ++ires ) {
01022     if( resfind[ires] > 0 ) {
01023       // First is left, second is right (returns (1,1) in the case of resonances, it could be improved avoiding the call in this case)
01024       // std::pair<double, double> windowFactor = backgroundHandler->windowFactors( useBackgroundWindow, ires );
01025       std::pair<double, double> windowBorder = backgroundHandler->windowBorders( useBackgroundWindow, ires );
01026       if( checkMassWindow(mass, windowBorder.first, windowBorder.second) ) {
01027         if (MuScleFitUtils::debug>1) std::cout << "massProb:resFound = " << ires << std::endl;
01028 
01029         // In this case the rapidity value is instead the resonance index again.
01030         PS[ires] = probability(mass, massResol, GLValue, GLNorm, ires, ires);
01031 
01032         std::pair<double, double> bgrResult = backgroundHandler->backgroundFunction( doBackgroundFit[loopCounter],
01033                                                                                      &(parval[bgrParShift]), MuScleFitUtils::totalResNum, ires,
01034                                                                                      // resConsidered, ResMass, ResHalfWidth, MuonType, mass, resEta );
01035                                                                                      resConsidered, ResMass, ResHalfWidth, MuonType, mass, eta1, eta2 );
01036         Bgrp1 = bgrResult.first;
01037         PB = bgrResult.second;
01038 
01039         if( PB != PB ) PB = 0;
01040         PStot[ires] = (1-Bgrp1)*PS[ires] + Bgrp1*PB;
01041         if( MuScleFitUtils::debug>0 ) std::cout << "PStot["<<ires<<"] = " << "(1-"<<Bgrp1<<")*"<<PS[ires]<<" + "<<Bgrp1<<"*"<<PB<<" = " << PStot[ires] << std::endl;
01042 
01043         PStot[ires] *= relativeCrossSections[ires];
01044       }
01045     }
01046   }
01047 
01048   for( int i=0; i<6; ++i ) {
01049     P += PStot[i];
01050   }
01051 
01052   if( MuScleFitUtils::signalProb_ != 0 && MuScleFitUtils::backgroundProb_ != 0 ) {
01053     double PStotTemp = 0.;
01054     for( int i=0; i<6; ++i ) {
01055       PStotTemp += PS[i]*relativeCrossSections[i];
01056     }
01057     if( PStotTemp != PStotTemp ) {
01058       std::cout << "ERROR: PStotTemp = nan!!!!!!!!!" << std::endl;
01059       int parnumber = (int)(parResol.size()+parScale.size()+crossSectionHandler->parNum()+parBgr.size());
01060       for( int i=0; i<6; ++i ) {
01061         std::cout << "PS[i] = " << PS[i] << std::endl;
01062         if( PS[i] != PS[i] ) {
01063           std::cout << "mass = " << mass << std::endl;
01064           std::cout << "massResol = " << massResol << std::endl;
01065           for( int j=0; j<parnumber; ++j ) {
01066             std::cout << "parval["<<j<<"] = " << parval[j] << std::endl;
01067           }
01068         }
01069       }
01070     }
01071     if( PStotTemp == PStotTemp ) {
01072       MuScleFitUtils::signalProb_->SetBinContent(MuScleFitUtils::minuitLoop_, MuScleFitUtils::signalProb_->GetBinContent(MuScleFitUtils::minuitLoop_) + PStotTemp);
01073     }
01074     if (debug>0) std::cout << "mass = " << mass << ", P = " << P << ", PStot = " << PStotTemp << ", PB = " << PB << ", bgrp1 = " << Bgrp1 << std::endl;
01075 
01076     MuScleFitUtils::backgroundProb_->SetBinContent(MuScleFitUtils::minuitLoop_, MuScleFitUtils::backgroundProb_->GetBinContent(MuScleFitUtils::minuitLoop_) + PB);
01077   }
01078   return P;
01079 }
01080 
01081 // Method to check if the mass value is within the mass window of the i-th resonance.
01082 // inline bool MuScleFitUtils::checkMassWindow( const double & mass, const int ires, const double & resMass, const double & leftFactor, const double & rightFactor )
01083 // {
01084 //   return( mass-resMass > -leftFactor*massWindowHalfWidth[MuonTypeForCheckMassWindow][ires]
01085 //           && mass-resMass < rightFactor*massWindowHalfWidth[MuonTypeForCheckMassWindow][ires] );
01086 // }
01087 inline bool MuScleFitUtils::checkMassWindow( const double & mass, const double & leftBorder, const double & rightBorder )
01088 {
01089   return( (mass > leftBorder) && (mass < rightBorder) );
01090 }
01091 
01092 // Function that returns the weight for a muon pair
01093 // ------------------------------------------------
01094 double MuScleFitUtils::computeWeight( const double & mass, const int iev, const bool doUseBkgrWindow )
01095 {
01096   // Compute weight for this event
01097   // -----------------------------
01098   double weight = 0.;
01099 
01100   // Take the highest-mass resonance within bounds
01101   // NB this must be revised once credible estimates of the relative xs of Y(1S), (2S), and (3S)
01102   // are made. Those are priors in the decision of which resonance to assign to an in-between event.
01103   // -----------------------------------------------------------------------------------------------
01104 
01105   if( doUseBkgrWindow && (debug > 0) ) std::cout << "using backgrond window for mass = " << mass << std::endl;
01106   // bool useBackgroundWindow = (doBackgroundFit[loopCounter] || doUseBkgrWindow);
01107   bool useBackgroundWindow = (doBackgroundFit[loopCounter]);
01108 
01109   for (int ires=0; ires<6; ires++) {
01110     if (resfind[ires]>0 && weight==0.) {
01111       // std::pair<double, double> windowFactor = backgroundHandler->windowFactors( useBackgroundWindow, ires );
01112       std::pair<double, double> windowBorder = backgroundHandler->windowBorders( useBackgroundWindow, ires );
01113       // if( checkMassWindow(mass, ires, backgroundHandler->resMass( useBackgroundWindow, ires ),
01114       //                     windowFactor.first, windowFactor.second) ) {
01115       if( checkMassWindow(mass, windowBorder.first, windowBorder.second) ) {
01116         weight = 1.0;
01117         if( doUseBkgrWindow && (debug > 0) ) std::cout << "setting weight to = " << weight << std::endl;
01118       }
01119     }
01120   }
01121 
01122   return weight;
01123 }
01124 
01125 // Likelihood minimization routine
01126 // -------------------------------
01127 void MuScleFitUtils::minimizeLikelihood()
01128 {
01129   // Output file with fit parameters resulting from minimization
01130   // -----------------------------------------------------------
01131   ofstream FitParametersFile;
01132   FitParametersFile.open ("FitParameters.txt", std::ios::app);
01133   FitParametersFile << "Fitting with resolution, scale, bgr function # "
01134                     << ResolFitType << " " << ScaleFitType << " " << BgrFitType
01135                     << " - Iteration " << loopCounter << std::endl;
01136 
01137   // Fill parvalue and other vectors needed for the fitting
01138   // ------------------------------------------------------
01139 
01140 
01141 
01142 
01143 
01144   // ----- //
01145   // FIXME //
01146   // ----- //
01147   // this was changed to verify the possibility that fixed parameters influence the errors.
01148   // It must be 0 otherwise the parameters for resonances will not be passed by minuit (will be always 0).
01149   // Should be removed.
01150   int parForResonanceWindows = 0;
01151   // int parnumber = (int)(parResol.size()+parScale.size()+parCrossSection.size()+parBgr.size() - parForResonanceWindows);
01152   int parnumber = (int)(parResol.size()+parScale.size()+crossSectionHandler->parNum()+parBgr.size() - parForResonanceWindows);
01153 
01154 
01155 
01156 
01157 
01158 
01159   int parnumberAll = (int)(parResol.size()+parScale.size()+crossSectionHandler->parNum()+parBgr.size());
01160 
01161   // parvalue is a std::vector<std::vector<double> > storing all the parameters from all the loops
01162   parvalue.push_back(parResol);
01163   std::vector<double> *tmpVec = &(parvalue.back());
01164 
01165   // If this is not the first loop we want to start from neutral values
01166   // Otherwise the scale will start with values correcting again a bias
01167   // that is already corrected.
01168   if( scaleFitNotDone_ ) {
01169     tmpVec->insert( tmpVec->end(), parScale.begin(), parScale.end() );
01170     std::cout << "scaleFitNotDone: tmpVec->size() = " << tmpVec->size() << std::endl;
01171   }
01172   else {
01173     scaleFunction->resetParameters(tmpVec);
01174     std::cout << "scaleFitDone: tmpVec->size() = " << tmpVec->size() << std::endl;
01175   }
01176   tmpVec->insert( tmpVec->end(), parCrossSection.begin(), parCrossSection.end() );
01177   tmpVec->insert( tmpVec->end(), parBgr.begin(), parBgr.end() );
01178   int i = 0;
01179   std::vector<double>::const_iterator it = tmpVec->begin();
01180   for( ; it != tmpVec->end(); ++it, ++i ) {
01181     std::cout << "tmpVec["<<i<<"] = " << *it << std::endl;
01182   }
01183 
01184   // Empty vector of size = number of cross section fitted parameters. Note that the cross section
01185   // fit works in a different way than the others and it uses ratios of the paramters passed via cfg.
01186   // We use this empty vector for compatibility with the rest of the structure.
01187   std::vector<int> crossSectionParNumSizeVec( MuScleFitUtils::crossSectionHandler->parNum(), 0 );
01188 
01189   std::vector<int> parfix(parResolFix);
01190   parfix.insert( parfix.end(), parScaleFix.begin(), parScaleFix.end() );
01191   parfix.insert( parfix.end(), crossSectionParNumSizeVec.begin(), crossSectionParNumSizeVec.end() );
01192   parfix.insert( parfix.end(), parBgrFix.begin(), parBgrFix.end() );
01193 
01194   std::vector<int> parorder(parResolOrder);
01195   parorder.insert( parorder.end(), parScaleOrder.begin(), parScaleOrder.end() );
01196   parorder.insert( parorder.end(), crossSectionParNumSizeVec.begin(), crossSectionParNumSizeVec.end() );
01197   parorder.insert( parorder.end(), parBgrOrder.begin(), parBgrOrder.end() );
01198 
01199   // This is filled later
01200   std::vector<double> parerr(3*parnumberAll,0.);
01201 
01202   if (debug>19) {
01203     std::cout << "[MuScleFitUtils-minimizeLikelihood]: Parameters before likelihood " << std::endl;
01204     for (unsigned int i=0; i<(unsigned int)parnumberAll; i++) {
01205       std::cout << "  Par # " << i << " = " << parvalue[loopCounter][i] << " : free = " << parfix[i] << "; order = "
01206            << parorder[i] << std::endl;
01207     }
01208   }
01209 
01210 //   // Background rescaling from regions to resonances
01211 //   // -----------------------------------------------
01212 //   // If we are in a loop > 0 and we are not fitting the background, but we have fitted it in the previous iteration
01213 //   if( loopCounter > 0 && !(doBackgroundFit[loopCounter]) && doBackgroundFit[loopCounter-1] ) {
01214 //     // This rescales from regions to resonances
01215 //     int localMuonType = MuonType;
01216 //     if( MuonType > 2 ) localMuonType = 2;
01217 //     backgroundHandler->rescale( parBgr, ResMass, massWindowHalfWidth[localMuonType],
01218 //                                 MuScleFitUtils::SavedPair);
01219 //   }
01220 
01221   // Init Minuit
01222   // -----------
01223   TMinuit rmin (parnumber);
01224   rminPtr_ = &rmin;
01225   rmin.SetFCN (likelihood);     // Unbinned likelihood
01226   // Standard initialization of minuit parameters:
01227   // sets input to be $stdin, output to be $stdout
01228   // and saving to a file.
01229   rmin.mninit (5, 6, 7);
01230   int ierror = 0;
01231   int istat;
01232   double arglis[4];
01233   arglis[0] = FitStrategy;      // Strategy 1 or 2
01234   // 1 standard
01235   // 2 try to improve minimum (slower)
01236   rmin.mnexcm ("SET STR", arglis, 1, ierror);
01237 
01238   arglis[0] = 10001;
01239   // Set the random seed for the generator used in SEEk to a fixed value for reproducibility
01240   rmin.mnexcm("SET RAN", arglis, 1, ierror);
01241 
01242   // Set fit parameters
01243   // ------------------
01244   double * Start = new double[parnumberAll];
01245   double * Step  = new double[parnumberAll];
01246   double * Mini  = new double[parnumberAll];
01247   double * Maxi  = new double[parnumberAll];
01248   int * ind = new int[parnumberAll]; // Order of release of parameters
01249   TString * parname = new TString[parnumberAll];
01250 
01251   if( !parResolStep.empty() && !parResolMin.empty() && !parResolMax.empty() ) {
01252     MuScleFitUtils::resolutionFunctionForVec->setParameters( Start, Step, Mini, Maxi, ind, parname, parResol, parResolOrder, parResolStep, parResolMin, parResolMax, MuonType );
01253   }
01254   else {
01255     MuScleFitUtils::resolutionFunctionForVec->setParameters( Start, Step, Mini, Maxi, ind, parname, parResol, parResolOrder, MuonType );
01256   }
01257 
01258   // Take the number of parameters in the resolutionFunction and displace the arrays passed to the scaleFunction
01259   int resParNum = MuScleFitUtils::resolutionFunctionForVec->parNum();
01260 
01261   if( !parScaleStep.empty() && !parScaleMin.empty() && !parScaleMax.empty() ) {
01262     MuScleFitUtils::scaleFunctionForVec->setParameters( &(Start[resParNum]), &(Step[resParNum]),
01263                                                         &(Mini[resParNum]), &(Maxi[resParNum]),
01264                                                         &(ind[resParNum]), &(parname[resParNum]),
01265                                                         parScale, parScaleOrder, parScaleStep,
01266                                                         parScaleMin, parScaleMax, MuonType );
01267   }
01268   else {
01269     MuScleFitUtils::scaleFunctionForVec->setParameters( &(Start[resParNum]), &(Step[resParNum]),
01270                                                         &(Mini[resParNum]), &(Maxi[resParNum]),
01271                                                         &(ind[resParNum]), &(parname[resParNum]),
01272                                                         parScale, parScaleOrder, MuonType );
01273   }
01274 
01275   // Initialize cross section parameters
01276   int crossSectionParShift = resParNum + MuScleFitUtils::scaleFunctionForVec->parNum();
01277   MuScleFitUtils::crossSectionHandler->setParameters( &(Start[crossSectionParShift]), &(Step[crossSectionParShift]), &(Mini[crossSectionParShift]),
01278                                                       &(Maxi[crossSectionParShift]), &(ind[crossSectionParShift]), &(parname[crossSectionParShift]),
01279                                                       parCrossSection, parCrossSectionOrder, resfind );
01280 
01281   // Initialize background parameters
01282   int bgrParShift = crossSectionParShift + crossSectionHandler->parNum();
01283   MuScleFitUtils::backgroundHandler->setParameters( &(Start[bgrParShift]), &(Step[bgrParShift]), &(Mini[bgrParShift]), &(Maxi[bgrParShift]),
01284                                                     &(ind[bgrParShift]), &(parname[bgrParShift]), parBgr, parBgrOrder, MuonType );
01285 
01286   for( int ipar=0; ipar<parnumber; ++ipar ) {
01287     std::cout << "parname["<<ipar<<"] = " << parname[ipar] << std::endl;
01288     std::cout << "Start["<<ipar<<"] = " << Start[ipar] << std::endl;
01289     std::cout << "Step["<<ipar<<"] = " << Step[ipar] << std::endl;
01290     std::cout << "Mini["<<ipar<<"] = " << Mini[ipar] << std::endl;
01291     std::cout << "Maxi["<<ipar<<"] = " << Maxi[ipar] << std::endl;
01292 
01293 
01294     rmin.mnparm( ipar, parname[ipar], Start[ipar], Step[ipar], Mini[ipar], Maxi[ipar], ierror );
01295 
01296 
01297     // Testing without limits
01298     // rmin.mnparm( ipar, parname[ipar], Start[ipar], Step[ipar], 0, 0, ierror );
01299 
01300 
01301   }
01302 
01303   // Do minimization
01304   // ---------------
01305   if (debug>19)
01306     std::cout << "[MuScleFitUtils-minimizeLikelihood]: Starting minimization" << std::endl;
01307   double fmin;
01308   double fdem;
01309   double errdef;
01310   int npari;
01311   int nparx;
01312   rmin.mnexcm ("CALL FCN", arglis, 1, ierror);
01313 
01314   // First, fix all parameters
01315   // -------------------------
01316   if (debug>19)
01317     std::cout << "[MuScleFitUtils-minimizeLikelihood]: First fix all parameters ...";
01318   for (int ipar=0; ipar<parnumber; ipar++) {
01319     rmin.FixParameter (ipar);
01320   }
01321 
01322   // Then release them in the specified order and refit
01323   // --------------------------------------------------
01324   if (debug>19) std::cout << " Then release them in order..." << std::endl;
01325 
01326   TString name;
01327   double pval;
01328   double pmin;
01329   double pmax;
01330   double errp;
01331   double errl;
01332   double errh;
01333   int ivar;
01334   double erro;
01335   double cglo;
01336   int n_times = 0;
01337   // n_times = number of loops required to unlock all parameters.
01338 
01339   if (debug>19) std::cout << "Before scale parNum" << std::endl;
01340   int scaleParNum = scaleFunction->parNum();
01341   if (debug>19) std::cout << "After scale parNum" << std::endl;
01342   // edm::LogInfo("minimizeLikelihood") << "number of parameters for scaleFunction = " << scaleParNum << std::endl;
01343   // edm::LogInfo("minimizeLikelihood") << "number of parameters for resolutionFunction = " << resParNum << std::endl;
01344   // edm::LogInfo("minimizeLikelihood") << "number of parameters for cross section = " << crossSectionHandler->parNum() << std::endl;
01345   // edm::LogInfo("minimizeLikelihood") << "number of parameters for backgroundFunction = " << parBgr.size() << std::endl;
01346   std::cout << "number of parameters for scaleFunction = " << scaleParNum << std::endl;
01347   std::cout << "number of parameters for resolutionFunction = " << resParNum << std::endl;
01348   std::cout << "number of parameters for cross section = " << crossSectionHandler->parNum() << std::endl;
01349   std::cout << "number of parameters for backgroundFunction = " << parBgr.size() << std::endl;
01350   // edm::LogInfo("minimizeLikelihood") << "number of parameters for backgroundFunction = " << backgroundFunction->parNum() << std::endl;
01351 
01352   for (int i=0; i<parnumber; i++) {
01353     // NB ind[] has been set as parorder[] previously
01354     if (n_times<ind[i]) {
01355       edm::LogInfo("minimizeLikelihood") << "n_times = " << n_times << ", ind["<<i<<"] = " << ind[i] << ", scaleParNum = " << scaleParNum << ", doScaleFit["<<loopCounter<<"] = " << doScaleFit[loopCounter] << std::endl;
01356       // Set the n_times only if we will do the fit
01357       if ( i<resParNum ) {
01358         if( doResolFit[loopCounter] ) n_times = ind[i];
01359       }
01360       else if( i<resParNum+scaleParNum ) {
01361         if( doScaleFit[loopCounter] ) n_times = ind[i];
01362       }
01363       else if( doBackgroundFit[loopCounter] ) n_times = ind[i];
01364     }
01365   }
01366   for (int iorder=0; iorder<n_times+1; iorder++) { // Repeat fit n_times times
01367     std::cout << "Starting minimization " << iorder << " of " << n_times << std::endl;
01368 
01369     bool somethingtodo = false;
01370 
01371     // Use parameters from cfg to select which fit to do
01372     // -------------------------------------------------
01373     if( doResolFit[loopCounter] ) {
01374       // Release resolution parameters and fit them
01375       // ------------------------------------------
01376       for( unsigned int ipar=0; ipar<parResol.size(); ++ipar ) {
01377         if( parfix[ipar]==0 && ind[ipar]==iorder ) {
01378           rmin.Release( ipar );
01379           somethingtodo = true;
01380         }
01381       }
01382     }
01383     if( doScaleFit[loopCounter] ) {
01384       // Release scale parameters and fit them
01385       // -------------------------------------
01386       for( unsigned int ipar=parResol.size(); ipar<parResol.size()+parScale.size(); ++ipar ) {
01387         if( parfix[ipar]==0 && ind[ipar]==iorder ) { // parfix=0 means parameter is free
01388           rmin.Release( ipar );
01389           somethingtodo = true;
01390         }
01391       }
01392       scaleFitNotDone_ = false;
01393     }
01394     unsigned int crossSectionParShift = parResol.size()+parScale.size();
01395     if( doCrossSectionFit[loopCounter] ) {
01396       // Release cross section parameters and fit them
01397       // ---------------------------------------------
01398       // Note that only cross sections of resonances that are being fitted are released
01399       bool doCrossSection = crossSectionHandler->releaseParameters( rmin, resfind, parfix, ind, iorder, crossSectionParShift );
01400       if( doCrossSection ) somethingtodo = true;
01401     }
01402     if( doBackgroundFit[loopCounter] ) {
01403       // Release background parameters and fit them
01404       // ------------------------------------------
01405       // for( int ipar=parResol.size()+parScale.size(); ipar<parnumber; ++ipar ) {
01406       // Free only the parameters for the regions, as the resonances intervals are never used to fit the background
01407       unsigned int bgrParShift = crossSectionParShift+crossSectionHandler->parNum();
01408       for( unsigned int ipar = bgrParShift; ipar < bgrParShift+backgroundHandler->regionsParNum(); ++ipar ) {
01409         // Release only those parameters for the resonances we are fitting
01410         if( parfix[ipar]==0 && ind[ipar]==iorder && backgroundHandler->unlockParameter(resfind, ipar - bgrParShift) ) {
01411           rmin.Release( ipar );
01412           somethingtodo = true;
01413         }
01414       }
01415     }
01416 
01417     // OK, now do minimization if some parameter has been released
01418     // -----------------------------------------------------------
01419     if( somethingtodo ) {
01420 // #ifdef DEBUG
01421 
01422       std::stringstream fileNum;
01423       fileNum << loopCounter;
01424 
01425       minuitLoop_ = 0;
01426       char name[50];
01427       sprintf(name, "likelihoodInLoop_%d_%d", loopCounter, iorder);
01428       TH1D * tempLikelihoodInLoop = new TH1D(name, "likelihood value in minuit loop", 10000, 0, 10000);
01429       likelihoodInLoop_ = tempLikelihoodInLoop;
01430       char signalProbName[50];
01431       sprintf(signalProbName, "signalProb_%d_%d", loopCounter, iorder);
01432       TH1D * tempSignalProb = new TH1D(signalProbName, "signal probability", 10000, 0, 10000);
01433       signalProb_ = tempSignalProb;
01434       char backgroundProbName[50];
01435       sprintf(backgroundProbName, "backgroundProb_%d_%d", loopCounter, iorder);
01436       TH1D * tempBackgroundProb = new TH1D(backgroundProbName, "background probability", 10000, 0, 10000);
01437       backgroundProb_ = tempBackgroundProb;
01438 // #endif
01439 
01440 
01441 
01442       // Before we start the minimization we create a list of events with only the events inside a smaller
01443       // window than the one in which the probability is != 0. We will compute the probability for all those
01444       // events and hopefully the margin will avoid them to get a probability = 0 (which causes a discontinuity
01445       // in the likelihood function). The width of this smaller window must be optimized, but we can start using
01446       // an 90% of the normalization window.
01447       double protectionFactor = 0.9;
01448 
01449       MuScleFitUtils::ReducedSavedPair.clear();
01450       for( unsigned int nev=0; nev<MuScleFitUtils::SavedPair.size(); ++nev ) {
01451         const lorentzVector * recMu1 = &(MuScleFitUtils::SavedPair[nev].first);
01452         const lorentzVector * recMu2 = &(MuScleFitUtils::SavedPair[nev].second);
01453         double mass = MuScleFitUtils::invDimuonMass( *recMu1, *recMu2 );
01454         // Test all resonances to see if the mass is inside at least one of the windows
01455         bool check = false;
01456         for( int ires = 0; ires < 6; ++ires ) {
01457           // std::pair<double, double> windowFactor = backgroundHandler->windowFactors( doBackgroundFit[loopCounter], ires );
01458           std::pair<double, double> windowBorder = backgroundHandler->windowBorders( doBackgroundFit[loopCounter], ires );
01459           // if( resfind[ires] && checkMassWindow( mass, ires, backgroundHandler->resMass( doBackgroundFit[loopCounter], ires ),
01460           //                                       0.9*windowFactor.first, 0.9*windowFactor.second ) ) {
01461           // double resMassValue = backgroundHandler->resMass( doBackgroundFit[loopCounter], ires );
01462           // double windowBorderLeft = resMassValue - protectionFactor*(resMassValue - windowBorder.first);
01463           // double windowBorderRight = resMassValue + protectionFactor*(windowBorder.second - resMassValue);
01464           double windowBorderShift = (windowBorder.second - windowBorder.first)*(1-protectionFactor)/2.;
01465           double windowBorderLeft = windowBorder.first + windowBorderShift;
01466           double windowBorderRight = windowBorder.second - windowBorderShift;
01467           if( resfind[ires] && checkMassWindow( mass, windowBorderLeft, windowBorderRight ) ) {
01468             check = true;
01469           }
01470         }
01471         if( check ) {
01472           MuScleFitUtils::ReducedSavedPair.push_back(std::make_pair(*recMu1, *recMu2));
01473         }
01474       }
01475       std::cout << "Fitting with " << MuScleFitUtils::ReducedSavedPair.size() << " events" << std::endl;
01476 
01477 
01478       // rmin.SetMaxIterations(500*parnumber);
01479 
01480       //Print some informations
01481       std::cout<<"MINUIT is starting the minimization for the iteration number "<<loopCounter<<std::endl;
01482 
01483       //Try to set iterations
01484       //      rmin.SetMaxIterations(100000);
01485 
01486       std::cout<<"maxNumberOfIterations (just set) = "<<rmin.GetMaxIterations()<<std::endl;
01487 
01488       MuScleFitUtils::normalizationChanged_ = 0;
01489 
01490       // Maximum number of iterations
01491       arglis[0] = 100000;
01492       // tolerance 
01493       arglis[1] = 0.1;
01494 
01495       // Run simplex first to get an initial estimate of the minimum
01496       if( startWithSimplex_ ) {
01497         rmin.mnexcm( "SIMPLEX", arglis, 0, ierror );
01498       }
01499 
01500       rmin.mnexcm( "MIGRAD", arglis, 2, ierror ); 
01501 
01502 
01503 
01504 
01505 // #ifdef DEBUG
01506       likelihoodInLoop_->Write();
01507       signalProb_->Write();
01508       backgroundProb_->Write();
01509       delete tempLikelihoodInLoop;
01510       delete tempSignalProb;
01511       delete tempBackgroundProb;
01512       likelihoodInLoop_ = 0;
01513       signalProb_ = 0;
01514       backgroundProb_ = 0;
01515 // #endif
01516 
01517 
01518       // Compute again the error matrix
01519       rmin.mnexcm( "HESSE", arglis, 0, ierror );
01520 
01521       // Peform minos error analysis.
01522       if( computeMinosErrors_ ) {
01523         duringMinos_ = true;
01524         rmin.mnexcm( "MINOS", arglis, 0, ierror );
01525         duringMinos_ = false;
01526       }
01527 
01528       if( normalizationChanged_ > 1 ) {
01529         std::cout << "WARNING: normalization changed during fit meaning that events exited from the mass window. This causes a discontinuity in the likelihood function. Please check the scan of the likelihood as a function of the parameters to see if there are discontinuities around the minimum." << std::endl;
01530       }
01531     }
01532 
01533     // bool notWritten = true;
01534     for (int ipar=0; ipar<parnumber; ipar++) {
01535 
01536       rmin.mnpout (ipar, name, pval, erro, pmin, pmax, ivar);
01537       // Save parameters in parvalue[] vector
01538       // ------------------------------------
01539       if (ierror!=0 && debug>0) {
01540         std::cout << "[MuScleFitUtils-minimizeLikelihood]: ierror!=0, bogus pars" << std::endl;
01541       }
01542       //     for (int ipar=0; ipar<parnumber; ipar++) {
01543       //       rmin.mnpout (ipar, name, pval, erro, pmin, pmax, ivar);
01544       parvalue[loopCounter][ipar] = pval;
01545       //     }
01546 
01547 
01548       // int ilax2 = 0;
01549       // Double_t val2pl, val2mi;
01550       // rmin.mnmnot (ipar+1, ilax2, val2pl, val2mi);
01551       rmin.mnerrs (ipar, errh, errl, errp, cglo);
01552 
01553 
01554       // Set error on params
01555       // -------------------
01556       if (errp!=0) {
01557         parerr[3*ipar] = errp;
01558       } else {
01559         parerr[3*ipar] = (((errh)>(fabs(errl)))?(errh):(fabs(errl)));
01560       }
01561       parerr[3*ipar+1] = errl;
01562       parerr[3*ipar+2] = errh;
01563 
01564       if( ipar == 0 ) {
01565         FitParametersFile << " Resolution fit parameters:" << std::endl;
01566       }
01567       if( ipar == int(parResol.size()) ) {
01568         FitParametersFile << " Scale fit parameters:" << std::endl;
01569       }
01570       if( ipar == int(parResol.size()+parScale.size()) ) {
01571         FitParametersFile << " Cross section fit parameters:" << std::endl;
01572       }
01573       if( ipar == int(parResol.size()+parScale.size()+crossSectionHandler->parNum()) ) {
01574         FitParametersFile << " Background fit parameters:" << std::endl;
01575       }
01576 //       if( ipar >= int(parResol.size()+parScale.size()) && ipar < int(parResol.size()+parScale.size()+crossSectionHandler->parNum()) && notWritted ) {
01577 
01578 //         std::vector<double> relativeCrossSections = crossSectionHandler->relativeCrossSections(&(parvalue[loopCounter][parResol.size()+parScale.size()]));
01579 //         std::vector<double>::const_iterator it = relativeCrossSections.begin();
01580 //         for( ; it != relativeCrossSections.end(); ++it ) {
01581 //           FitParametersFile << "  Results of the fit: parameter " << ipar << " has value "
01582 //                             << *it << "+-" << 0
01583 //                             << " + " << 0 << " - " << 0
01584 //                             << " /t/t (" << 0 << ")" << std::endl;
01585 //         }
01586 
01587 //         notWritten = false;
01588 //       }
01589 //       else {
01590         FitParametersFile << "  Results of the fit: parameter " << ipar << " has value "
01591                           << pval << "+-" << parerr[3*ipar]
01592                           << " + " << parerr[3*ipar+1] << " - " << parerr[3*ipar+2]
01593                           // << " \t\t (" << parname[ipar] << ")"
01594                           << std::endl;
01595 
01596 
01597 
01598     }
01599     rmin.mnstat (fmin, fdem, errdef, npari, nparx, istat); // NNBB Commented for a check!
01600     FitParametersFile << std::endl;
01601 
01602     if( minimumShapePlots_ ) {
01603       // Create plots of the minimum vs parameters
01604       // -----------------------------------------
01605       // Keep this after the parameters filling because it recomputes the values and it can compromise the fit results.
01606       if( somethingtodo ) {
01607         std::stringstream iorderString;
01608         iorderString << iorder;
01609         std::stringstream iLoopString;
01610         iLoopString << loopCounter;
01611 
01612         for (int ipar=0; ipar<parnumber; ipar++) {
01613           if( parfix[ipar] == 1 ) continue;
01614           std::cout << "plotting parameter = " << ipar+1 << std::endl;
01615           std::stringstream iparString;
01616           iparString << ipar+1;
01617           std::stringstream iparStringName;
01618           iparStringName << ipar;
01619           rmin.mncomd( ("scan "+iparString.str()).c_str(), ierror );
01620           if( ierror == 0 ) {
01621             TCanvas * canvas = new TCanvas(("likelihoodCanvas_loop_"+iLoopString.str()+"_oder_"+iorderString.str()+"_par_"+iparStringName.str()).c_str(), ("likelihood_"+iparStringName.str()).c_str(), 1000, 800);
01622             canvas->cd();
01623             // arglis[0] = ipar;
01624             // rmin.mnexcm( "SCA", arglis, 0, ierror );
01625             TGraph * graph = (TGraph*)rmin.GetPlot();
01626             graph->Draw("AP");
01627             // graph->SetTitle(("parvalue["+iparStringName.str()+"]").c_str());
01628             graph->SetTitle(parname[ipar]);
01629             // graph->Write();
01630 
01631             canvas->Write();
01632           }
01633         }
01634 
01635         //       // Draw contours of the fit
01636         //       TCanvas * canvas = new TCanvas(("contourCanvas_oder_"+iorderString.str()).c_str(), "contour", 1000, 800);
01637         //       canvas->cd();
01638         //       TGraph * contourGraph = (TGraph*)rmin.Contour(4, 2, 4);
01639         //       if( (rmin.GetStatus() == 0) || (rmin.GetStatus() >= 3) ) {
01640         //         contourGraph->Draw("AP");
01641         //       }
01642         //       else {
01643         //         std::cout << "Contour graph error: status = " << rmin.GetStatus() << std::endl;
01644         //       }
01645         //       canvas->Write();
01646       }
01647     }
01648 
01649   } // end loop on iorder
01650   FitParametersFile.close();
01651 
01652   std::cout << "[MuScleFitUtils-minimizeLikelihood]: Parameters after likelihood " << std::endl;
01653   for (unsigned int ipar=0; ipar<(unsigned int)parnumber; ipar++) {
01654     std::cout << ipar << " " << parvalue[loopCounter][ipar] << " : free = "
01655          << parfix[ipar] << "; order = " << parorder[ipar] << std::endl;
01656   }
01657 
01658   // Put back parvalue into parResol, parScale, parCrossSection, parBgr
01659   // ------------------------------------------------------------------
01660   for( int i=0; i<(int)(parResol.size()); ++i ) {
01661     parResol[i] = parvalue[loopCounter][i];
01662   }
01663   for( int i=0; i<(int)(parScale.size()); ++i ) {
01664     parScale[i] = parvalue[loopCounter][i+parResol.size()];
01665   }
01666   parCrossSection = crossSectionHandler->relativeCrossSections(&(parvalue[loopCounter][parResol.size()+parScale.size()]), resfind);
01667   for( unsigned int i=0; i<parCrossSection.size(); ++i ) {
01668   //   parCrossSection[i] = parvalue[loopCounter][i+parResol.size()+parScale.size()];
01669     std::cout << "relative cross section["<<i<<"] = " << parCrossSection[i] << std::endl;
01670   }
01671   // Save only the fitted background parameters
01672   for( unsigned int i = 0; i<(parBgr.size() - parForResonanceWindows); ++i ) {
01673     parBgr[i] = parvalue[loopCounter][i+parResol.size()+parScale.size()+crossSectionHandler->parNum()];
01674   }
01675 
01676   // Background rescaling from regions to resonances
01677   // -----------------------------------------------
01678   // Only if we fitted the background
01679   if( doBackgroundFit[loopCounter] ) {
01680     // This rescales from regions to resonances
01681     int localMuonType = MuonType;
01682     if( MuonType > 2 ) localMuonType = 2;
01683     backgroundHandler->rescale( parBgr, ResMass, massWindowHalfWidth[localMuonType],
01684                                 MuScleFitUtils::ReducedSavedPair );
01685   }
01686 
01687   // Delete the arrays used to set some parameters
01688   delete[] Start;
01689   delete[] Step;
01690   delete[] Mini;
01691   delete[] Maxi;
01692   delete[] ind;
01693   delete[] parname;
01694 }
01695 
01696 // Likelihood function
01697 // -------------------
01698 extern "C" void likelihood( int& npar, double* grad, double& fval, double* xval, int flag ) {
01699 
01700   if (MuScleFitUtils::debug>19) std::cout << "[MuScleFitUtils-likelihood]: In likelihood function" << std::endl;
01701 
01702   const lorentzVector * recMu1;
01703   const lorentzVector * recMu2;
01704   lorentzVector corrMu1;
01705   lorentzVector corrMu2;
01706 
01707   //   if (MuScleFitUtils::debug>19) {
01708   //     int parnumber = (int)(MuScleFitUtils::parResol.size()+MuScleFitUtils::parScale.size()+
01709   //                           MuScleFitUtils::parCrossSection.size()+MuScleFitUtils::parBgr.size());
01710   //     std::cout << "[MuScleFitUtils-likelihood]: Looping on tree with ";
01711   //     for (int ipar=0; ipar<parnumber; ipar++) {
01712   //       std::cout << "Parameter #" << ipar << " with value " << xval[ipar] << " ";
01713   //     }
01714   //     std::cout << std::endl;
01715   //   }
01716 
01717   // Loop on the tree
01718   // ----------------
01719   double flike = 0;
01720   int evtsinlik = 0;
01721   int evtsoutlik = 0;
01722   // std::cout << "SavedPair.size() = " << MuScleFitUtils::SavedPair.size() << std::endl;
01723   if( MuScleFitUtils::debug>0 ) {
01724     std::cout << "SavedPair.size() = " << MuScleFitUtils::SavedPair.size() << std::endl;
01725     std::cout << "ReducedSavedPair.size() = " << MuScleFitUtils::ReducedSavedPair.size() << std::endl;
01726   }
01727   // for( unsigned int nev=0; nev<MuScleFitUtils::SavedPair.size(); ++nev ) {
01728   for( unsigned int nev=0; nev<MuScleFitUtils::ReducedSavedPair.size(); ++nev ) {
01729 
01730     //     recMu1 = &(MuScleFitUtils::SavedPair[nev].first);
01731     //     recMu2 = &(MuScleFitUtils::SavedPair[nev].second);
01732     recMu1 = &(MuScleFitUtils::ReducedSavedPair[nev].first);
01733     recMu2 = &(MuScleFitUtils::ReducedSavedPair[nev].second);
01734 
01735     // Compute original mass
01736     // ---------------------
01737     double mass = MuScleFitUtils::invDimuonMass( *recMu1, *recMu2 );
01738 
01739     // Compute weight and reference mass (from original mass)
01740     // ------------------------------------------------------
01741     double weight = MuScleFitUtils::computeWeight(mass, MuScleFitUtils::iev_);
01742     if( weight!=0. ) {
01743       // Compute corrected mass (from previous biases) only if we are currently fitting the scale
01744       // ----------------------------------------------------------------------------------------
01745       if( MuScleFitUtils::doScaleFit[MuScleFitUtils::loopCounter] ) {
01746 //      std::cout << "Original pt1 = " << corrMu1.Pt() << std::endl;
01747 //      std::cout << "Original pt2 = " << corrMu2.Pt() << std::endl;
01748         corrMu1 = MuScleFitUtils::applyScale(*recMu1, xval, -1);
01749         corrMu2 = MuScleFitUtils::applyScale(*recMu2, xval,  1);
01750         
01751 //         if( (corrMu1.Pt() != corrMu1.Pt()) || (corrMu2.Pt() != corrMu2.Pt()) ) {
01752 //           std::cout << "Rescaled pt1 = " << corrMu1.Pt() << std::endl;
01753 //           std::cout << "Rescaled pt2 = " << corrMu2.Pt() << std::endl;
01754 //         }
01755 //      std::cout << "Rescaled pt1 = " << corrMu1.Pt() << std::endl;
01756 //      std::cout << "Rescaled pt2 = " << corrMu2.Pt() << std::endl;
01757       }
01758       else {
01759         corrMu1 = *recMu1;
01760         corrMu2 = *recMu2;
01761 
01762 //         if( (corrMu1.Pt() != corrMu1.Pt()) || (corrMu2.Pt() != corrMu2.Pt()) ) {
01763 //           std::cout << "Not rescaled pt1 = " << corrMu1.Pt() << std::endl;
01764 //           std::cout << "Not rescaled pt2 = " << corrMu2.Pt() << std::endl;
01765 //         }
01766       }
01767       double corrMass = MuScleFitUtils::invDimuonMass(corrMu1, corrMu2);
01768       double Y = (corrMu1+corrMu2).Rapidity();
01769       double resEta = (corrMu1+corrMu2).Eta();
01770       if( MuScleFitUtils::debug>19 ) {
01771         std::cout << "[MuScleFitUtils-likelihood]: Original/Corrected resonance mass = " << mass
01772              << " / " << corrMass << std::endl;
01773       }
01774 
01775       // Compute mass resolution
01776       // -----------------------
01777       double massResol = MuScleFitUtils::massResolution(corrMu1, corrMu2, xval);
01778       if (MuScleFitUtils::debug>19)
01779         std::cout << "[MuScleFitUtils-likelihood]: Resolution is " << massResol << std::endl;
01780 
01781       // Compute probability of this mass value including background modeling
01782       // --------------------------------------------------------------------
01783       if (MuScleFitUtils::debug>1) std::cout << "calling massProb inside likelihood function" << std::endl;
01784 
01785       // double prob = MuScleFitUtils::massProb( corrMass, resEta, Y, massResol, xval );
01786       double prob = MuScleFitUtils::massProb( corrMass, resEta, Y, massResol, xval, false, corrMu1.eta(), corrMu2.eta() );
01787       if (MuScleFitUtils::debug>1) std::cout << "likelihood:massProb = " << prob << std::endl;
01788 
01789       // Compute likelihood
01790       // ------------------
01791       if( prob>0 ) {
01792         // flike += log(prob*10000)*weight; // NNBB! x10000 to see if we can recover the problem of boundary
01793         flike += log(prob)*weight;
01794         evtsinlik += 1;  // NNBB test: see if likelihood per event is smarter (boundary problem)
01795       } else {
01796         if( MuScleFitUtils::debug > 0 ) {
01797           std::cout << "WARNING: corrMass = " << corrMass << " outside window, this will cause a discontinuity in the likelihood. Consider increasing the safety bands which are now set to 90% of the normalization window to avoid this problem" << std::endl;
01798           std::cout << "Original mass was = " << mass << std::endl;
01799           std::cout << "WARNING: massResol = " << massResol << " outside window" << std::endl;
01800         }
01801         evtsoutlik += 1;
01802       }
01803       if (MuScleFitUtils::debug>19)
01804         std::cout << "[MuScleFitUtils-likelihood]: Mass probability = " << prob << std::endl;
01805     } // weight!=0
01806 
01807   } // End of loop on tree events
01808 
01809 //   // Protection for low statistic. If the likelihood manages to throw out all the signal
01810 //   // events and stays with ~ 10 events in the resonance window it could have a better likelihood
01811 //   // because of ~ uniformly distributed events (a random combination could be good and spoil the fit).
01812 //   // We require that the number of events included in the fit does not change more than 5% in each minuit loop.
01813 //   bool lowStatPenalty = false;
01814 //   if( MuScleFitUtils::minuitLoop_ > 0 ) {
01815 //     double newEventsOutInRatio = double(evtsinlik);
01816 //     // double newEventsOutInRatio = double(evtsoutlik)/double(evtsinlik);
01817 //     double ratio = newEventsOutInRatio/MuScleFitUtils::oldEventsOutInRatio_;
01818 //     MuScleFitUtils::oldEventsOutInRatio_ = newEventsOutInRatio;
01819 //     if( ratio < 0.8 || ratio > 1.2 ) {
01820 //       std::cout << "Warning: too much change from oldEventsInLikelihood to newEventsInLikelihood, ratio is = " << ratio << std::endl;
01821 //       std::cout << "oldEventsInLikelihood = " << MuScleFitUtils::oldEventsOutInRatio_ << ", newEventsInLikelihood = " << newEventsOutInRatio << std::endl;
01822 //       lowStatPenalty = true;
01823 //     }
01824 //   }
01825 
01826   // It is a product of probabilities, we compare the sqrt_N of them. Thus N becomes a denominator of the logarithm.
01827   if( evtsinlik != 0 ) {
01828 
01829     if( MuScleFitUtils::normalizeLikelihoodByEventNumber_ ) {
01830       // && !(MuScleFitUtils::duringMinos_) ) {
01831       if( MuScleFitUtils::rminPtr_ == 0 ) {
01832         std::cout << "ERROR: rminPtr_ = " << MuScleFitUtils::rminPtr_ << ", code will crash" << std::endl;
01833       }
01834       double normalizationArg[] = {1/double(evtsinlik)};
01835       // Reset the normalizationArg only if it changed
01836       if( MuScleFitUtils::oldNormalization_ != normalizationArg[0] ) {
01837         int ierror = 0;
01838 //         if( MuScleFitUtils::likelihoodInLoop_ != 0 ) {
01839 //           // This condition is set only when minimizing. Later calls of hesse and minos will not change the value
01840 //           // This is done to avoid minos being confused by changing the UP parameter during its computation.
01841 //           MuScleFitUtils::rminPtr_->mnexcm("SET ERR", normalizationArg, 1, ierror);
01842 //         }
01843         MuScleFitUtils::rminPtr_->mnexcm("SET ERR", normalizationArg, 1, ierror);
01844         std::cout << "oldNormalization = " << MuScleFitUtils::oldNormalization_ << " new = " << normalizationArg[0] << std::endl;
01845         MuScleFitUtils::oldNormalization_ = normalizationArg[0];
01846         MuScleFitUtils::normalizationChanged_ += 1;
01847       }
01848       fval = -2.*flike/double(evtsinlik);
01849       // fval = -2.*flike;
01850       //     if( lowStatPenalty ) {
01851       //       fval *= 100;
01852       //     }
01853     }
01854     else {
01855       fval = -2.*flike;
01856     }
01857   }
01858   else {
01859     std::cout << "Problem: Events in likelihood = " << evtsinlik << std::endl;
01860     fval = 999999999.;
01861   }
01862   // fval = -2.*flike;
01863   if (MuScleFitUtils::debug>19)
01864     std::cout << "[MuScleFitUtils-likelihood]: End tree loop with likelihood value = " << fval << std::endl;
01865 
01866 //  #ifdef DEBUG
01867 
01868 //  if( MuScleFitUtils::minuitLoop_ < 10000 ) {
01869   if( MuScleFitUtils::likelihoodInLoop_ != 0 ) {
01870     ++MuScleFitUtils::minuitLoop_;
01871     MuScleFitUtils::likelihoodInLoop_->SetBinContent(MuScleFitUtils::minuitLoop_, fval);
01872   }
01873   //  }
01874   // else std::cout << "minuitLoop over 10000. Not filling histogram" << std::endl;
01875 
01876   std::cout<<"MINUIT loop number "<<MuScleFitUtils::minuitLoop_<<", likelihood = "<<fval<<std::endl;
01877 
01878   if( MuScleFitUtils::debug > 0 ) {
01879     //     if( MuScleFitUtils::duringMinos_ ) {
01880     //       int parnumber = (int)(MuScleFitUtils::parResol.size()+MuScleFitUtils::parScale.size()+
01881     //                             MuScleFitUtils::parCrossSection.size()+MuScleFitUtils::parBgr.size());
01882     //       std::cout << "[MuScleFitUtils-likelihood]: Looping on tree with ";
01883     //       for (int ipar=0; ipar<parnumber; ipar++) {
01884     //         std::cout << "Parameter #" << ipar << " with value " << xval[ipar] << " ";
01885     //       }
01886     //       std::cout << std::endl;
01887     //       std::cout << "[MuScleFitUtils-likelihood]: likelihood value = " << fval << std::endl;
01888     //     }
01889     std::cout << "Events in likelihood = " << evtsinlik << std::endl;
01890     std::cout << "Events out likelihood = " << evtsoutlik << std::endl;
01891   }
01892 
01893 //  #endif
01894 }
01895 
01896 // Mass fitting routine
01897 // --------------------
01898 std::vector<TGraphErrors*> MuScleFitUtils::fitMass (TH2F* histo) {
01899 
01900   if (MuScleFitUtils::debug>0) std::cout << "Fitting " << histo->GetName() << std::endl;
01901 
01902   std::vector<TGraphErrors *> results;
01903 
01904   // Results of the fit
01905   // ------------------
01906   std::vector<double> Ftop;
01907   std::vector<double> Fwidth;
01908   std::vector<double> Fmass;
01909   std::vector<double> Etop;
01910   std::vector<double> Ewidth;
01911   std::vector<double> Emass;
01912   std::vector<double> Fchi2;
01913   // X bin center and width
01914   // ----------------------
01915   std::vector<double> Xcenter;
01916   std::vector<double> Ex;
01917 
01918   // Fit with lorentzian peak
01919   // ------------------------
01920   TF1 *fitFcn = new TF1 ("fitFcn", lorentzianPeak, 70, 110, 3);
01921   fitFcn->SetParameters (100, 3, 91);
01922   fitFcn->SetParNames ("Ftop", "Fwidth", "Fmass");
01923   fitFcn->SetLineWidth (2);
01924 
01925   // Fit slices projected along Y from bins in X
01926   // -------------------------------------------
01927   double cont_min = 20;    // Minimum number of entries
01928   Int_t binx = histo->GetXaxis()->GetNbins();
01929   // TFile *f= new TFile("prova.root", "recreate");
01930   // histo->Write();
01931   for (int i=1; i<=binx; i++) {
01932     TH1 * histoY = histo->ProjectionY ("", i, i);
01933     // histoY->Write();
01934     double cont = histoY->GetEntries();
01935     if (cont>cont_min) {
01936       histoY->Fit ("fitFcn", "0", "", 70, 110);
01937       double *par = fitFcn->GetParameters();
01938       double *err = fitFcn->GetParErrors();
01939 
01940       Ftop.push_back(par[0]);
01941       Fwidth.push_back(par[1]);
01942       Fmass.push_back(par[2]);
01943       Etop.push_back(err[0]);
01944       Ewidth.push_back(err[1]);
01945       Emass.push_back(err[2]);
01946 
01947       double chi2 = fitFcn->GetChisquare();
01948       Fchi2.push_back(chi2);
01949 
01950       double xx = histo->GetXaxis()->GetBinCenter(i);
01951       Xcenter.push_back(xx);
01952       double ex = 0; // FIXME: you can use the bin width
01953       Ex.push_back(ex);
01954     }
01955   }
01956   // f->Close();
01957 
01958   // Put the fit results in arrays for TGraphErrors
01959   // ----------------------------------------------
01960   const int nn = Fmass.size();
01961   double *x    = new double[nn];
01962   double *ym   = new double[nn];
01963   double *e    = new double[nn];
01964   double *eym  = new double[nn];
01965   double *yw   = new double[nn];
01966   double *eyw  = new double[nn];
01967   double *yc   = new double[nn];
01968 
01969   for (int j=0; j<nn; j++) {
01970     x[j]   = Xcenter[j];
01971     ym[j]  = Fmass[j];
01972     eym[j] = Emass[j];
01973     yw[j]  = Fwidth[j];
01974     eyw[j] = Ewidth[j];
01975     yc[j]  = Fchi2[j];
01976     e[j]   = Ex[j];
01977   }
01978 
01979   // Create TGraphErrors
01980   // -------------------
01981   TString name = histo->GetName();
01982   TGraphErrors *grM = new TGraphErrors (nn, x, ym, e, eym);
01983   grM->SetTitle (name+"_M");
01984   grM->SetName (name+"_M");
01985   TGraphErrors *grW = new TGraphErrors (nn, x, yw, e, eyw);
01986   grW->SetTitle (name+"_W");
01987   grW->SetName (name+"_W");
01988   TGraphErrors *grC = new TGraphErrors (nn, x, yc, e, e);
01989   grC->SetTitle (name+"_chi2");
01990   grC->SetName (name+"_chi2");
01991 
01992   // Cleanup
01993   // -------
01994   delete x;
01995   delete ym;
01996   delete eym;
01997   delete yw;
01998   delete eyw;
01999   delete yc;
02000   delete e;
02001   delete fitFcn;
02002 
02003   results.push_back(grM);
02004   results.push_back(grW);
02005   results.push_back(grC);
02006   return results;
02007 }
02008 
02009 // Resolution fitting routine
02010 // --------------------------
02011 std::vector<TGraphErrors*> MuScleFitUtils::fitReso (TH2F* histo) {
02012   std::cout << "Fitting " << histo->GetName() << std::endl;
02013   std::vector<TGraphErrors *> results;
02014 
02015   // Results from fit
02016   // ----------------
02017   std::vector<double> maxs;
02018   std::vector<double> means;
02019   std::vector<double> sigmas;
02020   std::vector<double> chi2s;
02021   std::vector<double> Emaxs;
02022   std::vector<double> Emeans;
02023   std::vector<double> Esigmas;
02024 
02025   // X bin center and width
02026   // ----------------------
02027   std::vector<double> Xcenter;
02028   std::vector<double> Ex;
02029 
02030   // Fit with a gaussian
02031   // -------------------
02032   TF1 *fitFcn = new TF1 ("fitFunc", Gaussian, -0.2, 0.2, 3);
02033   fitFcn->SetParameters (100, 0, 0.02);
02034   fitFcn->SetParNames ("max", "mean", "sigma");
02035   fitFcn->SetLineWidth (2);
02036 
02037   // Fit slices projected along Y from bins in X
02038   // -------------------------------------------
02039   double cont_min = 20;    // Minimum number of entries
02040   Int_t binx = histo->GetXaxis()->GetNbins();
02041   for (int i=1; i<=binx; i++) {
02042     TH1 * histoY = histo->ProjectionY ("", i, i);
02043     double cont = histoY->GetEntries();
02044     if (cont>cont_min) {
02045       histoY->Fit ("fitFunc", "0", "", -0.2, 0.2);
02046       double *par = fitFcn->GetParameters();
02047       double *err = fitFcn->GetParErrors();
02048 
02049       maxs.push_back (par[0]);
02050       means.push_back (par[1]);
02051       sigmas.push_back (par[2]);
02052       Emaxs.push_back (err[0]);
02053       Emeans.push_back (err[1]);
02054       Esigmas.push_back (err[2]);
02055 
02056       double chi2 = fitFcn->GetChisquare();
02057       chi2s.push_back (chi2);
02058 
02059       double xx = histo->GetXaxis()->GetBinCenter(i);
02060       Xcenter.push_back (xx);
02061       double ex = 0; // FIXME: you can use the bin width
02062       Ex.push_back (ex);
02063     }
02064   }
02065 
02066   // Put the fit results in arrays for TGraphErrors
02067   // ----------------------------------------------
02068   const int nn = means.size();
02069   double *x    = new double[nn];
02070   double *ym   = new double[nn];
02071   double *e    = new double[nn];
02072   double *eym  = new double[nn];
02073   double *yw   = new double[nn];
02074   double *eyw  = new double[nn];
02075   double *yc   = new double[nn];
02076 
02077   for (int j=0; j<nn; j++) {
02078     x[j]   = Xcenter[j];
02079     ym[j]  = means[j];
02080     eym[j] = Emeans[j];
02081     // yw[j]  = maxs[j];
02082     // eyw[j] = Emaxs[j];
02083     yw[j]  = sigmas[j];
02084     eyw[j] = Esigmas[j];
02085     yc[j]  = chi2s[j];
02086     e[j]   = Ex[j];
02087   }
02088 
02089   // Create TGraphErrors
02090   // -------------------
02091   TString name = histo->GetName();
02092   TGraphErrors *grM = new TGraphErrors (nn, x, ym, e, eym);
02093   grM->SetTitle (name+"_mean");
02094   grM->SetName (name+"_mean");
02095   TGraphErrors *grW = new TGraphErrors (nn, x, yw, e, eyw);
02096   grW->SetTitle (name+"_sigma");
02097   grW->SetName (name+"_sigma");
02098   TGraphErrors *grC = new TGraphErrors (nn, x, yc, e, e);
02099   grC->SetTitle (name+"_chi2");
02100   grC->SetName (name+"_chi2");
02101 
02102   // Cleanup
02103   // -------
02104   delete x;
02105   delete ym;
02106   delete eym;
02107   delete yw;
02108   delete eyw;
02109   delete yc;
02110   delete e;
02111   delete fitFcn;
02112 
02113   results.push_back (grM);
02114   results.push_back (grW);
02115   results.push_back (grC);
02116   return results;
02117 }
02118 
02119 // Mass probability for likelihood computation - no-background version (not used anymore)
02120 // --------------------------------------------------------------------------------------
02121 double MuScleFitUtils::massProb( const double & mass, const double & rapidity, const int ires, const double & massResol )
02122 {
02123   // This routine computes the likelihood that a given measured mass "measMass" is
02124   // the result of resonance #ires if the resolution expected for the two muons is massResol
02125   // ---------------------------------------------------------------------------------------
02126 
02127   double P = 0.;
02128 
02129   // Return Lorentz value for zero resolution cases (like MC)
02130   // --------------------------------------------------------
02131   if (massResol==0.) {
02132     if (debug>9) std::cout <<  "Mass, gamma , mref, width, P: " << mass
02133                       << " " << ResGamma[ires] << " " << ResMass[ires]<< " " << massResol
02134                       << " : used Lorentz P-value" << std::endl;
02135     return (0.5*ResGamma[ires]/TMath::Pi())/((mass-ResMass[ires])*(mass-ResMass[ires])+
02136                                              .25*ResGamma[ires]*ResGamma[ires]);
02137   }
02138 
02139   // NB defined as below, P is not a "probability" but a likelihood that we observe
02140   // a dimuon mass "mass", given mRef, gamma, and massResol. It is what we need for the
02141   // fit which finds the best resolution parameters, though. A definition which is
02142   // more properly a probability is given below (in massProb2()), where the result
02143   // cannot be used to fit resolution parameters because the fit would always prefer
02144   // to set the res parameters to the minimum possible value (best resolution),
02145   // to have a probability close to one of observing any mass.
02146   // -------------------------------------------------------------------------------
02147   // NNBB the following two lines have been replaced with the six following them,
02148   // which provide an improvement of a factor 9 in speed of execution at a
02149   // negligible price in precision.
02150   // ----------------------------------------------------------------------------
02151   // GL->SetParameters(gamma,mRef,mass,massResol);
02152   // P = GL->Integral(mass-5*massResol, mass+5*massResol);
02153 
02154   Int_t np = 100;
02155   double * x = new double[np];
02156   double * w = new double[np];
02157   GL->SetParameters (ResGamma[ires], ResMass[ires], mass, massResol);
02158   GL->CalcGaussLegendreSamplingPoints (np, x, w, 0.1e-15);
02159   P = GL->IntegralFast (np, x, w, ResMass[ires]-10*ResGamma[ires], ResMass[ires]+10*ResGamma[ires]);
02160   delete[] x;
02161   delete[] w;
02162 
02163   // If we are too far away we set P to epsilon and forget about this event
02164   // ----------------------------------------------------------------------
02165   if (P<1.0e-12) {
02166     P = 1.0e-12;
02167     if (debug>9) std::cout << "Mass, gamma , mref, width, P: " << mass
02168                       << " " << ResGamma[ires] << " " << ResMass[ires] << " " << massResol
02169                       << ": used epsilon" << std::endl;
02170     return P;
02171   }
02172 
02173   if (debug>9) std::cout << "Mass, gamma , mref, width, P: " << mass
02174                     << " " << ResGamma[ires] << " " << ResMass[ires] << " " << massResol
02175                     << " " << P << std::endl;
02176   return P;
02177 }
02178 
02179 std::pair<lorentzVector, lorentzVector> MuScleFitUtils::findSimMuFromRes( const edm::Handle<edm::HepMCProduct> & evtMC,
02180                                                                           const edm::Handle<edm::SimTrackContainer> & simTracks )
02181 {
02182   //Loop on simulated tracks
02183   std::pair<lorentzVector, lorentzVector> simMuFromRes;
02184   for( edm::SimTrackContainer::const_iterator simTrack=simTracks->begin(); simTrack!=simTracks->end(); ++simTrack ) {
02185     //Chose muons
02186     if (fabs((*simTrack).type())==13) {
02187       //If tracks from IP than find mother
02188       if ((*simTrack).genpartIndex()>0) {
02189         HepMC::GenParticle* gp = evtMC->GetEvent()->barcode_to_particle ((*simTrack).genpartIndex());
02190         if( gp != 0 ) {
02191 
02192           for (HepMC::GenVertex::particle_iterator mother = gp->production_vertex()->particles_begin(HepMC::ancestors);
02193                mother!=gp->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
02194 
02195             bool fromRes = false;
02196             unsigned int motherPdgId = (*mother)->pdg_id();
02197             for( int ires = 0; ires < 6; ++ires ) {
02198               if( motherPdgId == motherPdgIdArray[ires] && resfind[ires] ) fromRes = true;
02199             }
02200             if( fromRes ) {
02201               if(gp->pdg_id() == 13)
02202                 simMuFromRes.first = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
02203                                                    simTrack->momentum().pz(),simTrack->momentum().e());
02204               else
02205                 simMuFromRes.second = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
02206                                                     simTrack->momentum().pz(),simTrack->momentum().e());
02207             }
02208           }
02209         }
02210         // else LogDebug("MuScleFitUtils") << "WARNING: no matching genParticle found for simTrack" << std::endl;
02211       }
02212     }
02213   }
02214   return simMuFromRes;
02215 }
02216 
02217 std::pair<lorentzVector, lorentzVector> MuScleFitUtils::findGenMuFromRes( const edm::HepMCProduct* evtMC )
02218 {
02219   const HepMC::GenEvent* Evt = evtMC->GetEvent();
02220   std::pair<lorentzVector,lorentzVector> muFromRes;
02221   //Loop on generated particles
02222   for (HepMC::GenEvent::particle_const_iterator part=Evt->particles_begin();
02223        part!=Evt->particles_end(); part++) {
02224     if (fabs((*part)->pdg_id())==13 && (*part)->status()==1) {
02225       bool fromRes = false;
02226       for (HepMC::GenVertex::particle_iterator mother = (*part)->production_vertex()->particles_begin(HepMC::ancestors);
02227            mother != (*part)->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
02228         unsigned int motherPdgId = (*mother)->pdg_id();
02229 
02230         // For sherpa the resonance is not saved. The muons from the resonance can be identified
02231         // by having as mother a muon of status 3.
02232         if( sherpa_ ) {
02233           if( motherPdgId == 13 && (*mother)->status() == 3 ) fromRes = true;
02234         }
02235         else {
02236           for( int ires = 0; ires < 6; ++ires ) {
02237             if( motherPdgId == motherPdgIdArray[ires] && resfind[ires] ) fromRes = true;
02238           }
02239         }
02240       }
02241       if(fromRes){
02242         if((*part)->pdg_id()==13)
02243           //   muFromRes.first = (*part)->momentum();
02244           muFromRes.first = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
02245                                            (*part)->momentum().pz(),(*part)->momentum().e()));
02246         else
02247           muFromRes.second = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
02248                                             (*part)->momentum().pz(),(*part)->momentum().e()));
02249       }
02250     }
02251   }
02252   return muFromRes;
02253 }
02254 
02255 std::pair<lorentzVector, lorentzVector> MuScleFitUtils::findGenMuFromRes( const reco::GenParticleCollection* genParticles)
02256 {
02257   std::pair<lorentzVector,lorentzVector> muFromRes;
02258 
02259   //Loop on generated particles
02260   if( debug>0 ) std::cout << "Starting loop on " << genParticles->size() << " genParticles" << std::endl;
02261   for( reco::GenParticleCollection::const_iterator part=genParticles->begin(); part!=genParticles->end(); ++part ) {
02262     if (fabs(part->pdgId())==13 && part->status()==1) {
02263       bool fromRes = false;
02264       unsigned int motherPdgId = part->mother()->pdgId();
02265       if( debug>0 ) {
02266         std::cout << "Found a muon with mother: " << motherPdgId << std::endl;
02267       }
02268       for( int ires = 0; ires < 6; ++ires ) {
02269         if( motherPdgId == motherPdgIdArray[ires] && resfind[ires] ) fromRes = true;
02270       }
02271       if(fromRes){
02272         if(part->pdgId()==13) {
02273           muFromRes.first = part->p4();
02274           if( debug>0 ) std::cout << "Found a genMuon + : " << muFromRes.first << std::endl;
02275           //      muFromRes.first = (lorentzVector(part->p4().px(),part->p4().py(),
02276           //                                       part->p4().pz(),part->p4().e()));
02277         }
02278         else {
02279           muFromRes.second = part->p4();
02280           if( debug>0 ) std::cout << "Found a genMuon - : " << muFromRes.second << std::endl;
02281           //      muFromRes.second = (lorentzVector(part->p4().px(),part->p4().py(),
02282           //                                        part->p4().pz(),part->p4().e()));
02283         }
02284       }
02285     }
02286   }
02287   return muFromRes;
02288 }