Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
26
27namespace RooStats {
28
29////////////////////////////////////////////////////////////////////////////////
30
32
33////////////////////////////////////////////////////////////////////////////////
34/// Populate xPrime with a new proposed point
35
37{
39 int n = xPrime.size();
40 int j = int( floor(RooRandom::uniform()*n) );
41 int i = 0;
42 for (auto *var : static_range_cast<RooRealVar *>(xPrime)) {
43 if (i == j) {
44 double val = var->getVal();
45 double max = var->getMax();
46 double min = var->getMin();
47 double len = max - min;
48 val += RooRandom::gaussian() * len * fDivisor;
49 while (val > max)
50 val -= len;
51 while (val < min)
52 val += len;
53 var->setVal(val);
54 // std::cout << "Proposing a step along " << var->GetName() << std::endl;
55 }
56 ++i;
57 }
58}
59
60////////////////////////////////////////////////////////////////////////////////
61/// Return the probability of proposing the point x1 given the starting
62/// point x2
63
65 return true;
66}
67
69 RooArgSet& )
70{
71 return 1.0; // should not be needed
72}
73
74}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:77
static double gaussian(TRandom *generator=randomGenerator())
Return a Gaussian random variable with mean 0 and variance 1.
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 CodegenImpl.h:58
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)