/*****************************************************************************
 * Project: RooFit                                                           *
 * Package: RooFitModels                                                     *
 *    File: $Id$
 * Authors:                                                                  *
 *   GR, Gerhard Raven,   Nikhef & VU, Gerhard.Raven@nikhef.nl
 *                                                                           *
 * Copyright (c) 2010, Nikhef & VU. All rights reserved.
 *                                                                           *
 * Redistribution and use in source and binary forms,                        *
 * with or without modification, are permitted according to the terms        *
 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
 *****************************************************************************/

//////////////////////////////////////////////////////////////////////////////
//
// BEGIN_HTML
// END_HTML
//

#include "RooFit.h"
#include "Riostream.h"
#include <math.h>
#include <string>
#include <algorithm>

#include "RooLegendre.h"
#include "RooAbsReal.h"
#include "Math/SpecFunc.h"
#include "TMath.h"

#include "TError.h"

using namespace std;

ClassImp(RooLegendre)
;


//_____________________________________________________________________________
namespace {
    inline double a(int p, int l, int m) {
        double r = TMath::Factorial(l+m)/TMath::Factorial(m+p)/TMath::Factorial(p)/TMath::Factorial(l-m-2*p);
        r /= pow(2.,m+2*p);
        return p%2==0 ? r : -r ;
    }
}

//_____________________________________________________________________________
RooLegendre::RooLegendre() :
  _l1(1),_m1(1),_l2(0),_m2(0)
{
}

//_____________________________________________________________________________
RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l, int m) 
 : RooAbsReal(name, title)
 , _ctheta("ctheta", "ctheta", this, ctheta)
 , _l1(l),_m1(m),_l2(0),_m2(0)
{
  //TODO: for now, we assume that ctheta has a range [-1,1]
  // should map the ctheta range onto this interval, and adjust integrals...

  //TODO: we assume m>=0
  //      should map m<0 back to m>=0...
}

//_____________________________________________________________________________
RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l1, int m1, int l2, int m2) 
 : RooAbsReal(name, title)
 , _ctheta("ctheta", "ctheta", this, ctheta)
 , _l1(l1),_m1(m1),_l2(l2),_m2(m2)
{
}

//_____________________________________________________________________________
RooLegendre::RooLegendre(const RooLegendre& other, const char* name) 
    : RooAbsReal(other, name)
    , _ctheta("ctheta", this, other._ctheta)
    , _l1(other._l1), _m1(other._m1)
    , _l2(other._l2), _m2(other._m2)
{
}

//_____________________________________________________________________________
Double_t RooLegendre::evaluate() const 
{
  // TODO: check that 0<=m_i<=l_i; on the other hand, assoc_legendre already does that ;-)
  // Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0
#ifdef R__HAS_MATHMORE  
  double r = 1;
  double ctheta = std::max(-1., std::min((double)_ctheta, +1.));
  if (_l1!=0||_m1!=0) r *= ROOT::Math::assoc_legendre(_l1,_m1,ctheta);
  if (_l2!=0||_m2!=0) r *= ROOT::Math::assoc_legendre(_l2,_m2,ctheta);
  if ((_m1+_m2)%2==1) r = -r;
  return r;
#else
  throw std::string("RooLegendre: ERROR: This class require installation of the MathMore library") ;
  return 0 ;
#endif
}

//_____________________________________________________________________________
namespace {
    bool fullRange(const RooRealProxy& x ,const char* range) 
    { return range==0 || strlen(range)==0 
          || ( x.min(range) == x.min() && x.max(range) == x.max() ) ; 
    }
}
Int_t RooLegendre::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
{
  // don't support indefinite integrals...
  if (fullRange(_ctheta,rangeName) && matchArgs(allVars, analVars, _ctheta)) return 1;
  return 0;
}

