Logo ROOT   6.10/09
Reference Guide
Roo2DMomentMorphFunction.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * author: Max Baak (mbaak@cern.ch) *
4  *****************************************************************************/
5 
6 /** \class Roo2DMomentMorphFunction
7  \ingroup Roofit
8 
9 2-dimensional morph function between a list of function-numbers as a function of two input parameters (m1 and m2).
10 The matrix mrefpoints assigns a function value to each m1,m2 coordinate.
11 Morphing can be set to be linear, non-linear or a mixture of the two.
12 
13 Usage example:
14 ~~~ {.cpp}
15  TMatrixD foo(9,3);
16 
17  foo(0,0) = 0; // coordinate of variable m1
18  foo(0,1) = 0; // coordinate of variable m2
19  foo(0,2) = 0**2; // function value at (m1,m2) = 0,0
20  foo(1,0) = 100;
21  foo(1,1) = 0;
22  foo(1,2) = 1**2;
23  foo(2,0) = 200;
24  foo(2,1) = 0;
25  foo(2,2) = 2**2;
26 
27  foo(3,0) = 0;
28  foo(3,1) = 150;
29  foo(3,2) = 4**2;
30  foo(4,0) = 100;
31  foo(4,1) = 150;
32  foo(4,2) = 5**2;
33  foo(5,0) = 200;
34  foo(5,1) = 150;
35  foo(5,2) = 8**2;
36 
37  foo(6,0) = 0;
38  foo(6,1) = 300;
39  foo(6,2) = 8**2;
40  foo(7,0) = 100;
41  foo(7,1) = 300;
42  foo(7,2) = 8.5**2;
43  foo(8,0) = 200;
44  foo(8,1) = 300;
45  foo(8,2) = 9**2;
46 
47  // need to provide at least 4 coordinates to Roo2DMomentMorphFunction for 2d extrapolation
48 
49  foo.Print();
50 
51  RooRealVar m1("m1","m1",50,0,200);
52  RooRealVar m2("m2","m2",50,0,300);
53  Roo2DMomentMorphFunction bar("bar","bar", m1, m2, foo );
54 
55  bar.Print();
56 ~~~
57 */
58 
59 
60 #include "Riostream.h"
61 
63 #include "RooAbsReal.h"
64 #include "RooAbsCategory.h"
65 #include <math.h>
66 #include "TMath.h"
67 #include "TTree.h"
68 
69 
70 using namespace std;
71 
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 /// Constructor.
76 /// Cross-check that we have enough reference points.
77 /// \param[in] name
78 /// \param[in] title
79 /// \param[in] _m1
80 /// \param[in] _m2
81 /// \param[in] mrefpoints
82 /// \param[in] setting
83 /// \param[in] verbose
84 
85 Roo2DMomentMorphFunction::Roo2DMomentMorphFunction(const char *name, const char *title,
86  RooAbsReal& _m1, RooAbsReal& _m2,
87  const TMatrixD& mrefpoints, const Setting& setting, const Bool_t& verbose ) :
88  RooAbsReal(name,title),
89  m1("m1","m1",this,_m1),
90  m2("m2","m2",this,_m2),
91  _setting(setting),
92  _verbose(verbose),
93  _npoints( mrefpoints.GetNrows() ),
94  _mref(mrefpoints)
95 {
96  if ( mrefpoints.GetNrows()<4 ) {
97  cerr << "Roo2DMomentMorphFunction::constructor(" << GetName() << ") ERROR: less than four reference points provided." << endl ;
98  assert(0);
99  }
100  // cross-check that we have enough reference points
101  if ( mrefpoints.GetNcols()!=3 ) {
102  cerr << "RooPolyMorph2D::constructor(" << GetName() << ") ERROR: no reference values provided." << endl ;
103  assert(0);
104  }
105 
106  // recast matrix into more useful form
107  _frac.ResizeTo( _npoints );
108 
109  initialize();
110 }
111 
112 ////////////////////////////////////////////////////////////////////////////////
113 /// Constructor.
114 /// Cross-check that we have enough reference points.
115 /// \param[in] name
116 /// \param[in] title
117 /// \param[in] _m1
118 /// \param[in] _m2
119 /// \param[in] nrows
120 /// \param[in] dm1arr
121 /// \param[in] dm2arr
122 /// \param[in] dvalarr
123 /// \param[in] setting
124 /// \param[in] verbose
125 
127  RooAbsReal& _m1, RooAbsReal& _m2,
128  const Int_t& nrows, const Double_t* dm1arr, const Double_t* dm2arr, const Double_t* dvalarr,
129  const Setting& setting, const Bool_t& verbose ) :
130  RooAbsReal(name,title),
131  m1("m1","m1",this,_m1),
132  m2("m2","m2",this,_m2),
133  _setting(setting),
134  _verbose( verbose ),
135  _npoints( nrows )
136 {
137  if ( nrows<4 ) {
138  cerr << "Roo2DMomentMorphFunction::constructor(" << GetName() << ") ERROR: less than four reference points provided." << endl ;
139  assert(0);
140  }
141 
142  // recast tree into more useful form
143  _mref.ResizeTo( _npoints,3 );
145 
146  for (int i=0; i<_npoints; ++i) {
147  _mref(i,0) = dm1arr[i] ;
148  _mref(i,1) = dm2arr[i] ;
149  _mref(i,2) = dvalarr[i] ;
150  }
151 
152  initialize();
153 }
154 
155 ////////////////////////////////////////////////////////////////////////////////
156 /// Copy constructor.
157 
159  RooAbsReal(other,name),
160  m1("m1",this,other.m1),
161  m2("m2",this,other.m2),
162  _setting(other._setting),
163  _verbose(other._verbose),
164  _npoints(other._npoints),
165  _mref(other._mref),
166  _frac(other._frac)
167 {
168  initialize();
169 }
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 /// Destructor.
173 
175 {
176 }
177 
178 ////////////////////////////////////////////////////////////////////////////////
179 
180 void
182 {
183  Double_t xmin(1e300), xmax(-1e300), ymin(1e300), ymax(-1e300);
184 
185  // transformation matrix for non-linear extrapolation, needed in evaluate()
186  for (Int_t k=0; k<_npoints; ++k) {
187  if (_mref(k,0)<xmin) { xmin=_mref(k,0); _ixmin=k; }
188  if (_mref(k,0)>xmax) { xmax=_mref(k,0); _ixmax=k; }
189  if (_mref(k,1)<ymin) { ymin=_mref(k,1); _iymin=k; }
190  if (_mref(k,1)>ymax) { ymax=_mref(k,1); _iymax=k; }
191  }
192 
193  // resize
194  _MSqr.ResizeTo(4,4);
195  _squareVec.ResizeTo(4,2);
196 }
197 
198 ////////////////////////////////////////////////////////////////////////////////
199 
200 Double_t
202 {
203  if (_verbose) { cout << "Now in evaluate." << endl; }
204  if (_verbose) { cout << "x = " << m1 << " ; y = " << m2 << endl; }
205 
206  calculateFractions(_verbose); // verbose turned off
207 
208  Double_t ret = 0.0;
209  for (Int_t i=0; i<4; ++i) { ret += ( _mref(_squareIdx[i],2) * _frac[_squareIdx[i]] ) ; }
210 
211  if (_verbose) { cout << "End of evaluate : " << ret << endl; }
212 
213  //return (ret>0 ? ret : 0) ;
214  return ret;
215 }
216 
217 ////////////////////////////////////////////////////////////////////////////////
218 
219 void
221 {
222  double sumfrac(0.);
223 
225  // reset all fractions
226  for (Int_t i=0; i<_npoints; ++i) { _frac[i]=0.0; }
227 
228  // this sets _squareVec and _squareIdx quantities and MSqr
229  (void) findSquare(m1,m2);
230 
231  std::vector<double> deltavec(4,1.0);
232  deltavec[1] = m1-_squareVec(0,0) ;
233  deltavec[2] = m2-_squareVec(0,1) ;
234  deltavec[3] = deltavec[1]*deltavec[2] ;
235 
236  for (Int_t i=0; i<4; ++i) {
237  double ffrac=0.;
238  for (Int_t j=0; j<4; ++j) { ffrac += _MSqr(j,i) * deltavec[j]; }
239 
240  // set fractions
241  _frac[_squareIdx[i]] = ffrac;
242  if (ffrac>0) sumfrac += ffrac;
243 
244  if (verbose) {
245  cout << _squareIdx[i] << " " << ffrac << " " << _squareVec(i,0) << " " << _squareVec(i,1) << endl;
246  }
247  }
248  }
249 
250  if (_setting == LinearPosFractions) {
251  for (Int_t i=0; i<4; ++i) {
252  if (_frac[_squareIdx[i]]<0) { _frac[_squareIdx[i]] = 0; }
253  _frac[_squareIdx[i]] *= (1.0/sumfrac) ;
254  }
255  }
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 
260 Bool_t
261 Roo2DMomentMorphFunction::findSquare(const double& x, const double& y) const
262 {
263  bool squareFound(false);
264 
265  std::vector< std::pair<int,double> > idvec;
266 
267  Double_t xspacing = (_mref(_ixmax,0)-_mref(_ixmin,0)) / TMath::Sqrt(_npoints) ;
268  Double_t yspacing = (_mref(_iymax,1)-_mref(_iymin,1)) / TMath::Sqrt(_npoints) ;
269 
270  Double_t dx(0), dy(0), delta(0);
271  for (Int_t k=0; k<_npoints; ++k) {
272  dx = (x-_mref(k,0))/xspacing ;
273  dy = (y-_mref(k,1))/yspacing ;
274  delta = TMath::Sqrt(dx*dx+dy*dy) ;
275  idvec.push_back( std::pair<int,double>(k,delta) );
276  }
277 
278  // order points closest to (x,y)
279  sort(idvec.begin(),idvec.end(),SorterL2H());
280  std::vector< std::pair<int,double> >::iterator itr = idvec.begin();
281  std::vector<int> cidx;
282  for(; itr!=idvec.end(); ++itr) {
283  cidx.push_back(itr->first);
284  }
285 
286  // 1. point falls outside available ref points: pick three square points.
287  // fall-back option
288  _squareVec(0,0) = _mref(cidx[0],0);
289  _squareVec(0,1) = _mref(cidx[0],1);
290  _squareIdx[0] = cidx[0];
291  _squareVec(1,0) = _mref(cidx[1],0);
292  _squareVec(1,1) = _mref(cidx[1],1);
293  _squareIdx[1] = cidx[1];
294  _squareVec(2,0) = _mref(cidx[2],0);
295  _squareVec(2,1) = _mref(cidx[2],1);
296  _squareIdx[2] = cidx[2];
297  _squareVec(3,0) = _mref(cidx[3],0);
298  _squareVec(3,1) = _mref(cidx[3],1);
299  _squareIdx[3] = cidx[3];
300 
301  // 2. try to find square enclosing square
302  if ( x>_mref(_ixmin,0) &&
303  x<_mref(_ixmax,0) &&
304  y>_mref(_iymin,1) &&
305  y<_mref(_iymax,1) )
306  {
307  for (unsigned int i=0; i<cidx.size() && !squareFound; ++i)
308  for (unsigned int j=i+1; j<cidx.size() && !squareFound; ++j)
309  for (unsigned int k=j+1; k<cidx.size() && !squareFound; ++k)
310  for (unsigned int l=k+1; l<cidx.size() && !squareFound; ++l) {
311  if ( isAcceptableSquare(_mref(cidx[i],0),_mref(cidx[i],1),_mref(cidx[j],0),_mref(cidx[j],1),_mref(cidx[k],0),_mref(cidx[k],1),_mref(cidx[l],0),_mref(cidx[l],1)) ) {
312  if ( pointInSquare(x,y,_mref(cidx[i],0),_mref(cidx[i],1),_mref(cidx[j],0),_mref(cidx[j],1),_mref(cidx[k],0),_mref(cidx[k],1),_mref(cidx[l],0),_mref(cidx[l],1)) ) {
313  _squareVec(0,0) = _mref(cidx[i],0);
314  _squareVec(0,1) = _mref(cidx[i],1);
315  _squareIdx[0] = cidx[i];
316  _squareVec(1,0) = _mref(cidx[j],0);
317  _squareVec(1,1) = _mref(cidx[j],1);
318  _squareIdx[1] = cidx[j];
319  _squareVec(2,0) = _mref(cidx[k],0);
320  _squareVec(2,1) = _mref(cidx[k],1);
321  _squareIdx[2] = cidx[k];
322  _squareVec(3,0) = _mref(cidx[l],0);
323  _squareVec(3,1) = _mref(cidx[l],1);
324  _squareIdx[3] = cidx[l];
325  squareFound=true;
326  }
327  }
328  }
329  }
330 
331  // construct transformation matrix for linear extrapolation
332  TMatrixD M(4,4);
333  for (Int_t k=0; k<4; ++k) {
334  dx = _squareVec(k,0) - _squareVec(0,0) ;
335  dy = _squareVec(k,1) - _squareVec(0,1) ;
336  M(k,0) = 1.0 ;
337  M(k,1) = dx ;
338  M(k,2) = dy ;
339  M(k,3) = dx*dy ;
340  }
341 
342  _MSqr = M.Invert();
343 
344  return squareFound;
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Returns true if p1 and p2 are on same side of line b-a.
349 /// \param[in] p1x
350 /// \param[in] p1y
351 /// \param[in] p2x
352 /// \param[in] p2y
353 /// \param[in] ax
354 /// \param[in] ay
355 /// \param[in] bx
356 /// \param[in] by
357 
358 Bool_t Roo2DMomentMorphFunction::onSameSide(const double& p1x, const double& p1y, const double& p2x, const double& p2y, const double& ax, const double& ay, const double& bx, const double& by) const
359 {
360  Double_t cp1 = myCrossProduct(bx-ax, by-ay, p1x-ax, p1y-ay);
361  Double_t cp2 = myCrossProduct(bx-ax, by-ay, p2x-ax, p2y-ay);
362  if (cp1*cp2 >= 0) return true;
363  else return false;
364 }
365 
366 ////////////////////////////////////////////////////////////////////////////////
367 
368 Bool_t
369 Roo2DMomentMorphFunction::pointInSquare(const double& px, const double& py, const double& ax, const double& ay, const double& bx, const double& by, const double& cx, const double& cy, const double& dx, const double& dy) const
370 {
371  bool insquare(false);
372 
373  int ntri(0);
374  if (ntri<2) ntri += static_cast<int>( pointInTriangle(px,py,ax,ay,bx,by,cx,cy) );
375  if (ntri<2) ntri += static_cast<int>( pointInTriangle(px,py,ax,ay,bx,by,dx,dy) );
376  if (ntri<2) ntri += static_cast<int>( pointInTriangle(px,py,ax,ay,cx,cy,dx,dy) );
377  if (ntri<2) ntri += static_cast<int>( pointInTriangle(px,py,bx,by,cx,cy,dx,dy) );
378 
379  if (ntri>=2) insquare=true;
380  else insquare=false;
381 
382  return insquare;
383 }
384 
385 ////////////////////////////////////////////////////////////////////////////////
386 
387 Bool_t
388 Roo2DMomentMorphFunction::pointInTriangle(const double& px, const double& py, const double& ax, const double& ay, const double& bx, const double& by, const double& cx, const double& cy) const
389 {
390  if (onSameSide(px,py,ax,ay,bx,by,cx,cy) && onSameSide(px,py,bx,by,ax,ay,cx,cy) && onSameSide(px,py,cx,cy,ax,ay,bx,by)) return true;
391  else return false;
392 }
393 
394 ////////////////////////////////////////////////////////////////////////////////
395 
396 Double_t
397 Roo2DMomentMorphFunction::myCrossProduct(const double& ax, const double& ay, const double& bx, const double& by) const
398 {
399  return ( ax*by - bx*ay );
400 }
401 
402 ////////////////////////////////////////////////////////////////////////////////
403 /// Reject kinked shapes.
404 /// \param[in] ax
405 /// \param[in] ay
406 /// \param[in] bx
407 /// \param[in] by
408 /// \param[in] cx
409 /// \param[in] cy
410 /// \param[in] dx
411 /// \param[in] dy
412 
413 Bool_t
414 Roo2DMomentMorphFunction::isAcceptableSquare(const double& ax, const double& ay, const double& bx, const double& by, const double& cx, const double& cy, const double& dx, const double& dy) const
415 {
416  if ( pointInTriangle(dx,dy,ax,ax,bx,by,cx,cy) ||
417  pointInTriangle(cx,cy,ax,ay,bx,by,dx,dy) ||
418  pointInTriangle(bx,by,ax,ay,cx,cy,dx,dy) ||
419  pointInTriangle(ax,ay,bx,by,cx,cy,dx,dy) ) return false;
420  else return true;
421 }
422 
423 
424 void
426 {
427  for( Int_t i=0; i<_npoints; i++ ){
428  cout << this << " " << i << " " << _mref(i,0) << " " << _mref(i,1) << " " << _mref(i,2) << endl;
429  }
430 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:145
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:292
float xmin
Definition: THbookFile.cxx:93
virtual ~Roo2DMomentMorphFunction()
Destructor.
Bool_t pointInTriangle(const double &px, const double &py, const double &ax, const double &ay, const double &bx, const double &by, const double &cx, const double &cy) const
float ymin
Definition: THbookFile.cxx:93
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
STL namespace.
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Definition: TMatrixT.cxx:1210
Float_t delta
Double_t x[n]
Definition: legend1.C:17
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
Definition: TMatrixT.cxx:1396
Double_t myCrossProduct(const double &ax, const double &ay, const double &bx, const double &by) const
2-dimensional morph function between a list of function-numbers as a function of two input parameters...
float ymax
Definition: THbookFile.cxx:93
Bool_t findSquare(const double &x, const double &y) const
Bool_t isAcceptableSquare(const double &ax, const double &ay, const double &bx, const double &by, const double &cx, const double &cy, const double &dx, const double &dy) const
Reject kinked shapes.
bool verbose
TLine * l
Definition: textangle.C:4
float xmax
Definition: THbookFile.cxx:93
void calculateFractions(Bool_t verbose=kTRUE) const
#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
Bool_t onSameSide(const double &p1x, const double &p1y, const double &p2x, const double &p2y, const double &ax, const double &ay, const double &bx, const double &by) const
Returns true if p1 and p2 are on same side of line b-a.
Double_t y[n]
Definition: legend1.C:17
typedef void((*Func_t)())
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
void initialize(typename Architecture_t::Matrix_t &A, EInitialization m)
Definition: Functions.h:257
Bool_t pointInSquare(const double &px, const double &py, const double &ax, const double &ay, const double &bx, const double &by, const double &cx, const double &cy, const double &dx, const double &dy) const