Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SequentialProposal.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Authors: Giovanni Petrucciani 4/21/2011
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11/** \class RooStats::SequentialProposal
12 \ingroup Roostats
13
14Class implementing a proposal function that samples the parameter space
15by moving only in one coordinate (chosen randomly) at each step
16*/
17
19#include <RooArgSet.h>
20#include <iostream>
21#include <memory>
22#include <RooRandom.h>
24
25
27
28namespace RooStats {
29
30////////////////////////////////////////////////////////////////////////////////
31
32SequentialProposal::SequentialProposal(double divisor) : fDivisor(1. / divisor) {}
33
34////////////////////////////////////////////////////////////////////////////////
35/// Populate xPrime with a new proposed point
36
38{
39 RooStats::SetParameters(&x, &xPrime);
40 int n = xPrime.size();
41 int j = int( floor(RooRandom::uniform()*n) );
42 int i = 0;
43 for (auto *var : static_range_cast<RooRealVar *>(xPrime)) {
44 if (i == j) {
45 double val = var->getVal();
46 double max = var->getMax();
47 double min = var->getMin();
48 double len = max - min;
49 val += RooRandom::gaussian() * len * fDivisor;
50 while (val > max)
51 val -= len;
52 while (val < min)
53 val += len;
54 var->setVal(val);
55 // std::cout << "Proposing a step along " << var->GetName() << std::endl;
56 }
57 ++i;
58 }
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Return the probability of proposing the point x1 given the starting
63/// point x2
64
66 return true;
67}
68
70 RooArgSet& )
71{
72 return 1.0; // should not be needed
73}
74
75}
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Storage_t::size_type size() const
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:78
static double gaussian(TRandom *generator=randomGenerator())
Return a Gaussian random variable with mean 0 and variance 1.
Class implementing a proposal function that samples the parameter space by moving only in one coordin...
double GetProposalDensity(RooArgSet &x1, RooArgSet &x2) override
Return the probability of proposing the point x1 given the starting point x2.
bool IsSymmetric(RooArgSet &x1, RooArgSet &x2) override
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is,...
void Propose(RooArgSet &xPrime, RooArgSet &x) override
Populate xPrime with a new proposed point.
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Namespace for the RooStats classes.
Definition Asimov.h:19
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)