//_____________________________________________________________________________
Double_t RooLegendre::analyticalIntegral(Int_t code, const char* ) const 
{
  // this was verified to match mathematica for 
  // l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
  R__ASSERT(code==1) ;
  if ( _m1==_m2 )                 return ( _l1 == _l2) ?  TMath::Factorial(_l1+_m2)/TMath::Factorial(_l1-_m1)*double(2)/(2*_l1+1) : 0.;
  if ( (_l1+_l2-_m1-_m2)%2 != 0 ) return 0; // these combinations are odd under x -> -x

  // from B.R. Wong, "On the overlap integral of associated Legendre Polynomials" 1998 J. Phys. A: Math. Gen. 31 1101
  // TODO: update to the result of 
  //       H. A. Mavromatis
  //       "A single-sum expression for the overlap integral of two associated Legendre polynomials"
  //       1999 J. Phys. A: Math. Gen. 32 2601
  //       http://iopscience.iop.org/0305-4470/32/13/011/pdf/0305-4470_32_13_011.pdf
  //       For that we need Wigner 3-j, which Lorenzo has added for Root 5.28... (note: check Condon-Shortly convention in this paper!)
  double r=0;
  for (int p1=0; 2*p1 <= _l1-_m1 ;++p1) {
    double a1 = a(p1,_l1,_m1);
    for (int p2=0; 2*p2 <= _l2-_m2 ; ++p2) {
       double a2 = a(p2,_l2,_m2);
       r+= a1*a2*TMath::Gamma( double(_l1+_l2-_m1-_m2-2*p1-2*p2+1)/2 )*TMath::Gamma( double(_m1+_m2+2*p1+2*p2+2)/2 );
    }
  }
  r /= TMath::Gamma( double(_l1+_l2+3)/2 );

  if ((_m1+_m2)%2==1) r = -r;
  return r;
}

Int_t RooLegendre::getMaxVal( const RooArgSet& /*vars*/) const {
    if (_m1==0&&_m2==0) return 1;
    // does anyone know the analytical expression for the  max values in case m!=0??
    if (_l1<3&&_l2<3) return 1;
    return 0;
}

