// @(#)root/roostats:$Id$
// Authors: Kevin Belasco        17/06/2009
// Authors: Kyle Cranmer         17/06/2009
/*************************************************************************
 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//_________________________________________________
/*
BEGIN_HTML
<p>
This class uses the Metropolis-Hastings algorithm to construct a Markov Chain
of data points using Monte Carlo. In the main algorithm, new points in the
parameter space are proposed and then visited based on their relative
likelihoods.  This class can use any implementation of the ProposalFunction,
including non-symmetric proposal functions, to propose parameter points and
still maintain detailed balance when constructing the chain.
</p>

<p>
The "Likelihood" function that is sampled when deciding what steps to take in
the chain has been given a very generic implementation.  The user can create
any RooAbsReal based on the parameters and pass it to a MetropolisHastings
object with the method SetFunction(RooAbsReal&).  Be sure to tell
MetropolisHastings whether your RooAbsReal is on a (+/-) regular or log scale,
so that it knows what logic to use when sampling your RooAbsReal.  For example,
a common use is to sample from a -log(Likelihood) distribution (NLL), for which
the appropriate configuration calls are SetType(MetropolisHastings::kLog);
SetSign(MetropolisHastings::kNegative);
If you're using a traditional likelihood function:
SetType(MetropolisHastings::kRegular);  SetSign(MetropolisHastings::kPositive);
You must set these type and sign flags or MetropolisHastings will not construct
a MarkovChain.
</p>

<p>
Also note that in ConstructChain(), the values of the variables are randomized
uniformly over their intervals before construction of the MarkovChain begins.
</p>
END_HTML
*/
//_________________________________________________

#ifndef RooStats_RooStatsUtils
#include "RooStats/RooStatsUtils.h"
#endif
#ifndef ROOT_Rtypes
#include "Rtypes.h"
#endif
#ifndef ROO_REAL_VAR
#include "RooRealVar.h"
#endif
#ifndef ROO_NLL_VAR
#include "RooNLLVar.h"
#endif
#ifndef ROO_GLOBAL_FUNC
#include "RooGlobalFunc.h"
#endif
#ifndef ROO_DATA_SET
#include "RooDataSet.h"
#endif
#ifndef ROO_ARG_SET
#include "RooArgSet.h"
#endif
#ifndef ROO_ARG_LIST
#include "RooArgList.h"
#endif
#ifndef ROO_MSG_SERVICE
#include "RooMsgService.h"
#endif
#ifndef ROO_RANDOM
#include "RooRandom.h"
#endif
#ifndef ROOT_TH1
#include "TH1.h"
#endif
#ifndef ROOT_TMath
#include "TMath.h"
#endif
#ifndef ROOT_TFile
#include "TFile.h"
#endif
#ifndef ROOSTATS_MetropolisHastings
#include "RooStats/MetropolisHastings.h"
#endif
#ifndef ROOSTATS_MarkovChain
#include "RooStats/MarkovChain.h"
#endif
#ifndef RooStats_MCMCInterval
#include "RooStats/MCMCInterval.h"
#endif

ClassImp(RooStats::MetropolisHastings);

using namespace RooFit;
using namespace RooStats;
using namespace std;

MetropolisHastings::MetropolisHastings()
{
   // default constructor
   fFunction = NULL;
   fPropFunc = NULL;
   fNumIters = 0;
   fNumBurnInSteps = 0;
   fSign = kSignUnset;
   fType = kTypeUnset;
}

MetropolisHastings::MetropolisHastings(RooAbsReal& function, const RooArgSet& paramsOfInterest,
      ProposalFunction& proposalFunction, Int_t numIters)
{
   fFunction = &function;
   SetParameters(paramsOfInterest);
   SetProposalFunction(proposalFunction);
   fNumIters = numIters;
   fNumBurnInSteps = 0;
   fSign = kSignUnset;
   fType = kTypeUnset;
}

