/***************************************************************************** 
  * Project: RooFit                                                           * 
  * @(#)root/roofit:$Id$ *
  *                                                                           * 
  * RooFit Lognormal PDF                                                      *
  *                                                                           * 
  * Author: Gregory Schott and Stefan Schmitz                                 *
  *                                                                           * 
  *****************************************************************************/ 

//////////////////////////////////////////////////////////////////////////////
//
// BEGIN_HTML
// RooFit Lognormal PDF. The two parameters are:
//  - m0: the median of the distribution
//  - k=exp(sigma): sigma is called the shape parameter in the TMath parametrization
//
//   Lognormal(x,m0,k) = 1/(sqrt(2*pi)*ln(k)*x)*exp(-ln^2(x/m0)/(2*ln^2(k)))
//
// The parametrization here is physics driven and differs from the ROOT::Math::lognormal_pdf(x,m,s,x0) with:
//  - m = log(m0)
//  - s = log(k)
//  - x0 = 0
// END_HTML


#include "RooFit.h"

#include "Riostream.h"
#include "Riostream.h"
#include <math.h>

#include "RooLognormal.h"
#include "RooAbsReal.h"
#include "RooRealVar.h"
#include "RooRandom.h"
#include "RooMath.h"
#include "TMath.h"

#include <Math/SpecFuncMathCore.h>
#include <Math/PdfFuncMathCore.h>
#include <Math/ProbFuncMathCore.h>

#include "TError.h"

using namespace std;

ClassImp(RooLognormal)


//_____________________________________________________________________________
RooLognormal::RooLognormal(const char *name, const char *title,
			 RooAbsReal& _x, RooAbsReal& _m0,
			 RooAbsReal& _k) :
  RooAbsPdf(name,title),
  x("x","Observable",this,_x),
  m0("m0","m0",this,_m0),
  k("k","k",this,_k)
{
}



//_____________________________________________________________________________
RooLognormal::RooLognormal(const RooLognormal& other, const char* name) : 
  RooAbsPdf(other,name), x("x",this,other.x), m0("m0",this,other.m0),
  k("k",this,other.k)
{
}



//_____________________________________________________________________________
Double_t RooLognormal::evaluate() const
{
  // ln(k)<1 would correspond to sigma < 0 in the parametrization
  // resulting by transforming a normal random variable in its
  // standard parametrization to a lognormal random variable
  // => treat ln(k) as -ln(k) for k<1

  Double_t xv = x;
  Double_t ln_k = TMath::Abs(TMath::Log(k));
  Double_t ln_m0 = TMath::Log(m0);
  Double_t x0 = 0;

  Double_t ret = ROOT::Math::lognormal_pdf(xv,ln_m0,ln_k,x0);
  return ret ;
}



//_____________________________________________________________________________
Int_t RooLognormal::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const 
{
  if (matchArgs(allVars,analVars,x)) return 1 ;
  return 0 ;
}



//_____________________________________________________________________________
Double_t RooLognormal::analyticalIntegral(Int_t code, const char* rangeName) const 
{
  R__ASSERT(code==1) ;

  static const Double_t root2 = sqrt(2.) ;

  Double_t ln_k = TMath::Abs(TMath::Log(k));
  Double_t ret = 0.5*( RooMath::erf( TMath::Log(x.max(rangeName)/m0)/(root2*ln_k) ) - RooMath::erf( TMath::Log(x.min(rangeName)/m0)/(root2*ln_k) ) ) ;

  return ret ;
}




//_____________________________________________________________________________
Int_t RooLognormal::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
{
  if (matchArgs(directVars,generateVars,x)) return 1 ;  
  return 0 ;
}