namespace {
    inline double maxSingle(int i, int j) {
        R__ASSERT(j<=i);
        //   x0 : 1 (ordinary Legendre)
        if (j==0) return 1;
        R__ASSERT(i<3);
        //   11: 1
        if (i<2) return 1;
        //   21: 3   22: 3
        static const double m2[3] = { 3,3 };
        return m2[j-1];
    }
}
Double_t RooLegendre::maxVal( Int_t /*code*/) const {
    return maxSingle(_l1,_m1)*maxSingle(_l2,_m2);
}
 RooLegendre.cxx:1
 RooLegendre.cxx:2
 RooLegendre.cxx:3
 RooLegendre.cxx:4
 RooLegendre.cxx:5
 RooLegendre.cxx:6
 RooLegendre.cxx:7
 RooLegendre.cxx:8
 RooLegendre.cxx:9
 RooLegendre.cxx:10
 RooLegendre.cxx:11
 RooLegendre.cxx:12
 RooLegendre.cxx:13
 RooLegendre.cxx:14
 RooLegendre.cxx:15
 RooLegendre.cxx:16
 RooLegendre.cxx:17
 RooLegendre.cxx:18
 RooLegendre.cxx:19
 RooLegendre.cxx:20
 RooLegendre.cxx:21
 RooLegendre.cxx:22
 RooLegendre.cxx:23
 RooLegendre.cxx:24
 RooLegendre.cxx:25
 RooLegendre.cxx:26
 RooLegendre.cxx:27
 RooLegendre.cxx:28
 RooLegendre.cxx:29
 RooLegendre.cxx:30
 RooLegendre.cxx:31
 RooLegendre.cxx:32
 RooLegendre.cxx:33
 RooLegendre.cxx:34
 RooLegendre.cxx:35
 RooLegendre.cxx:36
 RooLegendre.cxx:37
 RooLegendre.cxx:38
 RooLegendre.cxx:39
 RooLegendre.cxx:40
 RooLegendre.cxx:41
 RooLegendre.cxx:42
 RooLegendre.cxx:43
 RooLegendre.cxx:44
 RooLegendre.cxx:45
 RooLegendre.cxx:46
 RooLegendre.cxx:47
 RooLegendre.cxx:48
 RooLegendre.cxx:49
 RooLegendre.cxx:50
 RooLegendre.cxx:51
 RooLegendre.cxx:52
 RooLegendre.cxx:53
 RooLegendre.cxx:54
 RooLegendre.cxx:55
 RooLegendre.cxx:56
 RooLegendre.cxx:57
 RooLegendre.cxx:58
 RooLegendre.cxx:59
 RooLegendre.cxx:60
 RooLegendre.cxx:61
 RooLegendre.cxx:62
 RooLegendre.cxx:63
 RooLegendre.cxx:64
 RooLegendre.cxx:65
 RooLegendre.cxx:66
 RooLegendre.cxx:67
 RooLegendre.cxx:68
 RooLegendre.cxx:69
 RooLegendre.cxx:70
 RooLegendre.cxx:71
 RooLegendre.cxx:72
 RooLegendre.cxx:73
 RooLegendre.cxx:74
 RooLegendre.cxx:75
 RooLegendre.cxx:76
 RooLegendre.cxx:77
 RooLegendre.cxx:78
 RooLegendre.cxx:79
 RooLegendre.cxx:80
 RooLegendre.cxx:81
 RooLegendre.cxx:82
 RooLegendre.cxx:83
 RooLegendre.cxx:84
 RooLegendre.cxx:85
 RooLegendre.cxx:86
 RooLegendre.cxx:87
 RooLegendre.cxx:88
 RooLegendre.cxx:89
 RooLegendre.cxx:90
 RooLegendre.cxx:91
 RooLegendre.cxx:92
 RooLegendre.cxx:93
 RooLegendre.cxx:94
 RooLegendre.cxx:95
 RooLegendre.cxx:96
 RooLegendre.cxx:97
 RooLegendre.cxx:98
 RooLegendre.cxx:99
 RooLegendre.cxx:100
 RooLegendre.cxx:101
 RooLegendre.cxx:102
 RooLegendre.cxx:103
 RooLegendre.cxx:104
 RooLegendre.cxx:105
 RooLegendre.cxx:106
 RooLegendre.cxx:107
 RooLegendre.cxx:108
 RooLegendre.cxx:109
 RooLegendre.cxx:110
 RooLegendre.cxx:111
 RooLegendre.cxx:112
 RooLegendre.cxx:113
 RooLegendre.cxx:114
 RooLegendre.cxx:115
 RooLegendre.cxx:116
 RooLegendre.cxx:117
 RooLegendre.cxx:118
 RooLegendre.cxx:119
 RooLegendre.cxx:120
 RooLegendre.cxx:121
 RooLegendre.cxx:122
 RooLegendre.cxx:123
 RooLegendre.cxx:124
 RooLegendre.cxx:125
 RooLegendre.cxx:126
 RooLegendre.cxx:127
 RooLegendre.cxx:128
 RooLegendre.cxx:129
 RooLegendre.cxx:130
 RooLegendre.cxx:131
 RooLegendre.cxx:132
 RooLegendre.cxx:133
 RooLegendre.cxx:134
 RooLegendre.cxx:135
 RooLegendre.cxx:136
 RooLegendre.cxx:137
 RooLegendre.cxx:138
 RooLegendre.cxx:139
 RooLegendre.cxx:140
 RooLegendre.cxx:141
 RooLegendre.cxx:142
 RooLegendre.cxx:143
 RooLegendre.cxx:144
 RooLegendre.cxx:145
 RooLegendre.cxx:146
 RooLegendre.cxx:147
 RooLegendre.cxx:148
 RooLegendre.cxx:149
 RooLegendre.cxx:150
 RooLegendre.cxx:151
 RooLegendre.cxx:152
 RooLegendre.cxx:153
 RooLegendre.cxx:154
 RooLegendre.cxx:155
 RooLegendre.cxx:156
 RooLegendre.cxx:157
 RooLegendre.cxx:158
 RooLegendre.cxx:159
 RooLegendre.cxx:160
 RooLegendre.cxx:161
 RooLegendre.cxx:162
 RooLegendre.cxx:163
 RooLegendre.cxx:164
 RooLegendre.cxx:165
 RooLegendre.cxx:166
 RooLegendre.cxx:167
 RooLegendre.cxx:168
 RooLegendre.cxx:169