/// /// Perform 1d Prob scan. /// /// - Scan range defined through limit "scan". /// - Will fill the hCL histogram with the 1-CL curve. /// - Start at a scan value that is in the middle of the allowed /// range, preferably a solution, and scan up and down from there. /// - use the "probforce" command line flag to enable force minimum finding /// /// \param fast This will scan each scanpoint only once. /// \param reverse This will scan in reverse direction. /// When using the drag mode, this can sometimes make a difference. /// \return status: 2 = new global minimum found, 1 = error /// int MethodProbScan::scan1d(bool fast, bool reverse) { if ( arg->debug ) cout << "MethodProbScan::scan1d() : starting ... " << endl; nScansDone++; // The "improve" method doesn't need multiple scans. if ( arg->probforce || arg->probimprove ) fast = true; if ( arg->probforce ) scanDisableDragMode = true; // Save parameter values that were active at function call. if ( startPars ) delete startPars; startPars = new RooDataSet("startPars", "startPars", *w->set(parsName)); startPars->add(*w->set(parsName)); // // start scan from global minimum (not always a good idea as we need to set from other places as well) // setParameters(w, parsName, globalMin); // load scan parameter and scan range setLimit(w, scanVar1, "scan"); RooRealVar *par = w->var(scanVar1); assert(par); float min = hCL->GetXaxis()->GetXmin(); float max = hCL->GetXaxis()->GetXmax(); if ( fabs(par->getMin()-min)>1e-6 || fabs(par->getMax()-max)>1e-6 ){ cout << "MethodProbScan::scan1d() : WARNING : Scan range was changed after initScan()" << endl; cout << " was called so the old range will be used." << endl; } if ( arg->verbose ){ cout << "\nProb configuration:" << endl; cout << " combination : " << title << endl; cout << " scan variable : " << scanVar1 << endl; cout << " scan range : " << min << " ... " << max << endl; cout << " scan steps : " << nPoints1d << endl; cout << " fast mode : " << fast << endl; cout << endl; } // Set limit to all parameters. combiner->loadParameterLimits(); // fix scan parameter par->setConstant(true); // j = // 0 : start value -> upper limit // 1 : upper limit -> start value // 2 : start value -> lower limit // 3 : lower limit -> start value float startValue = par->getVal(); bool scanUp; // for the status bar float nTotalSteps = nPoints1d; nTotalSteps *= fast ? 1 : 2; float nStep = 0; float printFreq = nTotalSteps>15 ? 10 : nTotalSteps; // Report on the smallest new minimum we come across while scanning. // Sometimes the scan doesn't find the minimum // that was found before. Warn if this happens. double bestMinOld = chi2minGlobal; double bestMinFoundInScan = 100.; for ( int jj=0; jj<4; jj++ ) { int j = jj; if ( reverse ) switch(jj) { case 0: j = 2; break; case 1: j = 3; break; case 2: j = 0; break; case 3: j = 1; break; } float scanStart, scanStop; switch(j) { case 0: // UP setParameters(w, parsName, startPars->get(0)); scanStart = startValue; scanStop = par->getMax(); scanUp = true; break; case 1: // DOWN scanStart = par->getMax(); scanStop = startValue; scanUp = false; break; case 2: // DOWN setParameters(w, parsName, startPars->get(0)); scanStart = startValue; scanStop = par->getMin(); scanUp = false; break; case 3: // UP scanStart = par->getMin(); scanStop = startValue; scanUp = true; break; } if ( fast && ( j==1 || j==3 ) ) continue; for ( int i=0; i<nPoints1d; i++ ) { float scanvalue; if ( scanUp ) { scanvalue = min + (max-min)*(double)i/(double)nPoints1d + hCL->GetBinWidth(1)/2.; if ( scanvalue < scanStart ) continue; if ( scanvalue > scanStop ) break; } else { scanvalue = max - (max-min)*(double)(i+1)/(double)nPoints1d + hCL->GetBinWidth(1)/2.; if ( scanvalue > scanStart ) continue; if ( scanvalue < scanStop ) break; } // disable drag mode // (the improve method doesn't work with drag mode as parameter run // at their limits) if ( scanDisableDragMode ) setParameters(w, parsName, startPars->get(0)); // set the parameter of interest to the scan point par->setVal(scanvalue); // don't scan in unphysical region if ( scanvalue < par->getMin() || scanvalue > par->getMax() ) continue; // status bar if ( (((int)nStep % (int)(nTotalSteps/printFreq)) == 0)) cout << "MethodProbScan::scan1d() : scanning " << (float)nStep/(float)nTotalSteps*100. << "% \r" << flush; // fit! RooFitResult *fr = 0; if ( arg->probforce ) fr = fitToMinForce(w, combiner->getPdfName()); else if ( arg->probimprove ) fr = fitToMinImprove(w, combiner->getPdfName()); else fr = fitToMinBringBackAngles(w->pdf(pdfName), false, -1); double chi2minScan = fr->minNll(); if ( std::isinf(chi2minScan) ) chi2minScan=1e4; // else the toys in PDF_testConstraint don't work RooSlimFitResult *r = new RooSlimFitResult(fr); // try to save memory by using the slim fit result delete fr; allResults.push_back(r); bestMinFoundInScan = TMath::Min((double)chi2minScan, (double)bestMinFoundInScan); TString warningChi2Neg; if ( chi2minScan < 0 ){ float newChi2minScan = chi2minGlobal + 25.; // 5sigma more than best point warningChi2Neg = "MethodProbScan::scan1d() : WARNING : " + title; warningChi2Neg += TString(Form(" chi2 negative for scan point %i: %f",i,chi2minScan)); warningChi2Neg += " setting to: " + TString(Form("%f",newChi2minScan)); //cout << warningChi2Neg << "\r" << flush; cout << warningChi2Neg << endl; chi2minScan = newChi2minScan; } // If we find a minimum smaller than the old "global" minimum, this means that all // previous 1-CL values are too high. if ( chi2minScan<chi2minGlobal ){ if ( arg->verbose ) cout << "MethodProbScan::scan1d() : WARNING : '" << title << "' new global minimum found! " << " chi2minScan=" << chi2minScan << endl; chi2minGlobal = chi2minScan; // recompute previous 1-CL values for ( int k=1; k<=hCL->GetNbinsX(); k++ ){ hCL->SetBinContent(k, TMath::Prob(hChi2min->GetBinContent(k)-chi2minGlobal, 1)); } } double deltaChi2 = chi2minScan - chi2minGlobal; double oneMinusCL = TMath::Prob(deltaChi2, 1); // Save the 1-CL value and the corresponding fit result. // But only if better than before! if ( hCL->GetBinContent(hCL->FindBin(scanvalue)) <= oneMinusCL ){ hCL->SetBinContent(hCL->FindBin(scanvalue), oneMinusCL); hChi2min->SetBinContent(hCL->FindBin(scanvalue), chi2minScan); int iRes = hCL->FindBin(scanvalue)-1; curveResults[iRes] = r; } nStep++; } } cout << "MethodProbScan::scan1d() : scan done. " << endl; if ( bestMinFoundInScan-bestMinOld > 0.01 ){ cout << "MethodProbScan::scan1d() : WARNING: Scan didn't find similar minimum to what was found before!" << endl; cout << "MethodProbScan::scan1d() : Too strict parameter limits? Too coarse scan steps? Didn't load global minimum?" << endl; cout << "MethodProbScan::scan1d() : chi2 bestMinFoundInScan=" << bestMinFoundInScan << ", bestMinOld=" << bestMinOld << endl; } // attempt to correct for undercoverage if (pvalueCorrectorSet) { for ( int k=1; k<=hCL->GetNbinsX(); k++ ){ double pvalueProb = hCL->GetBinContent(k); pvalueProb = pvalueCorrector->transform(pvalueProb); hCL->SetBinContent(k, pvalueProb); } } setParameters(w, parsName, startPars->get(0)); saveSolutions(); if (arg->confirmsols) confirmSolutions(); if ( (bestMinFoundInScan-bestMinOld)/bestMinOld > 0.01 ) return 1; return 0; }
/// /// Find the global minimum in a more thorough way. /// First fit with external start parameters, then /// for each parameter that starts with "d" or "r" (typically angles and ratios): /// - at upper scan range, rest at start parameters /// - at lower scan range, rest at start parameters /// This amounts to a maximum of 1+2^n fits, where n is the number /// of parameters to be varied. /// /// \param w Workspace holding the pdf. /// \param name Name of the pdf without leading "pdf_". /// \param forceVariables Apply the force method for these variables only. Format /// "var1,var2,var3," (list must end with comma). Default is to apply for all angles, /// all ratios except rD_k3pi and rD_kpi, and the k3pi coherence factor. /// RooFitResult* Utils::fitToMinForce(RooWorkspace *w, TString name, TString forceVariables) { bool debug = true; TString parsName = "par_"+name; TString obsName = "obs_"+name; TString pdfName = "pdf_"+name; RooFitResult *r = 0; int printlevel = -1; RooMsgService::instance().setGlobalKillBelow(ERROR); // save start parameters if ( !w->set(parsName) ){ cout << "MethodProbScan::scan2d() : ERROR : parsName not found: " << parsName << endl; exit(1); } RooDataSet *startPars = new RooDataSet("startParsForce", "startParsForce", *w->set(parsName)); startPars->add(*w->set(parsName)); // set up parameters and ranges RooArgList *varyPars = new RooArgList(); TIterator* it = w->set(parsName)->createIterator(); while ( RooRealVar* p = (RooRealVar*)it->Next() ) { if ( p->isConstant() ) continue; if ( forceVariables=="" && ( false || TString(p->GetName()).BeginsWith("d") ///< use these variables // || TString(p->GetName()).BeginsWith("r") || TString(p->GetName()).BeginsWith("k") || TString(p->GetName()) == "g" ) && ! ( TString(p->GetName()) == "rD_k3pi" ///< don't use these || TString(p->GetName()) == "rD_kpi" // || TString(p->GetName()) == "dD_kpi" || TString(p->GetName()) == "d_dk" || TString(p->GetName()) == "d_dsk" )) { varyPars->add(*p); } else if ( forceVariables.Contains(TString(p->GetName())+",") ) { varyPars->add(*p); } } delete it; int nPars = varyPars->getSize(); if ( debug ) cout << "Utils::fitToMinForce() : nPars = " << nPars << " => " << pow(2.,nPars) << " fits" << endl; if ( debug ) cout << "Utils::fitToMinForce() : varying "; if ( debug ) varyPars->Print(); ////////// r = fitToMinBringBackAngles(w->pdf(pdfName), false, printlevel); ////////// int nErrors = 0; // We define a binary mask where each bit corresponds // to parameter at max or at min. for ( int i=0; i<pow(2.,nPars); i++ ) { if ( debug ) cout << "Utils::fitToMinForce() : fit " << i << " \r" << flush; setParameters(w, parsName, startPars->get(0)); for ( int ip=0; ip<nPars; ip++ ) { RooRealVar *p = (RooRealVar*)varyPars->at(ip); float oldMin = p->getMin(); float oldMax = p->getMax(); setLimit(w, p->GetName(), "force"); if ( i/(int)pow(2.,ip) % 2==0 ) { p->setVal(p->getMin()); } if ( i/(int)pow(2.,ip) % 2==1 ) { p->setVal(p->getMax()); } p->setRange(oldMin, oldMax); } // check if start parameters are sensible, skip if they're not double startParChi2 = getChi2(w->pdf(pdfName)); if ( startParChi2>2000 ){ nErrors += 1; continue; } // refit RooFitResult *r2 = fitToMinBringBackAngles(w->pdf(pdfName), false, printlevel); // In case the initial fit failed, accept the second one. // If both failed, still select the second one and hope the // next fit succeeds. if ( !(r->edm()<1 && r->covQual()==3) ){ delete r; r = r2; } else if ( r2->edm()<1 && r2->covQual()==3 && r2->minNll()<r->minNll() ){ // better minimum found! delete r; r = r2; } else{ delete r2; } } if ( debug ) cout << endl; if ( debug ) cout << "Utils::fitToMinForce() : nErrors = " << nErrors << endl; RooMsgService::instance().setGlobalKillBelow(INFO); // (re)set to best parameters setParameters(w, parsName, r); delete startPars; return r; }
/// /// Perform a 2d Prob scan. /// Scan range defined through limit "scan". /// Fills the hCL2d histogram with the 1-CL curve. /// Saves all encountered fit results to allResults. /// Saves the fit results that make it into the 1-CL curve into curveResults2d. /// Scan strategy: Spiral out! /// int MethodProbScan::scan2d() { if ( arg->debug ) cout << "MethodProbScan::scan2d() : starting ..." << endl; nScansDone++; sanityChecks(); if ( startPars ) delete startPars; // Define whether the 2d contours in hCL are "1D sigma" (ndof=1) or "2D sigma" (ndof=2). // Leave this at 1 for now, as the "2D sigma" contours are computed from hChi2min2d, not hCL. int ndof = 1; // Set up storage for fit results of this particular // scan. This is used for the drag start parameters. // We cannot use the curveResults2d member because that // only holds better results. vector<vector<RooSlimFitResult*> > mycurveResults2d; for ( int i=0; i<nPoints2dx; i++ ){ vector<RooSlimFitResult*> tmp; for ( int j=0; j<nPoints2dy; j++ ) tmp.push_back(0); mycurveResults2d.push_back(tmp); } // store start parameters so we can reset them later startPars = new RooDataSet("startPars", "startPars", *w->set(parsName)); startPars->add(*w->set(parsName)); // // start scan from global minimum (not always a good idea as we need to set from other places as well) // setParameters(w, parsName, globalMin); // Define scan parameters and scan range: RooRealVar *par1 = w->var(scanVar1); RooRealVar *par2 = w->var(scanVar2); // Set limit to all parameters. combiner->loadParameterLimits(); // fix scan parameters par1->setConstant(true); par2->setConstant(true); // Report on the smallest new minimum we come across while scanning. // Sometimes the scan doesn't find the minimum // that was found before. Warn if this happens. double bestMinOld = chi2minGlobal; double bestMinFoundInScan = 100.; // for the status bar int nSteps = 0; float nTotalSteps = nPoints2dx*nPoints2dy; float printFreq = nTotalSteps>100 && !arg->probforce ? 100 : nTotalSteps; ///< number of messages // initialize some control plots gStyle->SetOptTitle(1); TCanvas *cDbg = newNoWarnTCanvas(getUniqueRootName(), Form("DeltaChi2 for 2D scan %i",nScansDone)); cDbg->SetMargin(0.1,0.15,0.1,0.1); float hChi2min2dMin = hChi2min2d->GetMinimum(); bool firstScanDone = hChi2min2dMin<1e5; TH2F *hDbgChi2min2d = histHardCopy(hChi2min2d, firstScanDone); hDbgChi2min2d->SetTitle(Form("#Delta#chi^{2} for scan %i, %s",nScansDone,title.Data())); if ( firstScanDone ) hDbgChi2min2d->GetZaxis()->SetRangeUser(hChi2min2dMin,hChi2min2dMin+25); hDbgChi2min2d->GetXaxis()->SetTitle(par1->GetTitle()); hDbgChi2min2d->GetYaxis()->SetTitle(par2->GetTitle()); hDbgChi2min2d->GetZaxis()->SetTitle("#Delta#chi^{2}"); TH2F *hDbgStart = histHardCopy(hChi2min2d, false); // start coordinates // don't allow the under/overflow bins int iStart = min(hCL2d->GetXaxis()->FindBin(par1->getVal()), hCL2d->GetNbinsX()); int jStart = min(hCL2d->GetYaxis()->FindBin(par2->getVal()), hCL2d->GetNbinsY()); iStart = max(iStart, 1); jStart = max(jStart, 1); hDbgStart->SetBinContent(iStart, jStart, 500.); TMarker *startpointmark = new TMarker(par1->getVal(),par2->getVal(),3); // timer TStopwatch tFit; TStopwatch tSlimResult; TStopwatch tScan; TStopwatch tMemory; // set up the scan spiral int X = 2*nPoints2dx; int Y = 2*nPoints2dy; int x,y,dx,dy; x = y = dx = 0; dy = -1; int t = std::max(X,Y); int maxI = t*t; for ( int spiralstep=0; spiralstep<maxI; spiralstep++ ) { if ((-X/2 <= x) && (x <= X/2) && (-Y/2 <= y) && (y <= Y/2)) { int i = x+iStart; int j = y+jStart; if ( i>0 && i<=nPoints2dx && j>0 && j<=nPoints2dy ) { tScan.Start(false); // status bar if (((int)nSteps % (int)(nTotalSteps/printFreq)) == 0){ cout << Form("MethodProbScan::scan2d() : scanning %3.0f%%", (float)nSteps/(float)nTotalSteps*100.) << " \r" << flush; } nSteps++; // status histogram if ( spiralstep>0 ) hDbgStart->SetBinContent(i, j, 500./*firstScan ? 1. : hChi2min2dMin+36*/); // set start parameters from inner turn of the spiral int xStartPars, yStartPars; computeInnerTurnCoords(iStart, jStart, i, j, xStartPars, yStartPars, 1); RooSlimFitResult *rStartPars = mycurveResults2d[xStartPars-1][yStartPars-1]; if ( rStartPars ) setParameters(w, parsName, rStartPars); // memory management: tMemory.Start(false); // delete old, inner fit results, that we don't need for start parameters anymore // for this we take the second-inner-most turn. int iOld, jOld; bool innerTurnExists = computeInnerTurnCoords(iStart, jStart, i, j, iOld, jOld, 2); if ( innerTurnExists ){ deleteIfNotInCurveResults2d(mycurveResults2d[iOld-1][jOld-1]); mycurveResults2d[iOld-1][jOld-1] = 0; } tMemory.Stop(); // alternative choice for start parameters: always from what we found at function call // setParameters(w, parsName, startPars->get(0)); // set scan point float scanvalue1 = hCL2d->GetXaxis()->GetBinCenter(i); float scanvalue2 = hCL2d->GetYaxis()->GetBinCenter(j); par1->setVal(scanvalue1); par2->setVal(scanvalue2); // fit! tFit.Start(false); RooFitResult *fr; if ( !arg->probforce ) fr = fitToMinBringBackAngles(w->pdf(pdfName), false, -1); else fr = fitToMinForce(w, combiner->getPdfName()); double chi2minScan = fr->minNll(); tFit.Stop(); tSlimResult.Start(false); RooSlimFitResult *r = new RooSlimFitResult(fr); // try to save memory by using the slim fit result tSlimResult.Stop(); delete fr; allResults.push_back(r); bestMinFoundInScan = TMath::Min((double)chi2minScan, (double)bestMinFoundInScan); mycurveResults2d[i-1][j-1] = r; // If we find a new global minumum, this means that all // previous 1-CL values are too high. We'll save the new possible solution, adjust the global // minimum, return a status code, and stop. if ( chi2minScan > -500 && chi2minScan<chi2minGlobal ){ // warn only if there was a significant improvement if ( arg->debug || chi2minScan<chi2minGlobal-1e-2 ){ if ( arg->verbose ) cout << "MethodProbScan::scan2d() : WARNING : '" << title << "' new global minimum found! chi2minGlobal=" << chi2minGlobal << " chi2minScan=" << chi2minScan << endl; } chi2minGlobal = chi2minScan; // recompute previous 1-CL values for ( int k=1; k<=hCL2d->GetNbinsX(); k++ ) for ( int l=1; l<=hCL2d->GetNbinsY(); l++ ){ hCL2d->SetBinContent(k, l, TMath::Prob(hChi2min2d->GetBinContent(k,l)-chi2minGlobal, ndof)); } } double deltaChi2 = chi2minScan - chi2minGlobal; double oneMinusCL = TMath::Prob(deltaChi2, ndof); // Save the 1-CL value. But only if better than before! if ( hCL2d->GetBinContent(i, j) < oneMinusCL ){ hCL2d->SetBinContent(i, j, oneMinusCL); hChi2min2d->SetBinContent(i, j, chi2minScan); hDbgChi2min2d->SetBinContent(i, j, chi2minScan); curveResults2d[i-1][j-1] = r; } // draw/update histograms - doing only every 10th update saves // a lot of time for small combinations if ( ( arg->interactive && ((int)nSteps % 10 == 0) ) || nSteps==nTotalSteps ){ hDbgChi2min2d->Draw("colz"); hDbgStart->Draw("boxsame"); startpointmark->Draw(); cDbg->Update(); cDbg->Modified(); } tScan.Stop(); } } // spiral stuff: if( (x == y) || ((x < 0) && (x == -y)) || ((x > 0) && (x == 1-y))) { t = dx; dx = -dy; dy = t; } x += dx; y += dy; } cout << "MethodProbScan::scan2d() : scan done. " << endl; if ( arg->debug ){ cout << "MethodProbScan::scan2d() : full scan time: "; tScan.Print(); cout << "MethodProbScan::scan2d() : - fitting: "; tFit.Print(); cout << "MethodProbScan::scan2d() : - create RooSlimFitResults: "; tSlimResult.Print(); cout << "MethodProbScan::scan2d() : - memory management: "; tMemory.Print(); } setParameters(w, parsName, startPars->get(0)); saveSolutions2d(); if ( arg->debug ) printLocalMinima(); if ( arg->confirmsols ) confirmSolutions(); // clean all fit results that didn't make it into the final result for ( int i=0; i<allResults.size(); i++ ){ deleteIfNotInCurveResults2d(allResults[i]); } if ( bestMinFoundInScan-bestMinOld > 0.1 ) { cout << "MethodProbScan::scan2d() : WARNING: Scan didn't find minimum that was found before!" << endl; cout << "MethodProbScan::scan2d() : Are you using too strict parameter limits?" << endl; cout << "MethodProbScan::scan2d() : min chi2 found in scan: " << bestMinFoundInScan << ", old min chi2: " << bestMinOld << endl; return 1; } return 0; }