//_____________________________________________________________________________
void RooLognormal::generateEvent(Int_t code)
{
  R__ASSERT(code==1) ;

  Double_t xgen ;
  while(1) {    
    xgen = TMath::Exp(RooRandom::randomGenerator()->Gaus(TMath::Log(m0),TMath::Log(k)));
    if (xgen<=x.max() && xgen>=x.min()) {
      x = xgen ;
      break;
    }
  }

  return;
}


 RooLognormal.cxx:1
 RooLognormal.cxx:2
 RooLognormal.cxx:3
 RooLognormal.cxx:4
 RooLognormal.cxx:5
 RooLognormal.cxx:6
 RooLognormal.cxx:7
 RooLognormal.cxx:8
 RooLognormal.cxx:9
 RooLognormal.cxx:10
 RooLognormal.cxx:11
 RooLognormal.cxx:12
 RooLognormal.cxx:13
 RooLognormal.cxx:14
 RooLognormal.cxx:15
 RooLognormal.cxx:16
 RooLognormal.cxx:17
 RooLognormal.cxx:18
 RooLognormal.cxx:19
 RooLognormal.cxx:20
 RooLognormal.cxx:21
 RooLognormal.cxx:22
 RooLognormal.cxx:23
 RooLognormal.cxx:24
 RooLognormal.cxx:25
 RooLognormal.cxx:26
 RooLognormal.cxx:27
 RooLognormal.cxx:28
 RooLognormal.cxx:29
 RooLognormal.cxx:30
 RooLognormal.cxx:31
 RooLognormal.cxx:32
 RooLognormal.cxx:33
 RooLognormal.cxx:34
 RooLognormal.cxx:35
 RooLognormal.cxx:36
 RooLognormal.cxx:37
 RooLognormal.cxx:38
 RooLognormal.cxx:39
 RooLognormal.cxx:40
 RooLognormal.cxx:41
 RooLognormal.cxx:42
 RooLognormal.cxx:43
 RooLognormal.cxx:44
 RooLognormal.cxx:45
 RooLognormal.cxx:46
 RooLognormal.cxx:47
 RooLognormal.cxx:48
 RooLognormal.cxx:49
 RooLognormal.cxx:50
 RooLognormal.cxx:51
 RooLognormal.cxx:52
 RooLognormal.cxx:53
 RooLognormal.cxx:54
 RooLognormal.cxx:55
 RooLognormal.cxx:56
 RooLognormal.cxx:57
 RooLognormal.cxx:58
 RooLognormal.cxx:59
 RooLognormal.cxx:60
 RooLognormal.cxx:61
 RooLognormal.cxx:62
 RooLognormal.cxx:63
 RooLognormal.cxx:64
 RooLognormal.cxx:65
 RooLognormal.cxx:66
 RooLognormal.cxx:67
 RooLognormal.cxx:68
 RooLognormal.cxx:69
 RooLognormal.cxx:70
 RooLognormal.cxx:71
 RooLognormal.cxx:72
 RooLognormal.cxx:73
 RooLognormal.cxx:74
 RooLognormal.cxx:75
 RooLognormal.cxx:76
 RooLognormal.cxx:77
 RooLognormal.cxx:78
 RooLognormal.cxx:79
 RooLognormal.cxx:80
 RooLognormal.cxx:81
 RooLognormal.cxx:82
 RooLognormal.cxx:83
 RooLognormal.cxx:84
 RooLognormal.cxx:85
 RooLognormal.cxx:86
 RooLognormal.cxx:87
 RooLognormal.cxx:88
 RooLognormal.cxx:89
 RooLognormal.cxx:90
 RooLognormal.cxx:91
 RooLognormal.cxx:92
 RooLognormal.cxx:93
 RooLognormal.cxx:94
 RooLognormal.cxx:95
 RooLognormal.cxx:96
 RooLognormal.cxx:97
 RooLognormal.cxx:98
 RooLognormal.cxx:99
 RooLognormal.cxx:100
 RooLognormal.cxx:101
 RooLognormal.cxx:102
 RooLognormal.cxx:103
 RooLognormal.cxx:104
 RooLognormal.cxx:105
 RooLognormal.cxx:106
 RooLognormal.cxx:107
 RooLognormal.cxx:108
 RooLognormal.cxx:109
 RooLognormal.cxx:110
 RooLognormal.cxx:111
 RooLognormal.cxx:112
 RooLognormal.cxx:113
 RooLognormal.cxx:114
 RooLognormal.cxx:115
 RooLognormal.cxx:116
 RooLognormal.cxx:117
 RooLognormal.cxx:118
 RooLognormal.cxx:119
 RooLognormal.cxx:120
 RooLognormal.cxx:121
 RooLognormal.cxx:122
 RooLognormal.cxx:123
 RooLognormal.cxx:124
 RooLognormal.cxx:125
 RooLognormal.cxx:126
 RooLognormal.cxx:127
 RooLognormal.cxx:128
 RooLognormal.cxx:129
 RooLognormal.cxx:130
 RooLognormal.cxx:131
 RooLognormal.cxx:132
 RooLognormal.cxx:133
 RooLognormal.cxx:134
 RooLognormal.cxx:135
 RooLognormal.cxx:136
 RooLognormal.cxx:137
 RooLognormal.cxx:138
 RooLognormal.cxx:139
 RooLognormal.cxx:140
 RooLognormal.cxx:141
 RooLognormal.cxx:142
 RooLognormal.cxx:143