MarkovChain* MetropolisHastings::ConstructChain()
{
   if (fParameters.getSize() == 0 || !fPropFunc || !fFunction) {
      coutE(Eval) << "Critical members unintialized: parameters, proposal " <<
                     " function, or (log) likelihood function" << endl;
         return NULL;
   }
   if (fSign == kSignUnset || fType == kTypeUnset) {
      coutE(Eval) << "Please set type and sign of your function using "
         << "MetropolisHastings::SetType() and MetropolisHastings::SetSign()" <<
         endl;
      return NULL;
   }

   if (fChainParams.getSize() == 0) fChainParams.add(fParameters);

   RooArgSet x;
   RooArgSet xPrime;
   x.addClone(fParameters);
   RandomizeCollection(x);
   xPrime.addClone(fParameters);
   RandomizeCollection(xPrime);

   MarkovChain* chain = new MarkovChain();
   // only the POI will be added to the chain
   chain->SetParameters(fChainParams);

   Int_t weight = 0;
   Double_t xL = 0.0, xPrimeL = 0.0, a = 0.0;

   // ibucur: i think the user should have the possiblity to display all the message
   //    levels should they want to; maybe a setPrintLevel would be appropriate
   //    (maybe for the other classes that use this approach as well)?
   RooFit::MsgLevel oldMsgLevel = RooMsgService::instance().globalKillBelow();
   RooMsgService::instance().setGlobalKillBelow(RooFit::PROGRESS);

   // We will need to check if log-likelihood evaluation left an error status.
   // Now using faster eval error logging with CountErrors.
   if (fType == kLog) {
     RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CountErrors);
     //N.B: need to clear the count in case of previous errors !
     // the clear needs also to be done after calling setEvalErrorLoggingMode 
     RooAbsReal::clearEvalErrorLog();
   }

   bool hadEvalError = true;

   Int_t i = 0;
   // get a good starting point for x
   // for fType == kLog, this means that fFunction->getVal() did not cause
   // an eval error
   // for fType == kRegular this means fFunction->getVal() != 0
   //
   // kbelasco: i < 1000 is sort of arbitary, but way higher than the number of
   // steps we should have to take for any reasonable (log) likelihood function
   while (i < 1000 && hadEvalError) {
      RandomizeCollection(x);
      RooStats::SetParameters(&x, &fParameters);
      xL = fFunction->getVal();

      if (fType == kLog) {
         if (RooAbsReal::numEvalErrors() > 0) {
            RooAbsReal::clearEvalErrorLog();
            hadEvalError = true;
         } else
            hadEvalError = false;
      } else if (fType == kRegular) {
         if (xL == 0.0)
            hadEvalError = true;
         else
            hadEvalError = false;
      } else
         // for now the only 2 types are kLog and kRegular (won't get here)
         hadEvalError = false;
   }

   if(hadEvalError) {
      coutE(Eval) << "Problem finding a good starting point in " <<
                     "MetropolisHastings::ConstructChain() " << endl;
   }


   ooccoutP((TObject *)0, Generation) << "Metropolis-Hastings progress: ";

   // do main loop
   for (i = 0; i < fNumIters; i++) {
      // reset error handling flag
      hadEvalError = false;

      // print a dot every 1% of the chain construction
      if (i % (fNumIters / 100) == 0) ooccoutP((TObject*)0, Generation) << ".";

      fPropFunc->Propose(xPrime, x);

      RooStats::SetParameters(&xPrime, &fParameters);
      xPrimeL = fFunction->getVal();

      // check if log-likelihood for xprime had an error status
      if (fFunction->numEvalErrors() > 0 && fType == kLog) {
         xPrimeL = RooNumber::infinity();
         fFunction->clearEvalErrorLog();
         hadEvalError = true;
      }

      // why evaluate the last point again, can't we cache it?
      // kbelasco: commenting out lines below to add/test caching support
      //RooStats::SetParameters(&x, &fParameters);
      //xL = fFunction->getVal();

      if (fType == kLog) {
         if (fSign == kPositive)
            a = xL - xPrimeL;
         else
            a = xPrimeL - xL;
      }
      else
         a = xPrimeL / xL;
      //a = xL / xPrimeL;

      if (!hadEvalError && !fPropFunc->IsSymmetric(xPrime, x)) {
         Double_t xPrimePD = fPropFunc->GetProposalDensity(xPrime, x);
         Double_t xPD      = fPropFunc->GetProposalDensity(x, xPrime);
         if (fType == kRegular)
            a *= xPD / xPrimePD;
         else
            a += TMath::Log(xPrimePD) - TMath::Log(xPD);
      }

      if (!hadEvalError && ShouldTakeStep(a)) {
         // go to the proposed point xPrime

         // add the current point with the current weight
         if (weight != 0.0)
            chain->Add(x, CalcNLL(xL), (Double_t)weight);

         // reset the weight and go to xPrime
         weight = 1;
         RooStats::SetParameters(&xPrime, &x);
         xL = xPrimeL;
      } else {
         // stay at the current point
         weight++;
      }
   }

   // make sure to add the last point
   if (weight != 0.0)
      chain->Add(x, CalcNLL(xL), (Double_t)weight);
   ooccoutP((TObject *)0, Generation) << endl;

   RooMsgService::instance().setGlobalKillBelow(oldMsgLevel);

   Int_t numAccepted = chain->Size();
   coutI(Eval) << "Proposal acceptance rate: " <<
                   numAccepted/(Float_t)fNumIters * 100 << "%" << endl;
   coutI(Eval) << "Number of steps in chain: " << numAccepted << endl;

   //TFile chainDataFile("chainData.root", "recreate");
   //chain->GetDataSet()->Write();
   //chainDataFile.Close();

   return chain;
}

