Logo ROOT   6.10/09
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 /// Default constructor.
39 
40  Roo1DMomentMorphFunction::Roo1DMomentMorphFunction() : _mref(0), _frac(0), _M(0), _setting(Linear)
41 {
42  _varItr = _varList.createIterator() ;
43 }
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Constructor.
47 /// \param[in] name
48 /// \param[in] title
49 /// \param[in] _m
50 /// \param[in] varList
51 /// \param[in] mrefpoints
52 /// \param[in] setting
53 
55  RooAbsReal& _m,
56  const RooArgList& varList,
57  const TVectorD& mrefpoints,
58  const Setting& setting) :
59  RooAbsReal(name,title),
60  m("m","m",this,_m),
61  _varList("varList","List of variables",this),
62  _setting(setting)
63 {
64  // observables
65  TIterator* varItr = varList.createIterator() ;
66  RooAbsArg* var ;
67  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
68  if (!dynamic_cast<RooAbsReal*>(var)) {
69  coutE(InputArguments) << "Roo1DMomentMorphFunction::ctor(" << GetName() << ") ERROR: variable " << var->GetName() << " is not of type RooAbsReal" << endl ;
70  throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal") ;
71  }
72  _varList.add(*var) ;
73  }
74  delete varItr ;
75 
76  _mref = new TVectorD(mrefpoints);
77  _frac = 0 ;
79 
80  // initialization
81  initialize();
82 }
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 /// Copy constructor
86 /// \param[in] other
87 /// \param[in] name
88 
90  RooAbsReal(other,name),
91  m("m",this,other.m),
92  _varList("varList",this,other._varList),
93  _setting(other._setting)
94 {
95  _mref = new TVectorD(*other._mref) ;
96  _frac = 0 ;
98 
99  // initialization
100  initialize();
101 }
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 /// Destructor.
105 
107 {
108  if (_mref) delete _mref;
109  if (_frac) delete _frac;
110  if (_varItr) delete _varItr;
111  if (_M) delete _M;
112 }
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 
117 {
118  Int_t nVar = _varList.getSize();
119 
120  // other quantities needed
121  if (nVar!=_mref->GetNrows()) {
122  coutE(InputArguments) << "Roo1DMomentMorphFunction::initialize(" << GetName() << ") ERROR: nVar != nRefPoints" << endl ;
123  assert(0) ;
124  }
125 
126  _frac = new TVectorD(nVar);
127 
128  TVectorD* dm = new TVectorD(nVar);
129  _M = new TMatrixD(nVar,nVar);
130 
131  // transformation matrix for non-linear extrapolation, needed in evaluate()
132  TMatrixD M(nVar,nVar);
133  for (Int_t i=0; i<_mref->GetNrows(); ++i) {
134  (*dm)[i] = (*_mref)[i]-(*_mref)[0];
135  M(i,0) = 1.;
136  if (i>0) M(0,i) = 0.;
137  }
138  for (Int_t i=1; i<_mref->GetNrows(); ++i) {
139  for (Int_t j=1; j<_mref->GetNrows(); ++j) {
140  M(i,j) = TMath::Power((*dm)[i],(double)j);
141  }
142  }
143  (*_M) = M.Invert();
144 
145  delete dm ;
146 }
147 
148 ////////////////////////////////////////////////////////////////////////////////
149 
151 {
152  calculateFractions(); // this sets _frac vector, based on function of m
153 
154  _varItr->Reset() ;
155 
156  Double_t ret(0);
157  RooAbsReal* var(0) ;
158  for (Int_t i=0; (var = (RooAbsReal*)_varItr->Next()); ++i) {
159  ret += (*_frac)(i) * var->getVal();
160  }
161 
162  return ret ;
163 }
164 
165 ////////////////////////////////////////////////////////////////////////////////
166 
168 {
169  Int_t nVar = _varList.getSize();
170 
171  Double_t dm = m - (*_mref)[0];
172 
173  // fully non-linear
174  double sumposfrac=0.;
175  for (Int_t i=0; i<nVar; ++i) {
176  double ffrac=0.;
177  for (Int_t j=0; j<nVar; ++j) { ffrac += (*_M)(j,i) * (j==0?1.:TMath::Power(dm,(double)j)); }
178  if (ffrac>=0) sumposfrac+=ffrac;
179  (*_frac)(i) = ffrac;
180  }
181 
182  // various mode settings
183  int imin = idxmin(m);
184  int imax = idxmax(m);
185  double mfrac = (m-(*_mref)[imin])/((*_mref)[imax]-(*_mref)[imin]);
186  switch (_setting) {
187  case NonLinear:
188  // default already set above
189  break;
190  case Linear:
191  for (Int_t i=0; i<nVar; ++i)
192  (*_frac)(i) = 0.;
193  if (imax>imin) { // m in between mmin and mmax
194  (*_frac)(imin) = (1.-mfrac);
195  (*_frac)(imax) = (mfrac);
196  } else if (imax==imin) { // m outside mmin and mmax
197  (*_frac)(imin) = (1.);
198  }
199  break;
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  if ((*_frac)(i)<0) (*_frac)(i)=(0.);
213  (*_frac)(i) = (*_frac)(i)/sumposfrac;
214  }
215  break;
216  }
217 }
218 
219 ////////////////////////////////////////////////////////////////////////////////
220 
221 int Roo1DMomentMorphFunction::idxmin(const double& mval) const
222 {
223  int imin(0);
224  Int_t nVar = _varList.getSize();
225  double mmin=-DBL_MAX;
226  for (Int_t i=0; i<nVar; ++i)
227  if ( (*_mref)[i]>mmin && (*_mref)[i]<=mval ) { mmin=(*_mref)[i]; imin=i; }
228  return imin;
229 }
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 
233 int Roo1DMomentMorphFunction::idxmax(const double& mval) const
234 {
235  int imax(0);
236  Int_t nVar = _varList.getSize();
237  double mmax=DBL_MAX;
238  for (Int_t i=0; i<nVar; ++i)
239  if ( (*_mref)[i]<mmax && (*_mref)[i]>=mval ) { mmax=(*_mref)[i]; imax=i; }
240  return imax;
241 }
TMatrixD * _M
do not persist
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
TIterator * createIterator(Bool_t dir=kIterForward) const
#define coutE(a)
Definition: RooMsgService.h:34
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:75
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Iterator abstract base class.
Definition: TIterator.h:30
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:628
int idxmax(const double &m) const
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:22
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:22
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:336
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:58
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
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...