70 ClassImp(RooStats::NeymanConstruction); ;
72 using namespace RooFit;
73 using namespace RooStats;
80 NeymanConstruction::NeymanConstruction(RooAbsData& data, ModelConfig& model):
88 fAdaptiveSampling(false),
89 fAdditionalNToysFactor(1.),
90 fSaveBeltToFile(false),
105 NeymanConstruction::~NeymanConstruction() {
112 PointSetInterval* NeymanConstruction::GetInterval()
const {
117 oocoutI(f,Contents) <<
"NeymanConstruction saving ConfidenceBelt to file SamplingDistributions.root" << endl;
118 f =
new TFile(
"SamplingDistributions.root",
"recreate");
126 RooArgSet* fPOI =
new RooArgSet(*fModel.GetParametersOfInterest());
128 RooDataSet* pointsInInterval =
new RooDataSet(
"pointsInInterval",
129 "points in interval",
130 *(fPointsToTest->get(0)) );
133 for(Int_t i=0; i<fPointsToTest->numEntries(); ++i){
135 point = (RooArgSet*) fPointsToTest->get(i);
141 fTestStatSampler->SetParametersForTestStat(*fPOI);
144 Double_t thisTestStatistic = fTestStatSampler->EvaluateTestStatistic(fData, *fPOI );
156 SamplingDistribution* samplingDist=0;
158 Double_t upperEdgeOfAcceptance, upperEdgeMinusSigma, upperEdgePlusSigma;
159 Double_t lowerEdgeOfAcceptance, lowerEdgeMinusSigma, lowerEdgePlusSigma;
160 Int_t additionalMC=0;
172 Int_t totalMC = (Int_t) (2./fSize/TMath::Min(fLeftSideFraction,1.-fLeftSideFraction));
173 if(fLeftSideFraction==0. || fLeftSideFraction ==1.){
174 totalMC = (Int_t) (2./fSize);
177 Double_t tmc = Double_t(totalMC)*fAdditionalNToysFactor;
178 totalMC = (Int_t) tmc;
180 ToyMCSampler* toyMCSampler =
dynamic_cast<ToyMCSampler*
>(fTestStatSampler);
181 if(fAdaptiveSampling && toyMCSampler) {
188 additionalMC = 2*totalMC;
190 toyMCSampler->AppendSamplingDistribution(*point,
194 oocoutE((TObject*)0,Eval) <<
"Neyman Construction: error generating sampling distribution" << endl;
197 totalMC=samplingDist->GetSize();
203 upperEdgeOfAcceptance =
204 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) ,
205 sigma, upperEdgePlusSigma);
207 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) ,
208 sigma, upperEdgeMinusSigma);
211 lowerEdgeOfAcceptance =
212 samplingDist->InverseCDF( fLeftSideFraction * fSize ,
213 sigma, lowerEdgePlusSigma);
215 samplingDist->InverseCDF( fLeftSideFraction * fSize ,
216 sigma, lowerEdgeMinusSigma);
218 ooccoutD(samplingDist,Eval) <<
"NeymanConstruction: "
219 <<
"total MC = " << totalMC
220 <<
" this test stat = " << thisTestStatistic << endl
221 <<
" upper edge -1sigma = " << upperEdgeMinusSigma
222 <<
", upperEdge = "<<upperEdgeOfAcceptance
223 <<
", upper edge +1sigma = " << upperEdgePlusSigma << endl
224 <<
" lower edge -1sigma = " << lowerEdgeMinusSigma
225 <<
", lowerEdge = "<<lowerEdgeOfAcceptance
226 <<
", lower edge +1sigma = " << lowerEdgePlusSigma << endl;
228 (thisTestStatistic <= upperEdgeOfAcceptance &&
229 thisTestStatistic > upperEdgeMinusSigma)
230 || (thisTestStatistic >= upperEdgeOfAcceptance &&
231 thisTestStatistic < upperEdgePlusSigma)
232 || (thisTestStatistic <= lowerEdgeOfAcceptance &&
233 thisTestStatistic > lowerEdgeMinusSigma)
234 || (thisTestStatistic >= lowerEdgeOfAcceptance &&
235 thisTestStatistic < lowerEdgePlusSigma)
236 ) && (totalMC < 100./fSize)
241 samplingDist = fTestStatSampler->GetSamplingDistribution(*point);
243 oocoutE((TObject*)0,Eval) <<
"Neyman Construction: error generating sampling distribution" << endl;
247 lowerEdgeOfAcceptance =
248 samplingDist->InverseCDF( fLeftSideFraction * fSize );
249 upperEdgeOfAcceptance =
250 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) );
254 if(fConfBelt && fCreateBelt){
256 fConfBelt->AddAcceptanceRegion(*point, i,
257 lowerEdgeOfAcceptance,
258 upperEdgeOfAcceptance);
262 TIter itr = point->createIterator();
264 ooccoutP(samplingDist,Eval) <<
"NeymanConstruction: Prog: "<< i+1<<
"/"<<fPointsToTest->numEntries()
265 <<
" total MC = " << samplingDist->GetSize()
266 <<
" this test stat = " << thisTestStatistic << endl;
267 ooccoutP(samplingDist,Eval) <<
" ";
268 while ((myarg = (RooRealVar *)itr.Next())) {
269 ooccoutP(samplingDist,Eval) << myarg->GetName() <<
"=" << myarg->getVal() <<
" ";
271 ooccoutP(samplingDist,Eval) <<
"[" << lowerEdgeOfAcceptance <<
", "
272 << upperEdgeOfAcceptance <<
"] " <<
" in interval = " <<
273 (thisTestStatistic >= lowerEdgeOfAcceptance && thisTestStatistic <= upperEdgeOfAcceptance)
277 if(thisTestStatistic >= lowerEdgeOfAcceptance && thisTestStatistic <= upperEdgeOfAcceptance) {
280 pointsInInterval->add(*point);
286 samplingDist->Write();
287 string tmpName =
"hist_";
288 tmpName+=samplingDist->GetName();
289 TH1F* h =
new TH1F(tmpName.c_str(),
"",500,0.,5.);
290 for(
int ii=0; ii<samplingDist->GetSize(); ++ii){
291 h->Fill(samplingDist->GetSamplingDistribution().at(ii) );
300 oocoutI(pointsInInterval,Eval) << npass <<
" points in interval" << endl;
303 PointSetInterval* interval
304 =
new PointSetInterval(
"ClassicalConfidenceInterval", *pointsInInterval);