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 #ifndef RooStats_ConfidenceBelt
00036 #include "RooStats/ConfidenceBelt.h"
00037 #endif
00038
00039 #include "RooDataSet.h"
00040 #include "RooDataHist.h"
00041
00042 #include "RooStats/RooStatsUtils.h"
00043
00044 ClassImp(RooStats::ConfidenceBelt) ;
00045
00046 using namespace RooStats;
00047
00048
00049 ConfidenceBelt::ConfidenceBelt() :
00050 TNamed(), fParameterPoints(0)
00051 {
00052
00053 }
00054
00055
00056 ConfidenceBelt::ConfidenceBelt(const char* name) :
00057 TNamed(name,name), fParameterPoints(0)
00058 {
00059
00060 }
00061
00062
00063 ConfidenceBelt::ConfidenceBelt(const char* name, const char* title) :
00064 TNamed(name,title), fParameterPoints(0)
00065 {
00066
00067 }
00068
00069
00070 ConfidenceBelt::ConfidenceBelt(const char* name, RooAbsData& data) :
00071 TNamed(name,name), fParameterPoints((RooAbsData*)data.Clone("PointsToTestForBelt"))
00072 {
00073
00074 }
00075
00076
00077 ConfidenceBelt::ConfidenceBelt(const char* name, const char* title, RooAbsData& data) :
00078 TNamed(name,title), fParameterPoints((RooAbsData*)data.Clone("PointsToTestForBelt"))
00079 {
00080
00081 }
00082
00083
00084
00085
00086 ConfidenceBelt::~ConfidenceBelt()
00087 {
00088
00089
00090 }
00091
00092
00093
00094 Double_t ConfidenceBelt::GetAcceptanceRegionMin(RooArgSet& parameterPoint, Double_t cl, Double_t leftside) {
00095
00096 if(cl>0 || leftside > 0) cout <<"using default cl, leftside for now" <<endl;
00097 return GetAcceptanceRegion(parameterPoint, cl,leftside)->GetLowerLimit();
00098 }
00099
00100
00101 Double_t ConfidenceBelt::GetAcceptanceRegionMax(RooArgSet& parameterPoint, Double_t cl, Double_t leftside) {
00102 if(cl>0 || leftside > 0) cout <<"using default cl, leftside for now" <<endl;
00103 return GetAcceptanceRegion(parameterPoint, cl,leftside)->GetUpperLimit();
00104 }
00105
00106
00107 vector<Double_t> ConfidenceBelt::ConfidenceLevels() const {
00108 vector<Double_t> levels;
00109 return levels;
00110 }
00111
00112
00113 void ConfidenceBelt::AddAcceptanceRegion(RooArgSet& parameterPoint, Int_t dsIndex,
00114 Double_t lower, Double_t upper,
00115 Double_t cl, Double_t leftside){
00116
00117 if(cl>0 || leftside > 0) cout <<"using default cl, leftside for now" <<endl;
00118
00119 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
00120 RooDataHist* hist = dynamic_cast<RooDataHist*>( fParameterPoints );
00121
00122
00123
00124 if( !this->CheckParameters(parameterPoint) )
00125 std::cout << "problem with parameters" << std::endl;
00126
00127 Int_t luIndex = fSamplingSummaryLookup.GetLookupIndex(cl, leftside);
00128
00129 if(luIndex <0 ) {
00130 fSamplingSummaryLookup.Add(cl,leftside);
00131 luIndex = fSamplingSummaryLookup.GetLookupIndex(cl, leftside);
00132 cout << "lookup index = " << luIndex << endl;
00133 }
00134 AcceptanceRegion* thisRegion = new AcceptanceRegion(luIndex, lower, upper);
00135
00136 if( hist ) {
00137
00138
00139
00140
00141
00142 int index = hist->getIndex(parameterPoint);
00143
00144
00145
00146 if((Int_t)fSamplingSummaries.size() <= index) {
00147 fSamplingSummaries.reserve( hist->numEntries() );
00148 fSamplingSummaries.resize( hist->numEntries() );
00149 }
00150
00151
00152 fSamplingSummaries.at(index) = *thisRegion;
00153 }
00154 else if( tree ){
00155
00156 int index = dsIndex;
00157
00158
00159
00160 if((Int_t)fSamplingSummaries.size() <= index){
00161 fSamplingSummaries.reserve( tree->numEntries() );
00162 fSamplingSummaries.resize( tree->numEntries() );
00163 }
00164
00165
00166 fSamplingSummaries.at( index ) = *thisRegion;
00167 }
00168 }
00169
00170
00171 void ConfidenceBelt::AddAcceptanceRegion(RooArgSet& parameterPoint, AcceptanceRegion region,
00172 Double_t cl, Double_t leftside){
00173
00174 if(cl>0 || leftside > 0) cout <<"using default cl, leftside for now" <<endl;
00175
00176 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
00177 RooDataHist* hist = dynamic_cast<RooDataHist*>( fParameterPoints );
00178
00179 if( !this->CheckParameters(parameterPoint) )
00180 std::cout << "problem with parameters" << std::endl;
00181
00182
00183 if( hist ) {
00184
00185
00186
00187
00188
00189 int index = hist->getIndex(parameterPoint);
00190
00191
00192 if((Int_t)fSamplingSummaries.size() < index) fSamplingSummaries.reserve( hist->numEntries() );
00193
00194
00195 fSamplingSummaries.at(index) = region;
00196 }
00197 else if( tree ){
00198 tree->add( parameterPoint );
00199 int index = tree->numEntries() - 1;
00200
00201 if((Int_t)fSamplingSummaries.size() < index) fSamplingSummaries.reserve( tree->numEntries() );
00202
00203
00204 fSamplingSummaries.at( index ) = region;
00205 }
00206 }
00207
00208
00209 AcceptanceRegion* ConfidenceBelt::GetAcceptanceRegion(RooArgSet ¶meterPoint, Double_t cl, Double_t leftside)
00210 {
00211
00212
00213 if(cl>0 || leftside > 0) cout <<"using default cl, leftside for now" <<endl;
00214
00215 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
00216 RooDataHist* hist = dynamic_cast<RooDataHist*>( fParameterPoints );
00217
00218 if( !this->CheckParameters(parameterPoint) ){
00219 std::cout << "problem with parameters" << std::endl;
00220 return 0;
00221 }
00222
00223 if( hist ) {
00224
00225
00226
00227
00228
00229 int index = hist->getIndex(parameterPoint);
00230 return &(fSamplingSummaries.at(index).GetAcceptanceRegion());
00231 }
00232 else if( tree ){
00233
00234
00235 Int_t index = 0;
00236 const RooArgSet* thisPoint = 0;
00237 for(index=0; index<tree->numEntries(); ++index){
00238 thisPoint = tree->get(index);
00239 bool samePoint = true;
00240 TIter it = parameterPoint.createIterator();
00241 RooRealVar *myarg;
00242 while ( samePoint && (myarg = (RooRealVar *)it.Next())) {
00243 if(myarg->getVal() != thisPoint->getRealValue(myarg->GetName()))
00244 samePoint = false;
00245 }
00246 if(samePoint)
00247 break;
00248 }
00249
00250 return &(fSamplingSummaries.at(index).GetAcceptanceRegion());
00251 }
00252 else {
00253 std::cout << "dataset is not initialized properly" << std::endl;
00254 }
00255
00256 return 0;
00257
00258 }
00259
00260
00261 RooArgSet* ConfidenceBelt::GetParameters() const
00262 {
00263
00264 return new RooArgSet(*(fParameterPoints->get()));
00265 }
00266
00267
00268 Bool_t ConfidenceBelt::CheckParameters(RooArgSet ¶meterPoint) const
00269 {
00270
00271 if (parameterPoint.getSize() != fParameterPoints->get()->getSize() ) {
00272 std::cout << "size is wrong, parameters don't match" << std::endl;
00273 return false;
00274 }
00275 if ( ! parameterPoint.equals( *(fParameterPoints->get() ) ) ) {
00276 std::cout << "size is ok, but parameters don't match" << std::endl;
00277 return false;
00278 }
00279 return true;
00280 }