Logo ROOT   6.18/05
Reference Guide
RooJeffreysPrior.cxx
Go to the documentation of this file.
1/** \class RooJeffreysPrior
2\ingroup Roofit
3
4RooJeffreysPrior
5**/
6
7
8#include "RooFit.h"
9
10#include "Riostream.h"
11#include "Riostream.h"
12#include <math.h>
13
14#include "RooJeffreysPrior.h"
15#include "RooAbsReal.h"
16#include "RooAbsPdf.h"
17#include "RooErrorHandler.h"
18#include "RooArgSet.h"
19#include "RooMsgService.h"
20#include "RooFitResult.h"
21#include "TMatrixDSym.h"
22#include "RooDataHist.h"
23#include "RooFitResult.h"
24#include "RooNumIntConfig.h"
25#include "RooRealVar.h"
26
27#include "TError.h"
28
29using namespace std;
30
32
33using namespace RooFit;
34
35////////////////////////////////////////////////////////////////////////////////
36///_obsSet("!obsSet","obs-side variation",this),
37
38RooJeffreysPrior::RooJeffreysPrior(const char* name, const char* title,
39 RooAbsPdf& nominal,
40 const RooArgList& paramSet,
41 const RooArgList& obsSet) :
42 RooAbsPdf(name, title),
43 _nominal("nominal","nominal",this,nominal,kFALSE,kFALSE),
44 _obsSet("!obsSet","obs-side variation",this,kFALSE,kFALSE),
45 _paramSet("!paramSet","high-side variation",this)
46{
49
50
51 TIterator* inputIter1 = obsSet.createIterator() ;
52 RooAbsArg* comp ;
53 while((comp = (RooAbsArg*)inputIter1->Next())) {
54 if (!dynamic_cast<RooAbsReal*>(comp)) {
55 coutE(InputArguments) << "RooJeffreysPrior::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
56 << " in first list is not of type RooAbsReal" << endl ;
58 }
59 _obsSet.add(*comp) ;
60 // if (takeOwnership) {
61 // _ownedList.addOwned(*comp) ;
62 // }
63 }
64 delete inputIter1 ;
65
66 TIterator* inputIter3 = paramSet.createIterator() ;
67 while((comp = (RooAbsArg*)inputIter3->Next())) {
68 if (!dynamic_cast<RooAbsReal*>(comp)) {
69 coutE(InputArguments) << "RooJeffreysPrior::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
70 << " in first list is not of type RooAbsReal" << endl ;
72 }
73 _paramSet.add(*comp) ;
74 // if (takeOwnership) {
75 // _ownedList.addOwned(*comp) ;
76 // }
77 }
78 delete inputIter3 ;
79
80 // use a different integrator by default.
81 if(paramSet.getSize()==1)
82 this->specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
83
84}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Copy constructor
88
90 RooAbsPdf(other, name),
91 _nominal("!nominal",this,other._nominal),
92 _obsSet("!obsSet",this,other._obsSet),
93 _paramSet("!paramSet",this,other._paramSet)
94{
97
98 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// Default constructor
103
105{
106 _obsIter = NULL;
107 _paramIter = NULL;
108
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Destructor
113
115{
116 if (_obsIter) delete _obsIter ;
117 if (_paramIter) delete _paramIter ;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// Calculate and return current value of self
122
124{
127 // create Asimov dataset
128 // _paramSet.Print("v");
129 // return sqrt(_paramSet.getRealValue("mu"));
131 /*
132 cout << "_______________"<<endl;
133 _paramSet.Print("v");
134 _nominal->getVariables()->Print("v");
135 cout << "_______________"<<endl;
136 */
137 RooDataHist* data = ((RooAbsPdf&)(_nominal.arg())).generateBinned(_obsSet,ExpectedData());
138 // data->Print("v");
139 //RooFitResult* res = _nominal->fitTo(*data, Save(),PrintLevel(-1),Minos(kFALSE),SumW2Error(kTRUE));
140 RooFitResult* res = ((RooAbsPdf&)(_nominal.arg())).fitTo(*data, Save(),PrintLevel(-1),Minos(kFALSE),SumW2Error(kFALSE));
141 TMatrixDSym cov = res->covarianceMatrix();
142 cov.Invert();
143 double ret = sqrt(cov.Determinant());
144
145 /*
146 // for 1 parameter can avoid making TMatrix etc.
147 // but number of params may be > 1 with others held constant
148 if(_paramSet.getSize()==1){
149 RooRealVar* var = (RooRealVar*) _paramSet.first();
150 // also, the _paramSet proxy one does not pick up a different value
151 cout << "eval at "<< ret << " " << 1/(var->getError()) << endl;
152 // need to get the actual variable instance out of the pdf like below
153 var = (RooRealVar*) _nominal->getVariables()->find(var->GetName());
154 cout << "eval at "<< ret << " " << 1/(var->getError()) << endl;
155 }
156 */
157
158 // res->Print();
159 delete data;
160 delete res;
162
163 // cout << "eval at "<< ret << endl;
164 // _paramSet.Print("v");
165 return ret;
166
167}
168
169////////////////////////////////////////////////////////////////////////////////
170/// if (matchArgs(allVars,analVars,x)) return 1 ;
171/// if (matchArgs(allVars,analVars,mean)) return 2 ;
172/// return 1;
173
174Int_t RooJeffreysPrior::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
175{
176 return 0 ;
177}
178
179////////////////////////////////////////////////////////////////////////////////
180
181Double_t RooJeffreysPrior::analyticalIntegral(Int_t code, const char* /*rangeName*/) const
182{
183 R__ASSERT(code==1 );
184 //cout << "evaluating analytic integral" << endl;
185 return 1.;
186}
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:70
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2034
Int_t getSize() const
TIterator * createIterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
TIterator-style iteration over contained elements.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)
Set value by specifying the name of the desired state If printError is set, a message will be printed...
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
static void softAbort()
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
RooJeffreysPrior.
virtual ~RooJeffreysPrior()
Destructor.
Double_t evaluate() const
Iterator over lowSet.
RooRealProxy _nominal
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
TIterator * _paramIter
TIterator * _obsIter
Iterator over paramSet.
RooJeffreysPrior()
Default constructor.
RooListProxy _obsSet
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
if (matchArgs(allVars,analVars,x)) return 1 ; if (matchArgs(allVars,analVars,mean)) return 2 ; return...
RooListProxy _paramSet
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
RooFit::MsgLevel globalKillBelow() const
RooCategory & method1D()
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual Double_t Determinant() const
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Template specialisation used in RooAbsArg:
RooCmdArg SumW2Error(Bool_t flag)
@ InputArguments
Definition: RooGlobalFunc.h:58
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg ExpectedData(Bool_t flag=kTRUE)
RooCmdArg Minos(Bool_t flag=kTRUE)