Logo ROOT   6.08/07
Reference Guide
Roo1DMomentMorphFunction.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * author: Max Baak (mbaak@cern.ch) *
4  *****************************************************************************/
5 
6 /** \class Roo1DMomentMorphFunction
7  \ingroup Roofit
8 
9 1-dimensional morph function between a list of input functions (varlist) as a function of one input parameter (m).
10 The vector mrefpoints assigns an m-number to each function in the function list.
11 For example: varlist can contain MC histograms (or single numbers) of a reconstructed mass, for certain
12 true Higgs masses indicated in mrefpoints. the input parameter m is the true (continuous) Higgs mass.
13 Morphing can be set to be linear, non-linear or a mixture of the two.
14 */
15 
16 #include "Riostream.h"
17 
19 #include "RooAbsCategory.h"
20 #include "RooRealIntegral.h"
21 #include "RooRealConstant.h"
22 #include "RooRealVar.h"
23 #include "RooFormulaVar.h"
24 #include "RooCustomizer.h"
25 #include "RooAddPdf.h"
26 #include "RooAddition.h"
27 #include "RooMoment.h"
28 #include "RooLinearVar.h"
29 #include "RooChangeTracker.h"
30 
31 #include "TMath.h"
32 
33 using namespace std;
34 
36 
37 
38 ////////////////////////////////////////////////////////////////////////////////
39 /// Default constructor.
40 
41  Roo1DMomentMorphFunction::Roo1DMomentMorphFunction() : _mref(0), _frac(0), _M(0), _setting(Linear)
42 {
43  _varItr = _varList.createIterator() ;
44 }
45 
46 
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 /// Constructor.
50 /// \param[in] name
51 /// \param[in] title
52 /// \param[in] _m
53 /// \param[in] varList
54 /// \param[in] mrefpoints
55 /// \param[in] setting
56 
58  RooAbsReal& _m,
59  const RooArgList& varList,
60  const TVectorD& mrefpoints,
61  const Setting& setting) :
62  RooAbsReal(name,title),
63  m("m","m",this,_m),
64  _varList("varList","List of variables",this),
65  _setting(setting)
66 {
67  // observables
68  TIterator* varItr = varList.createIterator() ;
69  RooAbsArg* var ;
70  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
71  if (!dynamic_cast<RooAbsReal*>(var)) {
72  coutE(InputArguments) << "Roo1DMomentMorphFunction::ctor(" << GetName() << ") ERROR: variable " << var->GetName() << " is not of type RooAbsReal" << endl ;
73  throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal") ;
74  }
75  _varList.add(*var) ;
76  }
77  delete varItr ;
78 
79  _mref = new TVectorD(mrefpoints);
80  _frac = 0 ;
82 
83  // initialization
84  initialize();
85 }
86 
87 
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// Copy constructor
91 /// \param[in] other
92 /// \param[in] name
93 
95  RooAbsReal(other,name),
96  m("m",this,other.m),
97  _varList("varList",this,other._varList),
98  _setting(other._setting)
99 {
100  _mref = new TVectorD(*other._mref) ;
101  _frac = 0 ;
103 
104  // initialization
105  initialize();
106 }
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 /// Destructor.
110 
112 {
113  if (_mref) delete _mref;
114  if (_frac) delete _frac;
115  if (_varItr) delete _varItr;
116  if (_M) delete _M;
117 }
118 
119 
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 
123 
125 {
126  Int_t nVar = _varList.getSize();
127 
128  // other quantities needed
129  if (nVar!=_mref->GetNrows()) {
130  coutE(InputArguments) << "Roo1DMomentMorphFunction::initialize(" << GetName() << ") ERROR: nVar != nRefPoints" << endl ;
131  assert(0) ;
132  }
133 
134  _frac = new TVectorD(nVar);
135 
136  TVectorD* dm = new TVectorD(nVar);
137  _M = new TMatrixD(nVar,nVar);
138 
139  // transformation matrix for non-linear extrapolation, needed in evaluate()
140  TMatrixD M(nVar,nVar);
141  for (Int_t i=0; i<_mref->GetNrows(); ++i) {
142  (*dm)[i] = (*_mref)[i]-(*_mref)[0];
143  M(i,0) = 1.;
144  if (i>0) M(0,i) = 0.;
145  }
146  for (Int_t i=1; i<_mref->GetNrows(); ++i) {
147  for (Int_t j=1; j<_mref->GetNrows(); ++j) {
148  M(i,j) = TMath::Power((*dm)[i],(double)j);
149  }
150  }
151  (*_M) = M.Invert();
152 
153  delete dm ;
154 }
155 
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 
160 {
161  calculateFractions(); // this sets _frac vector, based on function of m
162 
163  _varItr->Reset() ;
164 
165  Double_t ret(0);
166  RooAbsReal* var(0) ;
167  for (Int_t i=0; (var = (RooAbsReal*)_varItr->Next()); ++i) {
168  ret += (*_frac)(i) * var->getVal();
169  }
170 
171  return ret ;
172 }
173 
174 
175 ////////////////////////////////////////////////////////////////////////////////
176 
178 {
179  Int_t nVar = _varList.getSize();
180 
181  Double_t dm = m - (*_mref)[0];
182 
183  // fully non-linear
184  double sumposfrac=0.;
185  for (Int_t i=0; i<nVar; ++i) {
186  double ffrac=0.;
187  for (Int_t j=0; j<nVar; ++j) { ffrac += (*_M)(j,i) * (j==0?1.:TMath::Power(dm,(double)j)); }
188  if (ffrac>=0) sumposfrac+=ffrac;
189  (*_frac)(i) = ffrac;
190  }
191 
192  // various mode settings
193  int imin = idxmin(m);
194  int imax = idxmax(m);
195  double mfrac = (m-(*_mref)[imin])/((*_mref)[imax]-(*_mref)[imin]);
196  switch (_setting) {
197  case NonLinear:
198  // default already set above
199  break;
200  case Linear:
201  for (Int_t i=0; i<nVar; ++i)
202  (*_frac)(i) = 0.;
203  if (imax>imin) { // m in between mmin and mmax
204  (*_frac)(imin) = (1.-mfrac);
205  (*_frac)(imax) = (mfrac);
206  } else if (imax==imin) { // m outside mmin and mmax
207  (*_frac)(imin) = (1.);
208  }
209  break;
211  for (Int_t i=0; i<nVar; ++i)
212  (*_frac)(i) = (0.);
213  if (imax>imin) { // m in between mmin and mmax
214  (*_frac)(imin) = (1.-mfrac);
215  (*_frac)(imax) = (mfrac);
216  } else if (imax==imin) { // m outside mmin and mmax
217  (*_frac)(imin) = (1.);
218  }
219  break;
221  for (Int_t i=0; i<nVar; ++i) {
222  if ((*_frac)(i)<0) (*_frac)(i)=(0.);
223  (*_frac)(i) = (*_frac)(i)/sumposfrac;
224  }
225  break;
226  }
227 }
228 
229 ////////////////////////////////////////////////////////////////////////////////
230 
231 int Roo1DMomentMorphFunction::idxmin(const double& mval) const
232 {
233  int imin(0);
234  Int_t nVar = _varList.getSize();
235  double mmin=-DBL_MAX;
236  for (Int_t i=0; i<nVar; ++i)
237  if ( (*_mref)[i]>mmin && (*_mref)[i]<=mval ) { mmin=(*_mref)[i]; imin=i; }
238  return imin;
239 }
240 
241 
242 ////////////////////////////////////////////////////////////////////////////////
243 
244 int Roo1DMomentMorphFunction::idxmax(const double& mval) const
245 {
246  int imax(0);
247  Int_t nVar = _varList.getSize();
248  double mmax=DBL_MAX;
249  for (Int_t i=0; i<nVar; ++i)
250  if ( (*_mref)[i]<mmax && (*_mref)[i]>=mval ) { mmax=(*_mref)[i]; imax=i; }
251  return imax;
252 }
TMatrixD * _M
do not persist
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
virtual void Reset()=0
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
int idxmin(const double &m) const
virtual ~Roo1DMomentMorphFunction()
Destructor.
Int_t GetNrows() const
Definition: TVectorT.h:81
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Iterator abstract base class.
Definition: TIterator.h:32
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
int idxmax(const double &m) const
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
Definition: TMatrixT.cxx:1396
Roo1DMomentMorphFunction()
Default constructor.
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
Int_t getSize() const
1-dimensional morph function between a list of input functions (varlist) as a function of one input p...
TMarker * m
Definition: textangle.C:8
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
static std::shared_ptr< std::function< double(double)> > Linear
Definition: NeuralNet.icc:53
virtual TObject * Next()=0
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
char name[80]
Definition: TGX11.cxx:109
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts &#39;var&#39; into set and registers &#39;var&#39; as server to owner with...