/*****************************************************************************
 * Project: RooFit                                                           *
 * Package: RooFitModels                                                     *
 * @(#)root/roofit:$Id$
 * Authors:                                                                  *
 *   JGS, Jim Smith    , University of Colorado, jgsmith@pizero.colorado.edu *
 * History:
 *   15-Aug-2002 JGS Created initial version
 *   11-Sep-2002 JGS Mods to introduce mu (Mirna van Hoek, JGS, Nick Danielson)
 *                                                                           *
 * Copyright (c) 2000-2005, Regents of the University of California,         *
 *                          University of Colorado                           *
 *                          and Stanford University. 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
// Implement standard CP physics model with S and C (no mention of lambda)
// Suitably stolen and modified from RooBCPEffDecay
// END_HTML
//

#include "RooFit.h"

#include "Riostream.h"
#include "Riostream.h"
#include "RooRealVar.h"
#include "RooRandom.h"
#include "RooBCPGenDecay.h"
#include "RooRealIntegral.h"

using namespace std;

ClassImp(RooBCPGenDecay) 
;



//_____________________________________________________________________________
RooBCPGenDecay::RooBCPGenDecay(const char *name, const char *title, 
			       RooRealVar& t, RooAbsCategory& tag,
			       RooAbsReal& tau, RooAbsReal& dm,
			       RooAbsReal& avgMistag, 
			       RooAbsReal& a, RooAbsReal& b,
			       RooAbsReal& delMistag,
                               RooAbsReal& mu,
			       const RooResolutionModel& model, DecayType type) :
  RooAbsAnaConvPdf(name,title,model,t), 
  _avgC("C","Coefficient of cos term",this,a),
  _avgS("S","Coefficient of cos term",this,b),
  _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
  _delMistag("delMistag","Delta mistag rate",this,delMistag),  
  _mu("mu","Tagg efficiency difference",this,mu),  
  _t("t","time",this,t),
  _tau("tau","decay time",this,tau),
  _dm("dm","mixing frequency",this,dm),
  _tag("tag","CP state",this,tag),
  _genB0Frac(0),
  _type(type)
{
  // Constructor
  switch(type) {
  case SingleSided:
    _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
    _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
    _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
    break ;
  case Flipped:
    _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
    _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
    _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
    break ;
  case DoubleSided:
    _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
    _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
    _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
    break ;
  }
}



//_____________________________________________________________________________
RooBCPGenDecay::RooBCPGenDecay(const RooBCPGenDecay& other, const char* name) : 
  RooAbsAnaConvPdf(other,name), 
  _avgC("C",this,other._avgC),
  _avgS("S",this,other._avgS),
  _avgMistag("avgMistag",this,other._avgMistag),
  _delMistag("delMistag",this,other._delMistag),
  _mu("mu",this,other._mu),
  _t("t",this,other._t),
  _tau("tau",this,other._tau),
  _dm("dm",this,other._dm),
  _tag("tag",this,other._tag),
  _genB0Frac(other._genB0Frac),
  _type(other._type),
  _basisExp(other._basisExp),
  _basisSin(other._basisSin),
  _basisCos(other._basisCos)
{
  // Copy constructor
}



//_____________________________________________________________________________
RooBCPGenDecay::~RooBCPGenDecay()
{
  // Destructor
}



//_____________________________________________________________________________
Double_t RooBCPGenDecay::coefficient(Int_t basisIndex) const 
{
  // B0    : _tag = +1 
  // B0bar : _tag = -1 

  if (basisIndex==_basisExp) {
    //exp term: (1 -/+ dw + mu*_tag*w)
    return (1 - _tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag)) ;
    // =    1 + _tag*deltaDil/2 + _mu*avgDil
  }

  if (basisIndex==_basisSin) {
    //sin term: (+/- (1-2w) + _mu*(1 -/+ delw))*S
    return (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS ;
    // =   (_tag*avgDil + _mu*(1 + tag*deltaDil/2)) * S
  }
  
  if (basisIndex==_basisCos) {
    //cos term: -(+/- (1-2w) + _mu*(1 -/+ delw))*C
    return -1.*(_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC ;
    // =   -(_tag*avgDil + _mu*(1 + _tag*deltaDil/2) )* C
  } 
  
  return 0 ;
}



//_____________________________________________________________________________
Int_t RooBCPGenDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
{
  if (rangeName) return 0 ;
  if (matchArgs(allVars,analVars,_tag)) return 1 ;
  return 0 ;
}



//_____________________________________________________________________________
Double_t RooBCPGenDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const 
{
  switch(code) {
    // No integration
  case 0: return coefficient(basisIndex) ;

    // Integration over 'tag'
  case 1:
    if (basisIndex==_basisExp) {
      return 2 ;
    }
    
    if (basisIndex==_basisSin) {
      return 2*_mu*_avgS ;
    }
    if (basisIndex==_basisCos) {
      return -2*_mu*_avgC ;
    }
    break ;
    
  default:
    assert(0) ;
  }
    
  return 0 ;
}



//_____________________________________________________________________________
Int_t RooBCPGenDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
{
  if (staticInitOK) {
    if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
  }
  if (matchArgs(directVars,generateVars,_t)) return 1 ;  
  return 0 ;
}



//_____________________________________________________________________________
void RooBCPGenDecay::initGenerator(Int_t code)
{
  if (code==2) {
    // Calculate the fraction of mixed events to generate
    Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
    _tag = 1 ;
    Double_t b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
    _genB0Frac = b0Int/sumInt ;
  }  
}



//_____________________________________________________________________________
void RooBCPGenDecay::generateEvent(Int_t code)
{
  // Generate mix-state dependent
  if (code==2) {
    Double_t rand = RooRandom::uniform() ;
    _tag = (rand<=_genB0Frac) ? 1 : -1 ;
  }

  // Generate delta-t dependent
  while(1) {
    Double_t rand = RooRandom::uniform() ;
    Double_t tval(0) ;

    switch(_type) {
    case SingleSided:
      tval = -_tau*log(rand);
      break ;
    case Flipped:
      tval= +_tau*log(rand);
      break ;
    case DoubleSided:
      tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
      break ;
    }

    // Accept event if T is in generated range
    Double_t maxDil = 1.0 ;
// 2 in next line is conservative and inefficient - allows for delMistag=1!
    Double_t maxAcceptProb = 2 + fabs(maxDil*_avgS) + fabs(maxDil*_avgC);        
    Double_t acceptProb    = (1-_tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag)) 
                           + (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS*sin(_dm*tval) 
                           - (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC*cos(_dm*tval);

    Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
    
    if (tval<_t.max() && tval>_t.min() && accept) {
      _t = tval ;
      break ;
    }
  }
  
}

 RooBCPGenDecay.cxx:1
 RooBCPGenDecay.cxx:2
 RooBCPGenDecay.cxx:3
 RooBCPGenDecay.cxx:4
 RooBCPGenDecay.cxx:5
 RooBCPGenDecay.cxx:6
 RooBCPGenDecay.cxx:7
 RooBCPGenDecay.cxx:8
 RooBCPGenDecay.cxx:9
 RooBCPGenDecay.cxx:10
 RooBCPGenDecay.cxx:11
 RooBCPGenDecay.cxx:12
 RooBCPGenDecay.cxx:13
 RooBCPGenDecay.cxx:14
 RooBCPGenDecay.cxx:15
 RooBCPGenDecay.cxx:16
 RooBCPGenDecay.cxx:17
 RooBCPGenDecay.cxx:18
 RooBCPGenDecay.cxx:19
 RooBCPGenDecay.cxx:20
 RooBCPGenDecay.cxx:21
 RooBCPGenDecay.cxx:22
 RooBCPGenDecay.cxx:23
 RooBCPGenDecay.cxx:24
 RooBCPGenDecay.cxx:25
 RooBCPGenDecay.cxx:26
 RooBCPGenDecay.cxx:27
 RooBCPGenDecay.cxx:28
 RooBCPGenDecay.cxx:29
 RooBCPGenDecay.cxx:30
 RooBCPGenDecay.cxx:31
 RooBCPGenDecay.cxx:32
 RooBCPGenDecay.cxx:33
 RooBCPGenDecay.cxx:34
 RooBCPGenDecay.cxx:35
 RooBCPGenDecay.cxx:36
 RooBCPGenDecay.cxx:37
 RooBCPGenDecay.cxx:38
 RooBCPGenDecay.cxx:39
 RooBCPGenDecay.cxx:40
 RooBCPGenDecay.cxx:41
 RooBCPGenDecay.cxx:42
 RooBCPGenDecay.cxx:43
 RooBCPGenDecay.cxx:44
 RooBCPGenDecay.cxx:45
 RooBCPGenDecay.cxx:46
 RooBCPGenDecay.cxx:47
 RooBCPGenDecay.cxx:48
 RooBCPGenDecay.cxx:49
 RooBCPGenDecay.cxx:50
 RooBCPGenDecay.cxx:51
 RooBCPGenDecay.cxx:52
 RooBCPGenDecay.cxx:53
 RooBCPGenDecay.cxx:54
 RooBCPGenDecay.cxx:55
 RooBCPGenDecay.cxx:56
 RooBCPGenDecay.cxx:57
 RooBCPGenDecay.cxx:58
 RooBCPGenDecay.cxx:59
 RooBCPGenDecay.cxx:60
 RooBCPGenDecay.cxx:61
 RooBCPGenDecay.cxx:62
 RooBCPGenDecay.cxx:63
 RooBCPGenDecay.cxx:64
 RooBCPGenDecay.cxx:65
 RooBCPGenDecay.cxx:66
 RooBCPGenDecay.cxx:67
 RooBCPGenDecay.cxx:68
 RooBCPGenDecay.cxx:69
 RooBCPGenDecay.cxx:70
 RooBCPGenDecay.cxx:71
 RooBCPGenDecay.cxx:72
 RooBCPGenDecay.cxx:73
 RooBCPGenDecay.cxx:74
 RooBCPGenDecay.cxx:75
 RooBCPGenDecay.cxx:76
 RooBCPGenDecay.cxx:77
 RooBCPGenDecay.cxx:78
 RooBCPGenDecay.cxx:79
 RooBCPGenDecay.cxx:80
 RooBCPGenDecay.cxx:81
 RooBCPGenDecay.cxx:82
 RooBCPGenDecay.cxx:83
 RooBCPGenDecay.cxx:84
 RooBCPGenDecay.cxx:85
 RooBCPGenDecay.cxx:86
 RooBCPGenDecay.cxx:87
 RooBCPGenDecay.cxx:88
 RooBCPGenDecay.cxx:89
 RooBCPGenDecay.cxx:90
 RooBCPGenDecay.cxx:91
 RooBCPGenDecay.cxx:92
 RooBCPGenDecay.cxx:93
 RooBCPGenDecay.cxx:94
 RooBCPGenDecay.cxx:95
 RooBCPGenDecay.cxx:96
 RooBCPGenDecay.cxx:97
 RooBCPGenDecay.cxx:98
 RooBCPGenDecay.cxx:99
 RooBCPGenDecay.cxx:100
 RooBCPGenDecay.cxx:101
 RooBCPGenDecay.cxx:102
 RooBCPGenDecay.cxx:103
 RooBCPGenDecay.cxx:104
 RooBCPGenDecay.cxx:105
 RooBCPGenDecay.cxx:106
 RooBCPGenDecay.cxx:107
 RooBCPGenDecay.cxx:108
 RooBCPGenDecay.cxx:109
 RooBCPGenDecay.cxx:110
 RooBCPGenDecay.cxx:111
 RooBCPGenDecay.cxx:112
 RooBCPGenDecay.cxx:113
 RooBCPGenDecay.cxx:114
 RooBCPGenDecay.cxx:115
 RooBCPGenDecay.cxx:116
 RooBCPGenDecay.cxx:117
 RooBCPGenDecay.cxx:118
 RooBCPGenDecay.cxx:119
 RooBCPGenDecay.cxx:120
 RooBCPGenDecay.cxx:121
 RooBCPGenDecay.cxx:122
 RooBCPGenDecay.cxx:123
 RooBCPGenDecay.cxx:124
 RooBCPGenDecay.cxx:125
 RooBCPGenDecay.cxx:126
 RooBCPGenDecay.cxx:127
 RooBCPGenDecay.cxx:128
 RooBCPGenDecay.cxx:129
 RooBCPGenDecay.cxx:130
 RooBCPGenDecay.cxx:131
 RooBCPGenDecay.cxx:132
 RooBCPGenDecay.cxx:133
 RooBCPGenDecay.cxx:134
 RooBCPGenDecay.cxx:135
 RooBCPGenDecay.cxx:136
 RooBCPGenDecay.cxx:137
 RooBCPGenDecay.cxx:138
 RooBCPGenDecay.cxx:139
 RooBCPGenDecay.cxx:140
 RooBCPGenDecay.cxx:141
 RooBCPGenDecay.cxx:142
 RooBCPGenDecay.cxx:143
 RooBCPGenDecay.cxx:144
 RooBCPGenDecay.cxx:145
 RooBCPGenDecay.cxx:146
 RooBCPGenDecay.cxx:147
 RooBCPGenDecay.cxx:148
 RooBCPGenDecay.cxx:149
 RooBCPGenDecay.cxx:150
 RooBCPGenDecay.cxx:151
 RooBCPGenDecay.cxx:152
 RooBCPGenDecay.cxx:153
 RooBCPGenDecay.cxx:154
 RooBCPGenDecay.cxx:155
 RooBCPGenDecay.cxx:156
 RooBCPGenDecay.cxx:157
 RooBCPGenDecay.cxx:158
 RooBCPGenDecay.cxx:159
 RooBCPGenDecay.cxx:160
 RooBCPGenDecay.cxx:161
 RooBCPGenDecay.cxx:162
 RooBCPGenDecay.cxx:163
 RooBCPGenDecay.cxx:164
 RooBCPGenDecay.cxx:165
 RooBCPGenDecay.cxx:166
 RooBCPGenDecay.cxx:167
 RooBCPGenDecay.cxx:168
 RooBCPGenDecay.cxx:169
 RooBCPGenDecay.cxx:170
 RooBCPGenDecay.cxx:171
 RooBCPGenDecay.cxx:172
 RooBCPGenDecay.cxx:173
 RooBCPGenDecay.cxx:174
 RooBCPGenDecay.cxx:175
 RooBCPGenDecay.cxx:176
 RooBCPGenDecay.cxx:177
 RooBCPGenDecay.cxx:178
 RooBCPGenDecay.cxx:179
 RooBCPGenDecay.cxx:180
 RooBCPGenDecay.cxx:181
 RooBCPGenDecay.cxx:182
 RooBCPGenDecay.cxx:183
 RooBCPGenDecay.cxx:184
 RooBCPGenDecay.cxx:185
 RooBCPGenDecay.cxx:186
 RooBCPGenDecay.cxx:187
 RooBCPGenDecay.cxx:188
 RooBCPGenDecay.cxx:189
 RooBCPGenDecay.cxx:190
 RooBCPGenDecay.cxx:191
 RooBCPGenDecay.cxx:192
 RooBCPGenDecay.cxx:193
 RooBCPGenDecay.cxx:194
 RooBCPGenDecay.cxx:195
 RooBCPGenDecay.cxx:196
 RooBCPGenDecay.cxx:197
 RooBCPGenDecay.cxx:198
 RooBCPGenDecay.cxx:199
 RooBCPGenDecay.cxx:200
 RooBCPGenDecay.cxx:201
 RooBCPGenDecay.cxx:202
 RooBCPGenDecay.cxx:203
 RooBCPGenDecay.cxx:204
 RooBCPGenDecay.cxx:205
 RooBCPGenDecay.cxx:206
 RooBCPGenDecay.cxx:207
 RooBCPGenDecay.cxx:208
 RooBCPGenDecay.cxx:209
 RooBCPGenDecay.cxx:210
 RooBCPGenDecay.cxx:211
 RooBCPGenDecay.cxx:212
 RooBCPGenDecay.cxx:213
 RooBCPGenDecay.cxx:214
 RooBCPGenDecay.cxx:215
 RooBCPGenDecay.cxx:216
 RooBCPGenDecay.cxx:217
 RooBCPGenDecay.cxx:218
 RooBCPGenDecay.cxx:219
 RooBCPGenDecay.cxx:220
 RooBCPGenDecay.cxx:221
 RooBCPGenDecay.cxx:222
 RooBCPGenDecay.cxx:223
 RooBCPGenDecay.cxx:224
 RooBCPGenDecay.cxx:225
 RooBCPGenDecay.cxx:226
 RooBCPGenDecay.cxx:227
 RooBCPGenDecay.cxx:228
 RooBCPGenDecay.cxx:229
 RooBCPGenDecay.cxx:230
 RooBCPGenDecay.cxx:231
 RooBCPGenDecay.cxx:232
 RooBCPGenDecay.cxx:233
 RooBCPGenDecay.cxx:234
 RooBCPGenDecay.cxx:235
 RooBCPGenDecay.cxx:236
 RooBCPGenDecay.cxx:237
 RooBCPGenDecay.cxx:238
 RooBCPGenDecay.cxx:239
 RooBCPGenDecay.cxx:240
 RooBCPGenDecay.cxx:241
 RooBCPGenDecay.cxx:242
 RooBCPGenDecay.cxx:243
 RooBCPGenDecay.cxx:244
 RooBCPGenDecay.cxx:245
 RooBCPGenDecay.cxx:246
 RooBCPGenDecay.cxx:247
 RooBCPGenDecay.cxx:248
 RooBCPGenDecay.cxx:249
 RooBCPGenDecay.cxx:250
 RooBCPGenDecay.cxx:251
 RooBCPGenDecay.cxx:252
 RooBCPGenDecay.cxx:253
 RooBCPGenDecay.cxx:254
 RooBCPGenDecay.cxx:255
 RooBCPGenDecay.cxx:256
 RooBCPGenDecay.cxx:257