CMS 3D CMS Logo

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