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