Bool_t MetropolisHastings::ShouldTakeStep(Double_t a)
{
   if ((fType == kLog && a <= 0.0) || (fType == kRegular && a >= 1.0)) {
      // The proposed point has a higher likelihood than the
      // current point, so we should go there
      return kTRUE;
   }
   else {
      // generate numbers on a log distribution to decide
      // whether to go to xPrime or stay at x
      //Double_t rand = fGen.Uniform(1.0);
      Double_t rand = RooRandom::uniform();
      if (fType == kLog) {
         rand = TMath::Log(rand);
         // kbelasco: should this be changed to just (-rand > a) for logical
         // consistency with below test when fType == kRegular?
         if (-1.0 * rand >= a)
            // we chose to go to the new proposed point
            // even though it has a lower likelihood than the current one
            return kTRUE;
      } else {
         // fType must be kRegular
         // kbelasco: ensure that we never visit a point where PDF == 0
         //if (rand <= a)
         if (rand < a)
            // we chose to go to the new proposed point
            // even though it has a lower likelihood than the current one
            return kTRUE;
      }
      return kFALSE;
   }
}

Double_t MetropolisHastings::CalcNLL(Double_t xL)
{
   if (fType == kLog) {
      if (fSign == kNegative)
         return xL;
      else
         return -xL;
   } else {
      if (fSign == kPositive)
         return -1.0 * TMath::Log(xL);
      else
         return -1.0 * TMath::Log(-xL);
   }
}
 MetropolisHastings.cxx:1
 MetropolisHastings.cxx:2
 MetropolisHastings.cxx:3
 MetropolisHastings.cxx:4
 MetropolisHastings.cxx:5
 MetropolisHastings.cxx:6
 MetropolisHastings.cxx:7
 MetropolisHastings.cxx:8
 MetropolisHastings.cxx:9
 MetropolisHastings.cxx:10
 MetropolisHastings.cxx:11
 MetropolisHastings.cxx:12
 MetropolisHastings.cxx:13
 MetropolisHastings.cxx:14
 MetropolisHastings.cxx:15
 MetropolisHastings.cxx:16
 MetropolisHastings.cxx:17
 MetropolisHastings.cxx:18
 MetropolisHastings.cxx:19
 MetropolisHastings.cxx:20
 MetropolisHastings.cxx:21
 MetropolisHastings.cxx:22
 MetropolisHastings.cxx:23
 MetropolisHastings.cxx:24
 MetropolisHastings.cxx:25
 MetropolisHastings.cxx:26
 MetropolisHastings.cxx:27
 MetropolisHastings.cxx:28
 MetropolisHastings.cxx:29
 MetropolisHastings.cxx:30
 MetropolisHastings.cxx:31
 MetropolisHastings.cxx:32
 MetropolisHastings.cxx:33
 MetropolisHastings.cxx:34
 MetropolisHastings.cxx:35
 MetropolisHastings.cxx:36
 MetropolisHastings.cxx:37
 MetropolisHastings.cxx:38
 MetropolisHastings.cxx:39
 MetropolisHastings.cxx:40
 MetropolisHastings.cxx:41
 MetropolisHastings.cxx:42
 MetropolisHastings.cxx:43
 MetropolisHastings.cxx:44
 MetropolisHastings.cxx:45
 MetropolisHastings.cxx:46
 MetropolisHastings.cxx:47
 MetropolisHastings.cxx:48
 MetropolisHastings.cxx:49
 MetropolisHastings.cxx:50
 MetropolisHastings.cxx:51
 MetropolisHastings.cxx:52
 MetropolisHastings.cxx:53
 MetropolisHastings.cxx:54
 MetropolisHastings.cxx:55
 MetropolisHastings.cxx:56
 MetropolisHastings.cxx:57
 MetropolisHastings.cxx:58
 MetropolisHastings.cxx:59
 MetropolisHastings.cxx:60
 MetropolisHastings.cxx:61
 MetropolisHastings.cxx:62
 MetropolisHastings.cxx:63
 MetropolisHastings.cxx:64
 MetropolisHastings.cxx:65
 MetropolisHastings.cxx:66
 MetropolisHastings.cxx:67
 MetropolisHastings.cxx:68
 MetropolisHastings.cxx:69
 MetropolisHastings.cxx:70
 MetropolisHastings.cxx:71
 MetropolisHastings.cxx:72
 MetropolisHastings.cxx:73
 MetropolisHastings.cxx:74
 MetropolisHastings.cxx:75
 MetropolisHastings.cxx:76
 MetropolisHastings.cxx:77
 MetropolisHastings.cxx:78
 MetropolisHastings.cxx:79
 MetropolisHastings.cxx:80
 MetropolisHastings.cxx:81
 MetropolisHastings.cxx:82
 MetropolisHastings.cxx:83
 MetropolisHastings.cxx:84
 MetropolisHastings.cxx:85
 MetropolisHastings.cxx:86
 MetropolisHastings.cxx:87
 MetropolisHastings.cxx:88
 MetropolisHastings.cxx:89
 MetropolisHastings.cxx:90
 MetropolisHastings.cxx:91
 MetropolisHastings.cxx:92
 MetropolisHastings.cxx:93
 MetropolisHastings.cxx:94
 MetropolisHastings.cxx:95
 MetropolisHastings.cxx:96
 MetropolisHastings.cxx:97
 MetropolisHastings.cxx:98
 MetropolisHastings.cxx:99
 MetropolisHastings.cxx:100
 MetropolisHastings.cxx:101
 MetropolisHastings.cxx:102
 MetropolisHastings.cxx:103
 MetropolisHastings.cxx:104
 MetropolisHastings.cxx:105
 MetropolisHastings.cxx:106
 MetropolisHastings.cxx:107
 MetropolisHastings.cxx:108
 MetropolisHastings.cxx:109
 MetropolisHastings.cxx:110
 MetropolisHastings.cxx:111
 MetropolisHastings.cxx:112
 MetropolisHastings.cxx:113
 MetropolisHastings.cxx:114
 MetropolisHastings.cxx:115
 MetropolisHastings.cxx:116
 MetropolisHastings.cxx:117
 MetropolisHastings.cxx:118
 MetropolisHastings.cxx:119
 MetropolisHastings.cxx:120
 MetropolisHastings.cxx:121
 MetropolisHastings.cxx:122
 MetropolisHastings.cxx:123
 MetropolisHastings.cxx:124
 MetropolisHastings.cxx:125
 MetropolisHastings.cxx:126
 MetropolisHastings.cxx:127
 MetropolisHastings.cxx:128
 MetropolisHastings.cxx:129
 MetropolisHastings.cxx:130
 MetropolisHastings.cxx:131
 MetropolisHastings.cxx:132
 MetropolisHastings.cxx:133
 MetropolisHastings.cxx:134
 MetropolisHastings.cxx:135
 MetropolisHastings.cxx:136
 MetropolisHastings.cxx:137
 MetropolisHastings.cxx:138
 MetropolisHastings.cxx:139
 MetropolisHastings.cxx:140
 MetropolisHastings.cxx:141
 MetropolisHastings.cxx:142
 MetropolisHastings.cxx:143
 MetropolisHastings.cxx:144
 MetropolisHastings.cxx:145
 MetropolisHastings.cxx:146
 MetropolisHastings.cxx:147
 MetropolisHastings.cxx:148
 MetropolisHastings.cxx:149
 MetropolisHastings.cxx:150
 MetropolisHastings.cxx:151
 MetropolisHastings.cxx:152
 MetropolisHastings.cxx:153
 MetropolisHastings.cxx:154
 MetropolisHastings.cxx:155
 MetropolisHastings.cxx:156
 MetropolisHastings.cxx:157
 MetropolisHastings.cxx:158
 MetropolisHastings.cxx:159
 MetropolisHastings.cxx:160
 MetropolisHastings.cxx:161
 MetropolisHastings.cxx:162
 MetropolisHastings.cxx:163
 MetropolisHastings.cxx:164
 MetropolisHastings.cxx:165
 MetropolisHastings.cxx:166
 MetropolisHastings.cxx:167
 MetropolisHastings.cxx:168
 MetropolisHastings.cxx:169
 MetropolisHastings.cxx:170
 MetropolisHastings.cxx:171
 MetropolisHastings.cxx:172
 MetropolisHastings.cxx:173
 MetropolisHastings.cxx:174
 MetropolisHastings.cxx:175
 MetropolisHastings.cxx:176
 MetropolisHastings.cxx:177
 MetropolisHastings.cxx:178
 MetropolisHastings.cxx:179
 MetropolisHastings.cxx:180
 MetropolisHastings.cxx:181
 MetropolisHastings.cxx:182
 MetropolisHastings.cxx:183
 MetropolisHastings.cxx:184
 MetropolisHastings.cxx:185
 MetropolisHastings.cxx:186
 MetropolisHastings.cxx:187
 MetropolisHastings.cxx:188
 MetropolisHastings.cxx:189
 MetropolisHastings.cxx:190
 MetropolisHastings.cxx:191
 MetropolisHastings.cxx:192
 MetropolisHastings.cxx:193
 MetropolisHastings.cxx:194
 MetropolisHastings.cxx:195
 MetropolisHastings.cxx:196
 MetropolisHastings.cxx:197
 MetropolisHastings.cxx:198
 MetropolisHastings.cxx:199
 MetropolisHastings.cxx:200
 MetropolisHastings.cxx:201
 MetropolisHastings.cxx:202
 MetropolisHastings.cxx:203
 MetropolisHastings.cxx:204
 MetropolisHastings.cxx:205
 MetropolisHastings.cxx:206
 MetropolisHastings.cxx:207
 MetropolisHastings.cxx:208
 MetropolisHastings.cxx:209
 MetropolisHastings.cxx:210
 MetropolisHastings.cxx:211
 MetropolisHastings.cxx:212
 MetropolisHastings.cxx:213
 MetropolisHastings.cxx:214
 MetropolisHastings.cxx:215
 MetropolisHastings.cxx:216
 MetropolisHastings.cxx:217
 MetropolisHastings.cxx:218
 MetropolisHastings.cxx:219
 MetropolisHastings.cxx:220
 MetropolisHastings.cxx:221
 MetropolisHastings.cxx:222
 MetropolisHastings.cxx:223
 MetropolisHastings.cxx:224
 MetropolisHastings.cxx:225
 MetropolisHastings.cxx:226
 MetropolisHastings.cxx:227
 MetropolisHastings.cxx:228
 MetropolisHastings.cxx:229
 MetropolisHastings.cxx:230
 MetropolisHastings.cxx:231
 MetropolisHastings.cxx:232
 MetropolisHastings.cxx:233
 MetropolisHastings.cxx:234
 MetropolisHastings.cxx:235
 MetropolisHastings.cxx:236
 MetropolisHastings.cxx:237
 MetropolisHastings.cxx:238
 MetropolisHastings.cxx:239
 MetropolisHastings.cxx:240
 MetropolisHastings.cxx:241
 MetropolisHastings.cxx:242
 MetropolisHastings.cxx:243
 MetropolisHastings.cxx:244
 MetropolisHastings.cxx:245
 MetropolisHastings.cxx:246
 MetropolisHastings.cxx:247
 MetropolisHastings.cxx:248
 MetropolisHastings.cxx:249
 MetropolisHastings.cxx:250
 MetropolisHastings.cxx:251
 MetropolisHastings.cxx:252
 MetropolisHastings.cxx:253
 MetropolisHastings.cxx:254
 MetropolisHastings.cxx:255
 MetropolisHastings.cxx:256
 MetropolisHastings.cxx:257
 MetropolisHastings.cxx:258
 MetropolisHastings.cxx:259
 MetropolisHastings.cxx:260
 MetropolisHastings.cxx:261
 MetropolisHastings.cxx:262
 MetropolisHastings.cxx:263
 MetropolisHastings.cxx:264
 MetropolisHastings.cxx:265
 MetropolisHastings.cxx:266
 MetropolisHastings.cxx:267
 MetropolisHastings.cxx:268
 MetropolisHastings.cxx:269
 MetropolisHastings.cxx:270
 MetropolisHastings.cxx:271
 MetropolisHastings.cxx:272
 MetropolisHastings.cxx:273
 MetropolisHastings.cxx:274
 MetropolisHastings.cxx:275
 MetropolisHastings.cxx:276
 MetropolisHastings.cxx:277
 MetropolisHastings.cxx:278
 MetropolisHastings.cxx:279
 MetropolisHastings.cxx:280
 MetropolisHastings.cxx:281
 MetropolisHastings.cxx:282
 MetropolisHastings.cxx:283
 MetropolisHastings.cxx:284
 MetropolisHastings.cxx:285
 MetropolisHastings.cxx:286
 MetropolisHastings.cxx:287
 MetropolisHastings.cxx:288
 MetropolisHastings.cxx:289
 MetropolisHastings.cxx:290
 MetropolisHastings.cxx:291
 MetropolisHastings.cxx:292
 MetropolisHastings.cxx:293
 MetropolisHastings.cxx:294
 MetropolisHastings.cxx:295
 MetropolisHastings.cxx:296
 MetropolisHastings.cxx:297
 MetropolisHastings.cxx:298
 MetropolisHastings.cxx:299
 MetropolisHastings.cxx:300
 MetropolisHastings.cxx:301
 MetropolisHastings.cxx:302
 MetropolisHastings.cxx:303
 MetropolisHastings.cxx:304
 MetropolisHastings.cxx:305
 MetropolisHastings.cxx:306
 MetropolisHastings.cxx:307
 MetropolisHastings.cxx:308
 MetropolisHastings.cxx:309
 MetropolisHastings.cxx:310
 MetropolisHastings.cxx:311
 MetropolisHastings.cxx:312
 MetropolisHastings.cxx:313
 MetropolisHastings.cxx:314
 MetropolisHastings.cxx:315
 MetropolisHastings.cxx:316
 MetropolisHastings.cxx:317
 MetropolisHastings.cxx:318
 MetropolisHastings.cxx:319
 MetropolisHastings.cxx:320
 MetropolisHastings.cxx:321
 MetropolisHastings.cxx:322
 MetropolisHastings.cxx:323
 MetropolisHastings.cxx:324
 MetropolisHastings.cxx:325
 MetropolisHastings.cxx:326
 MetropolisHastings.cxx:327
 MetropolisHastings.cxx:328
 MetropolisHastings.cxx:329
 MetropolisHastings.cxx:330
 MetropolisHastings.cxx:331
 MetropolisHastings.cxx:332
 MetropolisHastings.cxx:333
 MetropolisHastings.cxx:334
 MetropolisHastings.cxx:335
 MetropolisHastings.cxx:336