00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045 #ifndef RooStats_NeymanConstruction
00046 #include "RooStats/NeymanConstruction.h"
00047 #endif
00048
00049 #ifndef RooStats_RooStatsUtils
00050 #include "RooStats/RooStatsUtils.h"
00051 #endif
00052
00053 #ifndef RooStats_PointSetInterval
00054 #include "RooStats/PointSetInterval.h"
00055 #endif
00056
00057 #include "RooStats/SamplingDistribution.h"
00058 #include "RooStats/ToyMCSampler.h"
00059 #include "RooStats/ModelConfig.h"
00060
00061 #include "RooMsgService.h"
00062 #include "RooGlobalFunc.h"
00063
00064 #include "RooDataSet.h"
00065 #include "TFile.h"
00066 #include "TTree.h"
00067 #include "TMath.h"
00068 #include "TH1F.h"
00069
00070
00071
00072 ClassImp(RooStats::NeymanConstruction) ;
00073
00074 using namespace RooFit;
00075 using namespace RooStats;
00076
00077
00078
00079 NeymanConstruction::NeymanConstruction(RooAbsData& data, ModelConfig& model):
00080 fSize(0.05),
00081 fData(data),
00082 fModel(model),
00083 fTestStatSampler(0),
00084 fPointsToTest(0),
00085 fLeftSideFraction(0),
00086 fConfBelt(0),
00087 fAdaptiveSampling(false),
00088 fAdditionalNToysFactor(1.),
00089 fSaveBeltToFile(false),
00090 fCreateBelt(false)
00091
00092 {
00093
00094
00095
00096
00097
00098 }
00099
00100
00101 NeymanConstruction::~NeymanConstruction() {
00102
00103
00104
00105 }
00106
00107
00108 PointSetInterval* NeymanConstruction::GetInterval() const {
00109
00110
00111
00112
00113 TFile* f=0;
00114 if(fSaveBeltToFile){
00115 oocoutI(f,Contents) << "NeymanConstruction saving ConfidenceBelt to file SamplingDistributions.root" << endl;
00116 f = new TFile("SamplingDistributions.root","recreate");
00117 }
00118
00119 if(!&fData) {
00120 oocoutF((TObject*)0,Contents) << "Neyman Construction: data is not set, can't get interval" << endl;
00121 return 0;
00122 }
00123 Int_t npass = 0;
00124 RooArgSet* point;
00125
00126
00127
00128 RooArgSet* fPOI = new RooArgSet(*fModel.GetParametersOfInterest());
00129
00130 RooDataSet* pointsInInterval = new RooDataSet("pointsInInterval",
00131 "points in interval",
00132 *(fPointsToTest->get(0)) );
00133
00134
00135 for(Int_t i=0; i<fPointsToTest->numEntries(); ++i){
00136
00137 point = (RooArgSet*) fPointsToTest->get(i);
00138
00139
00140 *fPOI = *point;
00141
00142
00143 fTestStatSampler->SetParametersForTestStat(*fPOI);
00144
00145
00146 Double_t thisTestStatistic = fTestStatSampler->EvaluateTestStatistic(fData, *fPOI );
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158 SamplingDistribution* samplingDist=0;
00159 Double_t sigma;
00160 Double_t upperEdgeOfAcceptance, upperEdgeMinusSigma, upperEdgePlusSigma;
00161 Double_t lowerEdgeOfAcceptance, lowerEdgeMinusSigma, lowerEdgePlusSigma;
00162 Int_t additionalMC=0;
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174 Int_t totalMC = (Int_t) (2./fSize/TMath::Min(fLeftSideFraction,1.-fLeftSideFraction));
00175 if(fLeftSideFraction==0. || fLeftSideFraction ==1.){
00176 totalMC = (Int_t) (2./fSize);
00177 }
00178
00179 Double_t tmc = Double_t(totalMC)*fAdditionalNToysFactor;
00180 totalMC = (Int_t) tmc;
00181
00182 ToyMCSampler* toyMCSampler = dynamic_cast<ToyMCSampler*>(fTestStatSampler);
00183 if(fAdaptiveSampling && toyMCSampler) {
00184 do{
00185
00186
00187
00188
00189
00190 additionalMC = 2*totalMC;
00191 samplingDist =
00192 toyMCSampler->AppendSamplingDistribution(*point,
00193 samplingDist,
00194 additionalMC);
00195 totalMC=samplingDist->GetSize();
00196
00197
00198
00199
00200 sigma = 1;
00201 upperEdgeOfAcceptance =
00202 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) ,
00203 sigma, upperEdgePlusSigma);
00204 sigma = -1;
00205 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) ,
00206 sigma, upperEdgeMinusSigma);
00207
00208 sigma = 1;
00209 lowerEdgeOfAcceptance =
00210 samplingDist->InverseCDF( fLeftSideFraction * fSize ,
00211 sigma, lowerEdgePlusSigma);
00212 sigma = -1;
00213 samplingDist->InverseCDF( fLeftSideFraction * fSize ,
00214 sigma, lowerEdgeMinusSigma);
00215
00216 ooccoutD(samplingDist,Eval) << "NeymanConstruction: "
00217 << "total MC = " << totalMC
00218 << " this test stat = " << thisTestStatistic << endl
00219 << " upper edge -1sigma = " << upperEdgeMinusSigma
00220 << ", upperEdge = "<<upperEdgeOfAcceptance
00221 << ", upper edge +1sigma = " << upperEdgePlusSigma << endl
00222 << " lower edge -1sigma = " << lowerEdgeMinusSigma
00223 << ", lowerEdge = "<<lowerEdgeOfAcceptance
00224 << ", lower edge +1sigma = " << lowerEdgePlusSigma << endl;
00225 } while((
00226 (thisTestStatistic <= upperEdgeOfAcceptance &&
00227 thisTestStatistic > upperEdgeMinusSigma)
00228 || (thisTestStatistic >= upperEdgeOfAcceptance &&
00229 thisTestStatistic < upperEdgePlusSigma)
00230 || (thisTestStatistic <= lowerEdgeOfAcceptance &&
00231 thisTestStatistic > lowerEdgeMinusSigma)
00232 || (thisTestStatistic >= lowerEdgeOfAcceptance &&
00233 thisTestStatistic < lowerEdgePlusSigma)
00234 ) && (totalMC < 100./fSize)
00235 ) ;
00236 } else {
00237
00238
00239 samplingDist = fTestStatSampler->GetSamplingDistribution(*point);
00240
00241 lowerEdgeOfAcceptance =
00242 samplingDist->InverseCDF( fLeftSideFraction * fSize );
00243 upperEdgeOfAcceptance =
00244 samplingDist->InverseCDF( 1. - ((1.-fLeftSideFraction) * fSize) );
00245 }
00246
00247
00248 if(fConfBelt && fCreateBelt){
00249
00250 fConfBelt->AddAcceptanceRegion(*point, i,
00251 lowerEdgeOfAcceptance,
00252 upperEdgeOfAcceptance);
00253 }
00254
00255
00256 TIter itr = point->createIterator();
00257 RooRealVar* myarg;
00258 ooccoutP(samplingDist,Eval) << "NeymanConstruction: Prog: "<< i+1<<"/"<<fPointsToTest->numEntries()
00259 << " total MC = " << samplingDist->GetSize()
00260 << " this test stat = " << thisTestStatistic << endl;
00261 ooccoutP(samplingDist,Eval) << " ";
00262 while ((myarg = (RooRealVar *)itr.Next())) {
00263 ooccoutP(samplingDist,Eval) << myarg->GetName() << "=" << myarg->getVal() << " ";
00264 }
00265 ooccoutP(samplingDist,Eval) << "[" << lowerEdgeOfAcceptance << ", "
00266 << upperEdgeOfAcceptance << "] " << " in interval = " <<
00267 (thisTestStatistic >= lowerEdgeOfAcceptance && thisTestStatistic <= upperEdgeOfAcceptance)
00268 << endl << endl;
00269
00270
00271 if(thisTestStatistic >= lowerEdgeOfAcceptance && thisTestStatistic <= upperEdgeOfAcceptance) {
00272
00273
00274 pointsInInterval->add(*point);
00275 ++npass;
00276 }
00277
00278 if(fSaveBeltToFile){
00279
00280 samplingDist->Write();
00281 string tmpName = "hist_";
00282 tmpName+=samplingDist->GetName();
00283 TH1F* h = new TH1F(tmpName.c_str(),"",500,0.,5.);
00284 for(int ii=0; ii<samplingDist->GetSize(); ++ii){
00285 h->Fill(samplingDist->GetSamplingDistribution().at(ii) );
00286 }
00287 h->Write();
00288 delete h;
00289 }
00290
00291 delete samplingDist;
00292
00293 }
00294 oocoutI(pointsInInterval,Eval) << npass << " points in interval" << endl;
00295
00296
00297 PointSetInterval* interval
00298 = new PointSetInterval("ClassicalConfidenceInterval", *pointsInInterval);
00299
00300
00301 if(fSaveBeltToFile){
00302
00303 fConfBelt->Write();
00304
00305 f->Close();
00306 }
00307
00308 delete f;
00309
00310 return interval;
00311 }
00312
00313