44 _parItr = _parList.createIterator();
45 _obsItr = _obsList.createIterator();
54 :
RooAbsReal(name, title), _cacheMgr(this, 10,
kTRUE,
kTRUE), _parList(
"parList",
"List of morph parameters", this),
55 _obsList(
"obsList",
"List of observables", this), _referenceGrid(referenceGrid),
56 _pdfList(
"pdfList",
"List of pdfs", this), _setting(setting), _useHorizMorph(true)
84 for (
int i = 0; i < mrefpoints.
GetNrows(); ++i) {
85 for (
int j = 0; j < grid.numBoundaries(); ++j) {
86 if (mrefpoints[i] == grid.array()[j]) {
121 for (
int i = 0; (mref =
dynamic_cast<RooAbsReal *
>(mrefItr->Next())); ++i) {
124 <<
" is not of type RooAbsReal" << endl;
125 throw string(
"RooMomentMorphFuncND::ctor() ERROR mref is not of type RooAbsReal");
127 if (!dynamic_cast<RooConstVar *>(mref)) {
129 <<
" is not a constant, taking a snapshot of its value" << endl;
131 mrefpoints[i] = mref->
getVal();
135 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
138 for (
int i = 0; i < mrefpoints.GetNrows(); ++i) {
139 for (
int j = 0; j < grid.numBoundaries(); ++j) {
140 if (mrefpoints[i] == grid.array()[j]) {
200 if (!dynamic_cast<RooAbsReal *>(par)) {
202 <<
" is not of type RooAbsReal" << endl;
203 throw string(
"RooMomentMorphFuncND::initializeParameters() ERROR parameter is not of type RooAbsReal");
218 if (!dynamic_cast<RooAbsReal *>(var)) {
220 <<
" is not of type RooAbsReal" << endl;
221 throw string(
"RooMomentMorphFuncND::initializeObservables() ERROR variable is not of type RooAbsReal");
232 template <
typename T>
235 vector<Digits<T>> vd;
237 for (
typename vector<vector<T>>::const_iterator it = in.begin(); it != in.end(); ++it) {
244 for (
typename vector<
Digits<T>>::const_iterator it = vd.begin(); it != vd.end(); ++it) {
245 result.push_back(*(it->me));
247 out.push_back(result);
249 for (
typename vector<
Digits<T>>::iterator it = vd.begin();;) {
251 if (it->me == it->end) {
252 if (it + 1 == vd.end()) {
293 <<
": " << nPar <<
" !=" << nDim << endl;
299 <<
": " << nPdf <<
" !=" << nRef << endl;
309 vector<vector<double>> dm(nPdf);
310 for (
int k = 0; k < nPdf; ++k) {
312 for (
int idim = 0; idim < nPar; idim++) {
314 dm2.push_back(delta);
319 vector<vector<int>> powers;
320 for (
int idim = 0; idim < nPar; idim++) {
325 powers.push_back(xtmp);
328 vector<vector<int>>
output;
330 int nCombs = output.size();
332 for (
int k = 0; k < nPdf; ++k) {
334 for (
int i = 0; i < nCombs; i++) {
336 for (
int ix = 0; ix < nPar; ix++) {
338 tmpDm *=
TMath::Power(delta, static_cast<double>(output[i][ix]));
356 : _grid(other._grid),
_pdfList(other.
_pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
368 vector<int> thisBoundaries;
369 vector<double> thisBoundaryCoordinates;
370 thisBoundaries.push_back(bin_x);
371 thisBoundaryCoordinates.push_back(
_grid[0]->array()[bin_x]);
374 _nref.push_back(thisBoundaryCoordinates);
380 vector<int> thisBoundaries;
381 vector<double> thisBoundaryCoordinates;
382 thisBoundaries.push_back(bin_x);
383 thisBoundaryCoordinates.push_back(
_grid[0]->array()[bin_x]);
384 thisBoundaries.push_back(bin_y);
385 thisBoundaryCoordinates.push_back(
_grid[1]->array()[bin_y]);
388 _nref.push_back(thisBoundaryCoordinates);
394 vector<int> thisBoundaries;
395 vector<double> thisBoundaryCoordinates;
396 thisBoundaries.push_back(bin_x);
397 thisBoundaryCoordinates.push_back(
_grid[0]->array()[bin_x]);
398 thisBoundaries.push_back(bin_y);
399 thisBoundaryCoordinates.push_back(
_grid[1]->array()[bin_y]);
400 thisBoundaries.push_back(bin_z);
401 thisBoundaryCoordinates.push_back(
_grid[2]->array()[bin_z]);
404 _nref.push_back(thisBoundaryCoordinates);
410 vector<double> thisBoundaryCoordinates;
411 int nBins = bins.size();
412 for (
int i = 0; i < nBins; i++) {
413 thisBoundaryCoordinates.push_back(
_grid[i]->array()[bins[i]]);
417 _nref.push_back(thisBoundaryCoordinates);
434 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
435 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
436 vector<RooAbsReal *> myrms(nObs, null);
437 vector<RooAbsReal *> mypos(nObs, null);
438 vector<RooAbsReal *> slope(nPdf * nObs, null);
439 vector<RooAbsReal *> offsets(nPdf * nObs, null);
440 vector<RooAbsReal *> transVar(nPdf * nObs, null);
441 vector<RooAbsReal *> transPdf(nPdf, null);
451 for (
int i = 0; i < 3 * nPdf; ++i) {
452 string fracName =
Form(
"frac_%d", i);
453 double initval = 0.0;
459 else if (i < 2 * nPdf)
467 string sumfuncName =
Form(
"%s_sumfunc",
GetName());
472 for (
int i = 0; i < nPdf; ++i) {
473 for (
int j = 0; j < nObs; ++j) {
481 sigmarv[
sij(i, j)] = mom;
482 meanrv[
sij(i, j)] = mom->
mean();
484 ownedComps.
add(*sigmarv[
sij(i, j)]);
489 for (
int j = 0; j < nObs; ++j) {
492 for (
int i = 0; i < nPdf; ++i) {
493 meanList.
add(*meanrv[
sij(i, j)]);
494 rmsList.
add(*sigmarv[
sij(i, j)]);
496 string myrmsName =
Form(
"%s_rms_%d",
GetName(), j);
497 string myposName =
Form(
"%s_pos_%d",
GetName(), j);
498 mypos[j] =
new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
499 myrms[j] =
new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
508 for (
int i = 0; i < nPdf; ++i) {
513 string pdfName =
Form(
"pdf_%d", i);
516 for (
int j = 0; j < nObs; ++j) {
518 string slopeName =
Form(
"%s_slope_%d_%d",
GetName(), i, j);
519 string offsetName =
Form(
"%s_offset_%d_%d",
GetName(), i, j);
529 string transVarName =
Form(
"%s_transVar_%d_%d",
GetName(), i, j);
530 transVar[
sij(i, j)] =
new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[
sij(i, j)],
531 *offsets[
sij(i, j)]);
537 ownedComps.
add(*transVar[
sij(i, j)]);
538 cust.replaceArg(*var, *transVar[
sij(i, j)]);
541 transPdfList.
add(*transPdf[i]);
542 ownedComps.
add(*transPdf[i]);
546 theSumFunc =
new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), transPdfList, coefList);
557 string trackerName =
Form(
"%s_frac_tracker",
GetName());
561 cache =
new CacheElem(*theSumFunc, *tracker, fracl);
627 template <
typename Iterator>
630 if ((first == last) || (first == k) || (last == k)) {
633 Iterator itr1 = first;
634 Iterator itr2 = last;
643 while (first != itr1) {
644 if (*--itr1 < *itr2) {
646 while (!(*itr1 < *j)) ++j;
651 rotate(itr1, j, last);
656 rotate(k, itr2, last);
660 rotate(first, k, last);
667 int nPdf =
self._pdfList.getSize();
668 int nPar =
self._parList.getSize();
676 for (
int idim = 0; idim < nPar; idim++) {
678 dm2.push_back(delta);
681 vector<vector<int>> powers;
682 for (
int idim = 0; idim < nPar; idim++) {
684 for (
int ix = 0; ix <
self._referenceGrid._nnuis[idim]; ix++) {
687 powers.push_back(xtmp);
690 vector<vector<int>>
output;
692 int nCombs = output.size();
694 vector<double> deltavec(nPdf, 1.0);
697 for (
int i = 0; i < nCombs; i++) {
699 for (
int ix = 0; ix < nPar; ix++) {
701 tmpDm *=
TMath::Power(delta, static_cast<double>(output[i][ix]));
703 deltavec[nperm] = tmpDm;
707 double sumposfrac = 0.0;
708 for (
int i = 0; i < nPdf; ++i) {
711 for (
int j = 0; j < nPdf; ++j) {
712 ffrac += (*
self._M)(j, i) * deltavec[j] * fracNonLinear;
725 ((
RooRealVar *)frac(nPdf + i))->setVal(ffrac);
726 ((
RooRealVar *)frac(2 * nPdf + i))->setVal(ffrac);
729 cout <<
"NonLinear fraction " << ffrac << endl;
731 frac(nPdf + i)->Print();
732 frac(2 * nPdf + i)->Print();
737 for (
int i = 0; i < nPdf; ++i) {
747 self._parItr->Reset();
751 for (
int i = 0; i < nPdf; ++i) {
754 ((
RooRealVar *)frac(nPdf + i))->setVal(initval);
755 ((
RooRealVar *)frac(2 * nPdf + i))->setVal(initval);
760 for (
int j = 0; j < nPar; j++) {
762 mtmp.push_back(m->
getVal());
765 self.findShape(mtmp);
768 vector<double> deltavec(depth, 1.0);
773 for (
int ix = 0; ix < nPar; ix++) {
777 for (
int iperm = 1; iperm <= nPar; ++iperm) {
779 double dtmp = mtmp[xtmp[0]] -
self._squareVec[0][xtmp[0]];
780 for (
int itmp = 1; itmp < iperm; ++itmp) {
781 dtmp *= mtmp[xtmp[itmp]] -
self._squareVec[0][xtmp[itmp]];
783 deltavec[nperm + 1] = dtmp;
788 Double_t origFrac1(0.), origFrac2(0.);
789 for (
int i = 0; i < depth; ++i) {
791 for (
int j = 0; j < depth; ++j) {
792 ffrac += (*
self._MSqr)(j, i) * deltavec[j] * fracLinear;
808 cout <<
"Linear fraction " << ffrac << endl;
824 bool isEnclosed =
true;
825 for (
int i = 0; i < nPar; i++) {
836 vector<vector<double>> boundaries(nPar);
837 for (
int idim = 0; idim < nPar; idim++) {
841 boundaries[idim].push_back(lo);
842 boundaries[idim].push_back(hi);
845 vector<vector<double>>
output;
849 for (
int isq = 0; isq < depth; isq++) {
850 for (
int iref = 0; iref < nRef; iref++) {
873 for (
int ix = 0; ix < nPar; ix++) {
877 for (
int k = 0; k < depth; ++k) {
883 for (
int iperm = 1; iperm <= nPar; ++iperm) {
885 double dtmp =
_squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
886 for (
int itmp = 1; itmp < iperm; ++itmp) {
887 dtmp *=
_squareVec[k][xtmp[itmp]] - squareBase[xtmp[itmp]];
889 M(k, nperm + 1) = dtmp;
909 cout <<
"Currently BinIntegrator only knows how to deal with 1-d " << endl;
virtual const char * GetName() const
Returns name of object.
void findShape(const vector< double > &x) const
TIterator * _obsItr
do not persist
TIterator * createIterator(Bool_t dir=kIterForward) const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
vector< T >::const_iterator begin
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
RooAbsReal * sumFunc(const RooArgSet *nset)
void initializeObservables(const RooArgList &obsList)
Double_t getVal(const RooArgSet *set=0) const
CacheElem * getCache(const RooArgSet *nset) const
Bool_t setBinIntegrator(RooArgSet &allVars)
void initializeParameters(const RooArgList &parList)
void addPdf(const RooAbsReal &func, int bin_x)
Bool_t hasChanged(Bool_t clearState)
Returns true if state has changes since last call with clearState=kTRUE If clearState is true...
vector< vector< double > > _squareVec
int sij(const int &i, const int &j) const
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
Iterator abstract base class.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
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...
virtual ~RooMomentMorphFuncND()
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
RooAbsMoment * sigma(RooRealVar &obs)
vector< vector< double > > _nref
Class RooBinning is an implements RooAbsBinning in terms of an array of boundary values, posing no constraints on the choice of binning, thus allowing variable bin sizes.
virtual RooArgList containedArgs(Action)
virtual Double_t getVal(const RooArgSet *set=0) const
RooObjCacheManager _cacheMgr
void calculateFractions(const RooMomentMorphFuncND &self, Bool_t verbose=kTRUE) const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
const RooArgSet * nset() const
RooRealVar represents a fundamental (non-derived) real valued object.
Element * GetMatrixArray()
TMatrixT< Double_t > TMatrixD
void addServerList(RooAbsCollection &serverList, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE)
Register a list of RooAbsArg as servers to us by calls addServer() for each arg in the list...
static bool next_combination(const Iterator first, Iterator k, const Iterator last)
friend class RooRealSumFunc
RooAbsArg * at(Int_t idx) const
RooAbsArg * first() const
char * Form(const char *fmt,...)
Double_t evaluate() const
static void cartesian_product(vector< vector< T >> &out, vector< vector< T >> &in)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
RooChangeTracker * _tracker
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects...
static std::shared_ptr< std::function< double(double)> > Linear
typedef void((*Func_t)())
vector< RooAbsBinning * > _grid
virtual TObject * Next()=0
float type_of_call hi(const int &, const int &)
map< vector< int >, int > _pdfMap
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets...
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
void setLocalNoDirtyInhibit(Bool_t flag) const
void addBinning(const RooAbsBinning &binning)
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...