// @(#)root/splot:$Id$
// Author: Muriel Pivk, Anna Kreshuk    10/2005  

/**********************************************************************
 *                                                                    *
 * Copyright (c) 2005 ROOT Foundation,  CERN/PH-SFT                   *
 *                                                                    *
 **********************************************************************/

#include "TSPlot.h"
#include "TVirtualFitter.h"
#include "TH1.h"
#include "TTreePlayer.h"
#include "TTreeFormula.h"
#include "TTreeFormulaManager.h"
#include "TSelectorDraw.h"
#include "TBrowser.h"
#include "TClass.h"
#include "TMath.h"

extern void Yields(Int_t &, Double_t *, Double_t &f, Double_t *x, Int_t iflag);

ClassImp(TSPlot)

//____________________________________________________________________
//Begin_Html <!--
/* -->
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 3.2 Final//EN">
<p>
<b><font size="+2">Overview</font></b>

</p><p>
A common method used in High Energy Physics to perform measurements is
the maximum Likelihood method, exploiting discriminating variables to
disentangle signal from background. The crucial point for such an
analysis to be reliable is to use an exhaustive list of sources of
events combined with an accurate description of all the Probability
Density Functions (PDF).
</p><p>To assess the validity of the fit, a convincing quality check
is to explore further the data sample by examining the distributions of
control variables. A control variable can be obtained for instance by
removing one of the discriminating variables before performing again
the maximum Likelihood fit: this removed variable is a control
variable. The expected distribution of this control variable, for
signal, is to be compared to the one extracted, for signal, from the
data sample. In order to be able to do so, one must be able to unfold
from the distribution of the whole data sample.
</p><p>The TSPlot method allows to reconstruct the distributions for
the control variable, independently for each of the various sources of
events, without making use of any <em>a priori</em> knowledge on <u>this</u>
variable. The aim is thus to use the knowledge available for the
discriminating variables to infer the behaviour of the individual
sources of events with respect to the control variable.
</p><p>
TSPlot is optimal if the control variable is uncorrelated with the discriminating variables.

</p><p>
A detail description of the formalism itself, called <!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
<img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48">, is given in&nbsp;[<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/node1.html#bib:sNIM">1</a>].

</p><p>
<b><font size="+2">The method</font></b> 

</p><p>
The <!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
<img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48"> technique is developped in the above context of a maximum Likelihood method making use of discriminating variables.

</p><p>One considers a data sample in which are merged several species
of events. These species represent various signal components and
background components which all together account for the data sample.
The different terms of the log-Likelihood are:
</p><ul>
<li><img src="gif/sPlot_img6.png" alt="$N$" align="bottom" border="0" height="17" width="22">: the total number of events in the data sample,
</li>
<li><!-- MATH
 ${\rm N}_{\rm s}$
 -->
<img src="gif/sPlot_img7.png" alt="${\rm N}_{\rm s}$" align="middle" border="0" height="34" width="25">: the number of species of events populating the data sample,
</li>
<li><img src="gif/sPlot_img8.png" alt="$N_i$" align="middle" border="0" height="34" width="25">: the number of events expected on the average for the <img src="gif/sPlot_img9.png" alt="$i^{\rm th}$" align="bottom" border="0" height="20" width="25"> species,
</li>
<li><!-- MATH
 ${\rm f}_i(y_e)$
 -->
<img src="gif/sPlot_img10.png" alt="${\rm f}_i(y_e)$" align="middle" border="0" height="37" width="47">: the value of the PDFs of the discriminating variables <img src="gif/sPlot_img11.png" alt="$y$" align="middle" border="0" height="33" width="15"> for the <img src="gif/sPlot_img12.png" alt="$i^{th}$" align="bottom" border="0" height="20" width="25"> species and for event <img src="gif/sPlot_img13.png" alt="$e$" align="bottom" border="0" height="17" width="13">, 
</li>
<li><img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">: the set of control variables which, by definition, do not appear in the expression of the Likelihood function <img src="gif/sPlot_img15.png" alt="${\cal L}$" align="bottom" border="0" height="18" width="18">.
</li>
</ul>
The extended log-Likelihood reads:
<br>
<div align="right">

<!-- MATH
 \begin{equation}
{\cal L}=\sum_{e=1}^{N}\ln \Big\{ \sum_{i=1}^{{\rm N}_{\rm s}}N_i{\rm f}_i(y_e) \Big\} -\sum_{i=1}^{{\rm N}_{\rm s}}N_i ~.
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:eLik"></a><img src="gif/sPlot_img16.png" alt="\begin{displaymath}
{\cal L}=\sum_{e=1}^{N}\ln \Big\{ \sum_{i=1}^{{\rm N}_{\rm s}}N_i{\rm f}_i(y_e) \Big\} -\sum_{i=1}^{{\rm N}_{\rm s}}N_i ~.
\end{displaymath}" border="0" height="59" width="276"></td>
<td align="right" width="10">
(1)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
From this expression, after maximization of <img src="gif/sPlot_img15.png" alt="${\cal L}$" align="bottom" border="0" height="18" width="18"> with respect to the <img src="gif/sPlot_img8.png" alt="$N_i$" align="middle" border="0" height="34" width="25"> parameters, a weight can be computed for every event and each species, in order to obtain later the true distribution <!-- MATH
 ${\hbox{\bf {M}}}_i(x)$
 -->
<img src="gif/sPlot_img17.png" alt="${\hbox{\bf {M}}}_i(x)$" align="middle" border="0" height="37" width="56"> of variable <img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">. If <img src="gif/sPlot_img18.png" alt="${\rm n}$" align="bottom" border="0" height="17" width="15"> is one of the <!-- MATH
 ${\rm N}_{\rm s}$
 -->
<img src="gif/sPlot_img7.png" alt="${\rm N}_{\rm s}$" align="middle" border="0" height="34" width="25"> species present in the data sample, the weight for this species is defined by:
<br>
<div align="right">

<!-- MATH
 \begin{equation}
\begin{Large}\fbox{$
{_s{\cal P}}_{\rm n}(y_e)={\sum_{j=1}^{{\rm N}_{\rm s}} \hbox{\bf V}_{{\rm n}j}{\rm f}_j(y_e)\over\sum_{k=1}^{{\rm N}_{\rm s}}N_k{\rm f}_k(y_e) } $}\end{Large} ~,
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:weightxnotiny"></a><img src="gif/sPlot_img19.png" alt="\begin{displaymath}
\begin{Large}
\fbox{$
{_s{\cal P}}_{\rm n}(y_e)={\sum_{j=1}^...
...um_{k=1}^{{\rm N}_{\rm s}}N_k{\rm f}_k(y_e) } $}\end{Large} ~,
\end{displaymath}" border="0" height="76" width="279"></td>
<td align="right" width="10">
(2)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
where <!-- MATH
 $\hbox{\bf V}_{{\rm n}j}$
 -->
<img src="gif/sPlot_img20.png" alt="$\hbox{\bf V}_{{\rm n}j}$" align="middle" border="0" height="34" width="35">
is the covariance matrix resulting from the Likelihood maximization.
This matrix can be used directly from the fit, but this is numerically
less accurate than the direct computation:
<br>
<div align="right">

<!-- MATH
 \begin{equation}
\hbox{\bf V}^{-1}_{{\rm n}j}~=~
{\partial^2(-{\cal L})\over\partial N_{\rm n}\partial N_j}~=~
\sum_{e=1}^N {{\rm f}_{\rm n}(y_e){\rm f}_j(y_e)\over(\sum_{k=1}^{{\rm N}_{\rm s}}N_k{\rm f}_k(y_e))^2} ~.
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:VarianceMatrixDirect"></a><img src="gif/sPlot_img21.png" alt="\begin{displaymath}
\hbox{\bf V}^{-1}_{{\rm n}j}~=~
{\partial^2(-{\cal L})\over\...
...y_e)\over(\sum_{k=1}^{{\rm N}_{\rm s}}N_k{\rm f}_k(y_e))^2} ~.
\end{displaymath}" border="0" height="58" width="360"></td>
<td align="right" width="10">
(3)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
The distribution of the control variable&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15"> obtained by histogramming the weighted events reproduces, on average, the true distribution <!-- MATH
 ${\hbox{\bf {M}}}_{\rm n}(x)$
 -->
<img src="gif/sPlot_img22.png" alt="${\hbox{\bf {M}}}_{\rm n}(x)$" align="middle" border="0" height="37" width="59">.

<p>
The class TSPlot allows to reconstruct the true distribution <!-- MATH
 ${\hbox{\bf {M}}}_{\rm n}(x)$
 -->
<img src="gif/sPlot_img22.png" alt="${\hbox{\bf {M}}}_{\rm n}(x)$" align="middle" border="0" height="37" width="59"> of a control variable&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15"> for each of the <!-- MATH
 ${\rm N}_{\rm s}$
 -->
<img src="gif/sPlot_img7.png" alt="${\rm N}_{\rm s}$" align="middle" border="0" height="34" width="25"> species from the sole knowledge of the PDFs of the discriminating variables <img src="gif/sPlot_img23.png" alt="${\rm f}_i(y)$" align="middle" border="0" height="37" width="40">. The plots obtained thanks to the TSPlot class are called <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57">.

</p><p>
<b><font size="+2">Some properties and checks</font></b>

</p><p>
Beside reproducing the true distribution, <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> bear remarkable properties:

</p><ul>
<li>
Each <img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">-distribution is properly normalized:
<br>
<div align="right">

<!-- MATH
 \begin{equation}
\sum_{e=1}^{N} {_s{\cal P}}_{\rm n}(y_e)~=~N_{\rm n}~.
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:NormalizationOK"></a><img src="gif/sPlot_img24.png" alt="\begin{displaymath}
\sum_{e=1}^{N} {_s{\cal P}}_{\rm n}(y_e)~=~N_{\rm n}~.
\end{displaymath}" border="0" height="58" width="158"></td>
<td align="right" width="10">
(4)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
</li>
<li>
For any event:
<br>
<div align="right">

<!-- MATH
 \begin{equation}
\sum_{l=1}^{{\rm N}_{\rm s}} {_s{\cal P}}_l(y_e) ~=~1 ~.
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:numberconservation"></a><img src="gif/sPlot_img25.png" alt="\begin{displaymath}
\sum_{l=1}^{{\rm N}_{\rm s}} {_s{\cal P}}_l(y_e) ~=~1 ~.
\end{displaymath}" border="0" height="59" width="140"></td>
<td align="right" width="10">
(5)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
That is to say that, summing up the <!-- MATH
 ${\rm N}_{\rm s}$
 -->
<img src="gif/sPlot_img7.png" alt="${\rm N}_{\rm s}$" align="middle" border="0" height="34" width="25"> <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57">, one recovers the data sample distribution in&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">, and summing up the number of events entering in a <!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
<img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48"> for a given species, one recovers the yield of the species, as provided by the fit. The property&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#eq:NormalizationOK">4</a> is implemented in the TSPlot class as a check.
</li>
<li>the sum of the statistical uncertainties per bin
<br>
<div align="right">

<!-- MATH
 \begin{equation}
\sigma[N_{\rm n}\  _s\tilde{\rm M}_{\rm n}(x) {\delta x}]~=~\sqrt{\sum_{e \subset {\delta x}} ({_s{\cal P}}_{\rm n})^2} ~.
\end{equation}
 -->
<table align="center" width="100%">
<tbody><tr valign="middle"><td align="center" nowrap="nowrap"><a name="eq:ErrorPerBin"></a><img src="gif/sPlot_img26.png" alt="\begin{displaymath}
\sigma[N_{\rm n}\ _s\tilde{\rm M}_{\rm n}(x) {\delta x}]~=~\sqrt{\sum_{e \subset {\delta x}} ({_s{\cal P}}_{\rm n})^2} ~.
\end{displaymath}" border="0" height="55" width="276"></td>
<td align="right" width="10">
(6)</td></tr>
</tbody></table>
<br clear="all"></div><p></p>
reproduces the statistical uncertainty on the yield <img src="gif/sPlot_img27.png" alt="$N_{\rm n}$" align="middle" border="0" height="34" width="28">, as provided by the fit: <!-- MATH
 $\sigma[N_{\rm n}]\equiv\sqrt{\hbox{\bf V}_{{\rm n}{\rm n}}}$
 -->
<img src="gif/sPlot_img28.png" alt="$\sigma[N_{\rm n}]\equiv\sqrt{\hbox{\bf V}_{{\rm n}{\rm n}}}$" align="middle" border="0" height="40" width="123">.
Because of that and since the determination of the yields is optimal
when obtained using a Likelihood fit, one can conclude that the<!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
 <img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48"> technique is itself an optimal method to reconstruct distributions of control variables.
</li>
</ul>

<p>
<b><font size="+2">Different steps followed by TSPlot</font></b>

</p><p>

</p><ol>
<li>A maximum Likelihood fit is performed to obtain the yields <img src="gif/sPlot_img8.png" alt="$N_i$" align="middle" border="0" height="34" width="25"> of the various species. 
The fit relies on discriminating variables&nbsp;<img src="gif/sPlot_img11.png" alt="$y$" align="middle" border="0" height="33" width="15"> uncorrelated with a control variable&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">:
the later is therefore totally absent from the fit. 
</li>
<li>The weights <img src="gif/sPlot_img29.png" alt="${_s{\cal P}}$" align="middle" border="0" height="34" width="27"> are calculated using Eq.&nbsp;(<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#eq:weightxnotiny">2</a>) where the covariance matrix is taken from Minuit.
</li>
<li>Histograms of&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15"> are filled by weighting the events with <img src="gif/sPlot_img29.png" alt="${_s{\cal P}}$" align="middle" border="0" height="34" width="27">. 
</li>
<li>Error bars per bin are given by Eq.&nbsp;(<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#eq:ErrorPerBin">6</a>). 
</li>
</ol>
The <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> reproduce the true distributions of the species in the control variable&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">, within the above defined statistical uncertainties.

<p>
<b><font size="+2">Illustrations</font></b>

</p><p>
To illustrate the technique, one considers an example derived from the analysis where <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57">
have been first used (charmless B decays). One is dealing with a data
sample in which two species are present: the first is termed signal and
the second background. A maximum Likelihood fit is performed to obtain
the two yields <img src="gif/sPlot_img30.png" alt="$N_1$" align="middle" border="0" height="34" width="27"> and <img src="gif/sPlot_img31.png" alt="$N_2$" align="middle" border="0" height="34" width="27">. The fit relies on two discriminating variables collectively denoted&nbsp;<img src="gif/sPlot_img11.png" alt="$y$" align="middle" border="0" height="33" width="15"> which are chosen within three possible variables denoted <img src="gif/sPlot_img1.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39">, <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35"> and <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20">.
The variable which is not incorporated in&nbsp;<img src="gif/sPlot_img11.png" alt="$y$" align="middle" border="0" height="33" width="15"> is used as the control variable&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15">. The six distributions of the three variables are assumed to be the ones depicted in Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:pdfs">1</a>.

</p><p>

</p><div align="center"><a name="fig:pdfs"></a><a name="106"></a>
<table>
<caption align="bottom"><strong>Figure 1:</strong>
Distributions of the three discriminating variables available to perform the Likelihood fit: 
<img src="gif/sPlot_img32.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39">, <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35">, <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20">.
Among the three variables, two are used to perform the fit while one is
kept out of the fit to serve the purpose of a control variable. The
three distributions on the top (resp. bottom) of the figure correspond
to the signal (resp. background). The unit of the vertical axis is
chosen such that it indicates the number of entries per bin, if one
slices the histograms in 25 bins.</caption>
<tbody><tr><td><img src="gif/sPlot_img33.png" alt="\begin{figure}\begin{center}
\mbox{{\psfig{file=pdfmesNIM.eps,width=0.33\linewi...
...th}}
{\psfig{file=pdffiNIM.eps,width=0.33\linewidth}}}
\end{center}\end{figure}" border="0" height="162" width="544"></td></tr>
</tbody></table>
</div>

<p>
A data sample being built through a Monte Carlo simulation based on the distributions shown in Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:pdfs">1</a>, one obtains the three distributions of Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:pdfstot">2</a>. Whereas the distribution of&nbsp;<img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35"> clearly indicates the presence of the signal, the distribution of <img src="gif/sPlot_img1.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39"> and <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20"> are less obviously populated by signal.

</p><p>

</p><div align="center"><a name="fig:pdfstot"></a><a name="169"></a>
<table>
<caption align="bottom"><strong>Figure 2:</strong>
Distributions of the three discriminating variables for signal plus
background. The three distributions are the ones obtained from a data
sample obtained through a Monte Carlo simulation based on the
distributions shown in Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:pdfs">1</a>.  The data sample consists of 500 signal events and 5000 background events.</caption>
<tbody><tr><td><img src="gif/sPlot_img34.png" alt="\begin{figure}\begin{center}
\mbox{{\psfig{file=genmesTOTNIM.eps,width=0.33\lin...
...}
{\psfig{file=genfiTOTNIM.eps,width=0.33\linewidth}}}
\end{center}\end{figure}" border="0" height="160" width="545"></td></tr>
</tbody></table>
</div>

<p>
Chosing <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35"> and <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20"> as discriminating variables to determine <img src="gif/sPlot_img30.png" alt="$N_1$" align="middle" border="0" height="34" width="27"> and <img src="gif/sPlot_img31.png" alt="$N_2$" align="middle" border="0" height="34" width="27"> through a maximum Likelihood fit, one builds, for the control variable <img src="gif/sPlot_img1.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39"> which is unknown to the fit, the two <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> for signal and background shown in Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:messPlots">3</a>. One observes that the <!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
<img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48">
for signal reproduces correctly the PDF even where the latter vanishes,
although the error bars remain sizeable. This results from the almost
complete cancellation between positive and negative weights: the sum of
weights is close to zero while the sum of weights squared is not. The
occurence of negative weights occurs through the appearance of the
covariance matrix, and its negative components, in the definition of
Eq.&nbsp;(<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#eq:weightxnotiny">2</a>).

</p><p>
A word of caution is in order with respect to the error bars. Whereas
their sum in quadrature is identical to the statistical uncertainties
of the yields determined by the fit, and if, in addition, they are
asymptotically correct, the error bars should be handled with care for
low statistics and/or for too fine binning. This is because the error
bars do not incorporate two known properties of the PDFs: PDFs are
positive definite and can be non-zero in a given x-bin, even if in the
particular data sample at hand, no event is observed in this bin. The
latter limitation is not specific to<!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
 <img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57">,
rather it is always present when one is willing to infer the PDF at the
origin of an histogram, when, for some bins, the number of entries does
not guaranty the applicability of the Gaussian regime. In such
situations, a satisfactory practice is to attach allowed ranges to the
histogram to indicate the upper and lower limits of the PDF value which
are consistent with the actual observation, at a given confidence
level.
</p><p>

</p><div align="center"><a name="fig:messPlots"></a><a name="127"></a>
<table>
<caption align="bottom"><strong>Figure 3:</strong>
The <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> (signal on the left, background on the right) obtained for <img src="gif/sPlot_img32.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39"> are represented as dots with error bars. They are obtained from a fit using only information from <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35"> and <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20">.</caption>
<tbody><tr><td><img src="gif/sPlot_img35.png" alt="\begin{figure}\begin{center}
\mbox{\psfig{file=mass-sig-sPlot.eps,width=0.48\li...
... \psfig{file=mass-bkg-sPlot.eps,width=0.48\linewidth}}
\end{center}\end{figure}" border="0" height="181" width="539"></td></tr>
</tbody></table>
</div>

<p>
Chosing <img src="gif/sPlot_img1.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39"> and <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35"> as discriminating variables to determine <img src="gif/sPlot_img30.png" alt="$N_1$" align="middle" border="0" height="34" width="27"> and <img src="gif/sPlot_img31.png" alt="$N_2$" align="middle" border="0" height="34" width="27"> through a maximum Likelihood fit, one builds, for the control variable <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20"> which is unknown to the fit, the two <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> for signal and background shown in Fig.&nbsp;<a href="http://www.slac.stanford.edu/%7Epivk/sPlot/sPlot_ROOT/sPlot_ROOT.html#fig:FisPlots">4</a>. In the <!-- MATH
 $\hbox{$_s$}{\cal P}lot$
 -->
<img src="gif/sPlot_img5.png" alt="$\hbox{$_s$}{\cal P}lot$" align="middle" border="0" height="34" width="48"> for signal one observes that error bars are the largest in the&nbsp;<img src="gif/sPlot_img14.png" alt="$x$" align="bottom" border="0" height="17" width="15"> regions where the background is the largest.

</p><p>

</p><div align="center"><a name="fig:FisPlots"></a><a name="136"></a>
<table>
<caption align="bottom"><strong>Figure 4:</strong>
The <!-- MATH
 $\hbox{$_s$}{\cal P}lots$
 -->
<img src="gif/sPlot_img4.png" alt="$\hbox {$_s$}{\cal P}lots$" align="middle" border="0" height="34" width="57"> (signal on the left, background on the right) obtained for <img src="gif/sPlot_img3.png" alt="${\cal F}$" align="bottom" border="0" height="18" width="20"> are represented as dots with error bars. They are obtained from a fit using only information from <img src="gif/sPlot_img32.png" alt="${m_{\rm ES}}$" align="middle" border="0" height="33" width="39"> and <img src="gif/sPlot_img2.png" alt="$\Delta E$" align="bottom" border="0" height="17" width="35">.</caption>
<tbody><tr><td><img src="gif/sPlot_img36.png" alt="\begin{figure}\begin{center}
\mbox{\psfig{file=fisher-sig-sPlot.eps,width=0.48\...
...psfig{file=fisher-bkg-sPlot.eps,width=0.48\linewidth}}
\end{center}\end{figure}" border="0" height="180" width="539"></td></tr>
</tbody></table>
</div>

<p>
The results above can be obtained by running the tutorial TestSPlot.C
</p>
<!--*/
//-->End_Html


//____________________________________________________________________
TSPlot::TSPlot() : 
 fTree(0),
 fTreename(0),
 fVarexp(0),
 fSelection(0)
{
   // default constructor (used by I/O only)
   fNx = 0;
   fNy=0;
   fNevents = 0;
   fNSpecies=0;
   fNumbersOfEvents=0;
}

//____________________________________________________________________
TSPlot::TSPlot(Int_t nx, Int_t ny, Int_t ne, Int_t ns, TTree *tree) : 
 fTreename(0),
 fVarexp(0),
 fSelection(0)

{
   //normal TSPlot constructor
   // nx :  number of control variables
   // ny :  number of discriminating variables
   // ne :  total number of events
   // ns :  number of species
   // tree: input data
   
   fNx = nx;
   fNy=ny;
   fNevents = ne;
   fNSpecies=ns;

   fXvar.ResizeTo(fNevents, fNx);
   fYvar.ResizeTo(fNevents, fNy);
   fYpdf.ResizeTo(fNevents, fNSpecies*fNy);
   fSWeights.ResizeTo(fNevents, fNSpecies*(fNy+1));
   fTree = tree;
   fNumbersOfEvents = 0;
}

//____________________________________________________________________
TSPlot::~TSPlot()
{
   // destructor
   
   if (fNumbersOfEvents)
      delete [] fNumbersOfEvents;
   if (!fXvarHists.IsEmpty())
      fXvarHists.Delete();
   if (!fYvarHists.IsEmpty())
      fYvarHists.Delete();
   if (!fYpdfHists.IsEmpty())
      fYpdfHists.Delete();
}

//____________________________________________________________________
void TSPlot::Browse(TBrowser *b)
{
   //To browse the histograms

   if (!fSWeightsHists.IsEmpty()) {
      TIter next(&fSWeightsHists);
      TH1D* h = 0;
      while ((h = (TH1D*)next()))
         b->Add(h,h->GetName());
   }

   if (!fYpdfHists.IsEmpty()) {
      TIter next(&fYpdfHists);
      TH1D* h = 0;
      while ((h = (TH1D*)next()))
         b->Add(h,h->GetName());
   }
   if (!fYvarHists.IsEmpty()) {
      TIter next(&fYvarHists);
      TH1D* h = 0;
      while ((h = (TH1D*)next()))
         b->Add(h,h->GetName());
   }
   if (!fXvarHists.IsEmpty()) {
      TIter next(&fXvarHists);
      TH1D* h = 0;
      while ((h = (TH1D*)next()))
         b->Add(h,h->GetName());
   }
   b->Add(&fSWeights, "sWeights");
}


//____________________________________________________________________
void TSPlot::SetInitialNumbersOfSpecies(Int_t *numbers)
{
//Set the initial number of events of each species - used
//as initial estimates in minuit

   if (!fNumbersOfEvents)
      fNumbersOfEvents = new Double_t[fNSpecies];
   for (Int_t i=0; i<fNSpecies; i++)
      fNumbersOfEvents[i]=numbers[i];
}

//____________________________________________________________________
void TSPlot::MakeSPlot(Option_t *option)
{
//Calculates the sWeights
//The option controls the print level
//"Q" - no print out
//"V" - prints the estimated #of events in species - default
//"VV" - as "V" + the minuit printing + sums of weights for control


   if (!fNumbersOfEvents){
      Error("MakeSPlot","Initial numbers of events in species have not been set");
      return;
   }
   Int_t i, j, ispecies;

   TString opt = option;
   opt.ToUpper();
   opt.ReplaceAll("VV", "W");

   //make sure that global fitter is minuit
   char s[]="TFitter";
   if (TVirtualFitter::GetFitter()){
      Int_t strdiff=strcmp(TVirtualFitter::GetFitter()->IsA()->GetName(), s);
      if (strdiff!=0)
         delete TVirtualFitter::GetFitter();
   }


   TVirtualFitter *minuit = TVirtualFitter::Fitter(0, 2);
   fPdfTot.ResizeTo(fNevents, fNSpecies);

   //now let's do it, excluding different yvars
   //for iplot = -1 none is excluded
   for (Int_t iplot=-1; iplot<fNy; iplot++){
      for (i=0; i<fNevents; i++){
         for (ispecies=0; ispecies<fNSpecies; ispecies++){
            fPdfTot(i, ispecies)=1;
            for (j=0; j<fNy; j++){
               if (j!=iplot)
                  fPdfTot(i, ispecies)*=fYpdf(i, ispecies*fNy+j);
            }
         }
      }
      minuit->Clear();
      minuit->SetFCN(Yields);
      Double_t arglist[10];
      //set the print level
      if (opt.Contains("Q")||opt.Contains("V")){
         arglist[0]=-1;
      }
      if (opt.Contains("W"))
         arglist[0]=0;
      minuit->ExecuteCommand("SET PRINT", arglist, 1);

      minuit->SetObjectFit(&fPdfTot); //a tricky way to get fPdfTot matrix to fcn
      for (ispecies=0; ispecies<fNSpecies; ispecies++)
         minuit->SetParameter(ispecies, "", fNumbersOfEvents[ispecies], 1, 0, 0);

      minuit->ExecuteCommand("MIGRAD", arglist, 0);
      for (ispecies=0; ispecies<fNSpecies; ispecies++){
         fNumbersOfEvents[ispecies]=minuit->GetParameter(ispecies);
         if (!opt.Contains("Q"))
            printf("estimated #of events in species %d = %f\n", ispecies, fNumbersOfEvents[ispecies]);
      }
      if (!opt.Contains("Q"))
         printf("\n");
      Double_t *covmat = minuit->GetCovarianceMatrix();
      SPlots(covmat, iplot);

      if (opt.Contains("W")){
         Double_t *sumweight = new Double_t[fNSpecies];
         for (i=0; i<fNSpecies; i++){
            sumweight[i]=0;
            for (j=0; j<fNevents; j++)
               sumweight[i]+=fSWeights(j, (iplot+1)*fNSpecies + i);
            printf("checking sum of weights[%d]=%f\n", i, sumweight[i]);
         }
         printf("\n");
         delete [] sumweight;
      }
   }
}

//____________________________________________________________________
void TSPlot::SPlots(Double_t *covmat, Int_t i_excl)
{
//Computes the sWeights from the covariance matrix

   Int_t i, ispecies, k;
   Double_t numerator, denominator;
   for (i=0; i<fNevents; i++){
      denominator=0;
      for (ispecies=0; ispecies<fNSpecies; ispecies++)
         denominator+=fNumbersOfEvents[ispecies]*fPdfTot(i, ispecies);
      for (ispecies=0; ispecies<fNSpecies; ispecies++){
         numerator=0;
         for (k=0; k<fNSpecies; k++)
            numerator+=covmat[ispecies*fNSpecies+k]*fPdfTot(i, k);
         fSWeights(i, (i_excl+1)*fNSpecies + ispecies)=numerator/denominator;
      }
   }

}

//____________________________________________________________________
void TSPlot::GetSWeights(TMatrixD &weights)
{
//Returns the matrix of sweights

   if (weights.GetNcols()!=fNSpecies*(fNy+1) || weights.GetNrows()!=fNevents)
      weights.ResizeTo(fNevents, fNSpecies*(fNy+1));
   weights = fSWeights;
}

//____________________________________________________________________
void TSPlot::GetSWeights(Double_t *weights)
{
//Returns the matrix of sweights. It is assumed that the array passed in the argurment is big enough

   for (Int_t i=0; i<fNevents; i++){
      for (Int_t j=0; j<fNSpecies; j++){
         weights[i*fNSpecies+j]=fSWeights(i, j);
      }
   }
}

//____________________________________________________________________
void TSPlot::FillXvarHists(Int_t nbins)
{
//Fills the histograms of x variables (not weighted) with nbins

   Int_t i, j;

   if (!fXvarHists.IsEmpty()){
      if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
         fXvarHists.Delete();
      else
         return;
   }

   //make the histograms
   char name[10];
   for (i=0; i<fNx; i++){
      snprintf(name,10, "x%d", i);
      TH1D *h = new TH1D(name, name, nbins, fMinmax(0, i), fMinmax(1, i));
      for (j=0; j<fNevents; j++)
         h->Fill(fXvar(j, i));
      fXvarHists.Add(h);
   }
                    
}

//____________________________________________________________________
TObjArray* TSPlot::GetXvarHists()
{
//Returns the array of histograms of x variables (not weighted)
//If histograms have not already
//been filled, they are filled with default binning 100.

   Int_t nbins = 100;
   if (fXvarHists.IsEmpty())
      FillXvarHists(nbins);
   else if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
      FillXvarHists(nbins);
   return &fXvarHists;
}

//____________________________________________________________________
TH1D *TSPlot::GetXvarHist(Int_t ixvar)
{
//Returns the histogram of variable #ixvar
//If histograms have not already
//been filled, they are filled with default binning 100.

   Int_t nbins = 100;
   if (fXvarHists.IsEmpty())
      FillXvarHists(nbins);
   else if (((TH1D*)fXvarHists.First())->GetNbinsX()!=nbins)
      FillXvarHists(nbins);

   return (TH1D*)fXvarHists.UncheckedAt(ixvar);
}

//____________________________________________________________________
void TSPlot::FillYvarHists(Int_t nbins)
{
//Fill the histograms of y variables

   Int_t i, j;

   if (!fYvarHists.IsEmpty()){
      if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
         fYvarHists.Delete();
      else
         return;
   }

   //make the histograms
   char name[10];
   for (i=0; i<fNy; i++){
      snprintf(name,10, "y%d", i);
      TH1D *h=new TH1D(name, name, nbins, fMinmax(0, fNx+i), fMinmax(1, fNx+i));
      for (j=0; j<fNevents; j++)
         h->Fill(fYvar(j, i));
      fYvarHists.Add(h);
   }
}

//____________________________________________________________________
TObjArray* TSPlot::GetYvarHists()
{
//Returns the array of histograms of y variables. If histograms have not already
//been filled, they are filled with default binning 100.

   Int_t nbins = 100;
   if (fYvarHists.IsEmpty())
      FillYvarHists(nbins);
   else if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
      FillYvarHists(nbins);
   return &fYvarHists;
}

//____________________________________________________________________
TH1D *TSPlot::GetYvarHist(Int_t iyvar)
{
//Returns the histogram of variable iyvar.If histograms have not already
//been filled, they are filled with default binning 100.
 
   Int_t nbins = 100;
   if (fYvarHists.IsEmpty())
      FillYvarHists(nbins);
   else if (((TH1D*)fYvarHists.First())->GetNbinsX()!=nbins)
      FillYvarHists(nbins);
   return (TH1D*)fYvarHists.UncheckedAt(iyvar);
}

//____________________________________________________________________
void TSPlot::FillYpdfHists(Int_t nbins)
{
//Fills the histograms of pdf-s of y variables with binning nbins

   Int_t i, j, ispecies;

   if (!fYpdfHists.IsEmpty()){
      if (((TH1D*)fYpdfHists.First())->GetNbinsX()!=nbins)
         fYpdfHists.Delete();
      else
         return;
   }

   char name[30];
   for (ispecies=0; ispecies<fNSpecies; ispecies++){
      for (i=0; i<fNy; i++){
         snprintf(name,30, "pdf_species%d_y%d", ispecies, i);
         //TH1D *h = new TH1D(name, name, nbins, ypdfmin[ispecies*fNy+i], ypdfmax[ispecies*fNy+i]);
         TH1D *h = new TH1D(name, name, nbins, fMinmax(0, fNx+fNy+ispecies*fNy+i), fMinmax(1, fNx+fNy+ispecies*fNy+i)); 
         for (j=0; j<fNevents; j++)
            h->Fill(fYpdf(j, ispecies*fNy+i));
         fYpdfHists.Add(h);
      }
   }
}

//____________________________________________________________________
TObjArray* TSPlot::GetYpdfHists()
{
//Returns the array of histograms of pdf's of y variables with binning nbins
//If histograms have not already
//been filled, they are filled with default binning 100.

   Int_t nbins = 100;
   if (fYpdfHists.IsEmpty())
      FillYpdfHists(nbins);

   return &fYpdfHists;
}

//____________________________________________________________________
TH1D *TSPlot::GetYpdfHist(Int_t iyvar, Int_t ispecies)
{
//Returns the histogram of the pdf of variable #iyvar for species #ispecies, binning nbins
//If histograms have not already
//been filled, they are filled with default binning 100.

   Int_t nbins = 100;
   if (fYpdfHists.IsEmpty())
      FillYpdfHists(nbins);

   return (TH1D*)fYpdfHists.UncheckedAt(fNy*ispecies+iyvar);
}

//____________________________________________________________________
void TSPlot::FillSWeightsHists(Int_t nbins)
{
   //The order of histograms in the array:
   //x0_species0, x0_species1,..., x1_species0, x1_species1,..., y0_species0, y0_species1,...
   //If the histograms have already been filled with a different binning, they are refilled
   //and all histograms are deleted

   if (fSWeights.GetNoElements()==0){
      Error("GetSWeightsHists", "SWeights were not computed");
      return;
   }

   if (!fSWeightsHists.IsEmpty()){
      if (((TH1D*)fSWeightsHists.First())->GetNbinsX()!=nbins)
         fSWeightsHists.Delete();
      else
         return;
   }

   char name[30];

   //Fill histograms of x-variables weighted with sWeights
   for (Int_t ivar=0; ivar<fNx; ivar++){
      for (Int_t ispecies=0; ispecies<fNSpecies; ispecies++){
            snprintf(name,30, "x%d_species%d", ivar, ispecies);
            TH1D *h = new TH1D(name, name, nbins, fMinmax(0, ivar), fMinmax(1, ivar));
            h->Sumw2();
            for (Int_t ievent=0; ievent<fNevents; ievent++)
               h->Fill(fXvar(ievent, ivar), fSWeights(ievent, ispecies));
            fSWeightsHists.AddLast(h); 
         }
   }

   //Fill histograms of y-variables (exluded from the fit), weighted with sWeights
   for (Int_t iexcl=0; iexcl<fNy; iexcl++){
      for(Int_t ispecies=0; ispecies<fNSpecies; ispecies++){
            snprintf(name,30, "y%d_species%d", iexcl, ispecies);
            TH1D *h = new TH1D(name, name, nbins, fMinmax(0, fNx+iexcl), fMinmax(1, fNx+iexcl));
            h->Sumw2();
            for (Int_t ievent=0; ievent<fNevents; ievent++)
               h->Fill(fYvar(ievent, iexcl), fSWeights(ievent, fNSpecies*(iexcl+1)+ispecies));
            fSWeightsHists.AddLast(h);
      }
   }
}

//____________________________________________________________________
TObjArray *TSPlot::GetSWeightsHists()
{
   //Returns an array of all histograms of variables, weighted with sWeights
   //If histograms have not been already filled, they are filled with default binning 50
   //The order of histograms in the array:
   //x0_species0, x0_species1,..., x1_species0, x1_species1,..., y0_species0, y0_species1,...

   Int_t nbins = 50; //default binning
   if (fSWeightsHists.IsEmpty())
      FillSWeightsHists(nbins);

   return &fSWeightsHists;
}

//____________________________________________________________________
void TSPlot::RefillHist(Int_t type, Int_t nvar, Int_t nbins, Double_t min, Double_t max, Int_t nspecies)
{
   //The Fill...Hist() methods fill the histograms with the real limits on the variables
   //This method allows to refill the specified histogram with user-set boundaries min and max
   //Parameters:
   //type = 1 - histogram of x variable #nvar
   //     = 2 - histogram of y variable #nvar
   //     = 3 - histogram of y_pdf for y #nvar and species #nspecies
   //     = 4 - histogram of x variable #nvar, species #nspecies, WITH sWeights
   //     = 5 - histogram of y variable #nvar, species #nspecies, WITH sWeights

   if (type<1 || type>5){
      Error("RefillHist", "type must lie between 1 and 5");
      return;
   }
   char name[20];
   Int_t j;
   TH1D *hremove;
   if (type==1){
      hremove = (TH1D*)fXvarHists.RemoveAt(nvar);
      delete hremove;
      snprintf(name,20,"x%d",nvar);
      TH1D *h = new TH1D(name, name, nbins, min, max);
      for (j=0; j<fNevents;j++)
         h->Fill(fXvar(j, nvar));
      fXvarHists.AddAt(h, nvar);
   }
   if (type==2){
      hremove = (TH1D*)fYvarHists.RemoveAt(nvar);
      delete hremove;
      snprintf(name,20, "y%d", nvar);
      TH1D *h = new TH1D(name, name, nbins, min, max);
      for (j=0; j<fNevents;j++)
         h->Fill(fYvar(j, nvar));
      fXvarHists.AddAt(h, nvar); 
   }
   if (type==3){
      hremove = (TH1D*)fYpdfHists.RemoveAt(nspecies*fNy+nvar);
      delete hremove;
      snprintf(name,20, "pdf_species%d_y%d", nspecies, nvar);
      TH1D *h=new TH1D(name, name, nbins, min, max);
      for (j=0; j<fNevents; j++)
         h->Fill(fYpdf(j, nspecies*fNy+nvar));
      fYpdfHists.AddAt(h, nspecies*fNy+nvar);
   }
   if (type==4){
      hremove = (TH1D*)fSWeightsHists.RemoveAt(fNSpecies*nvar+nspecies);
      delete hremove;
      snprintf(name,20, "x%d_species%d", nvar, nspecies);
      TH1D *h = new TH1D(name, name, nbins, min, max);
      h->Sumw2();
      for (Int_t ievent=0; ievent<fNevents; ievent++)
         h->Fill(fXvar(ievent, nvar), fSWeights(ievent, nspecies));
      fSWeightsHists.AddAt(h, fNSpecies*nvar+nspecies);
   }
   if (type==5){
      hremove = (TH1D*)fSWeightsHists.RemoveAt(fNx*fNSpecies + fNSpecies*nvar+nspecies);
      delete hremove;
      snprintf(name,20, "y%d_species%d", nvar, nspecies);
      TH1D *h = new TH1D(name, name, nbins, min, max);
      h->Sumw2();
      for (Int_t ievent=0; ievent<fNevents; ievent++)
         h->Fill(fYvar(ievent, nvar), fSWeights(ievent, nspecies));
      fSWeightsHists.AddAt(h, fNx*fNSpecies + fNSpecies*nvar+nspecies); 
   }  
}
//____________________________________________________________________
TH1D *TSPlot::GetSWeightsHist(Int_t ixvar, Int_t ispecies,Int_t iyexcl)
{
   //Returns the histogram of a variable, weithed with sWeights
   //If histograms have not been already filled, they are filled with default binning 50
   //If parameter ixvar!=-1, the histogram of x-variable #ixvar is returned for species ispecies
   //If parameter ixvar==-1, the histogram of y-variable #iyexcl is returned for species ispecies
   //If the histogram has already been filled and the binning is different from the parameter nbins
   //all histograms with old binning will be deleted and refilled.


   Int_t nbins = 50; //default binning
   if (fSWeightsHists.IsEmpty())
      FillSWeightsHists(nbins);

   if (ixvar==-1)
      return (TH1D*)fSWeightsHists.UncheckedAt(fNx*fNSpecies + fNSpecies*iyexcl+ispecies);
   else
      return (TH1D*)fSWeightsHists.UncheckedAt(fNSpecies*ixvar + ispecies);

}


//____________________________________________________________________
void TSPlot::SetTree(TTree *tree)
{
   // Set the input Tree
   fTree = tree;
}

//____________________________________________________________________
void TSPlot::SetTreeSelection(const char* varexp, const char *selection, Long64_t firstentry)
{
   //Specifies the variables from the tree to be used for splot
   //
   //Variables fNx, fNy, fNSpecies and fNEvents should already be set!
   //
   //In the 1st parameter it is assumed that first fNx variables are x(control variables), 
   //then fNy y variables (discriminating variables), 
   //then fNy*fNSpecies ypdf variables (probability distribution functions of dicriminating 
   //variables for different species). The order of pdfs should be: species0_y0, species0_y1,...
   //species1_y0, species1_y1,...species[fNSpecies-1]_y0...
   //The 2nd parameter allows to make a cut
   //TTree::Draw method description contains more details on specifying expression and selection 

   TTreeFormula **var;
   std::vector<TString> cnames;
   TList *formulaList = new TList();
   TSelectorDraw *selector = (TSelectorDraw*)(((TTreePlayer*)fTree->GetPlayer())->GetSelector());

   Long64_t entry, entryNumber;
   Int_t i,nch;
   Int_t ncols;
   TObjArray *leaves = fTree->GetListOfLeaves();

   fTreename= new TString(fTree->GetName());
   if (varexp)
      fVarexp = new TString(varexp);
   if (selection)
      fSelection= new TString(selection);

   nch = varexp ? strlen(varexp) : 0;


//*-*- Compile selection expression if there is one
   TTreeFormula *select = 0;
   if (selection && strlen(selection)) {
      select = new TTreeFormula("Selection",selection,fTree);
      if (!select) return;
      if (!select->GetNdim()) { delete select; return; }
      formulaList->Add(select);
   }
//*-*- if varexp is empty, take first nx + ny + ny*nspecies columns by default

   if (nch == 0) {
      ncols = fNx + fNy + fNy*fNSpecies;
      for (i=0;i<ncols;i++) {
         cnames.push_back( leaves->At(i)->GetName() );
      }
//*-*- otherwise select only the specified columns
   } else {
      ncols = selector->SplitNames(varexp,cnames);
   }
   var = new TTreeFormula* [ncols];
   Double_t *xvars = new Double_t[ncols];

   fMinmax.ResizeTo(2, ncols);
   for (i=0; i<ncols; i++){
      fMinmax(0, i)=1e30; 
      fMinmax(1, i)=-1e30;
   }

//*-*- Create the TreeFormula objects corresponding to each column
   for (i=0;i<ncols;i++) {
      var[i] = new TTreeFormula("Var1",cnames[i].Data(),fTree);
      formulaList->Add(var[i]);
   }

//*-*- Create a TreeFormulaManager to coordinate the formulas
   TTreeFormulaManager *manager=0;
   if (formulaList->LastIndex()>=0) {
      manager = new TTreeFormulaManager;
      for(i=0;i<=formulaList->LastIndex();i++) {
         manager->Add((TTreeFormula*)formulaList->At(i));
      }
      manager->Sync();
   }
//*-*- loop on all selected entries
   // fSelectedRows = 0;
   Int_t tnumber = -1;
   Long64_t selectedrows=0;
   for (entry=firstentry;entry<firstentry+fNevents;entry++) {
      entryNumber = fTree->GetEntryNumber(entry);
      if (entryNumber < 0) break;
      Long64_t localEntry = fTree->LoadTree(entryNumber);
      if (localEntry < 0) break;
      if (tnumber != fTree->GetTreeNumber()) {
         tnumber = fTree->GetTreeNumber();
         if (manager) manager->UpdateFormulaLeaves();
      }
      Int_t ndata = 1;
      if (manager && manager->GetMultiplicity()) {
         ndata = manager->GetNdata();
      }

      for(Int_t inst=0;inst<ndata;inst++) {
         Bool_t loaded = kFALSE;
         if (select) {
            if (select->EvalInstance(inst) == 0) {
               continue;
            }
         }

         if (inst==0) loaded = kTRUE;
         else if (!loaded) {
            // EvalInstance(0) always needs to be called so that
            // the proper branches are loaded.
            for (i=0;i<ncols;i++) {
               var[i]->EvalInstance(0);
            }
            loaded = kTRUE;
         }

         for (i=0;i<ncols;i++) {
            xvars[i] = var[i]->EvalInstance(inst);
         }
         
         // curentry = entry-firstentry;
         //printf("event#%d\n", curentry);
         //for (i=0; i<ncols; i++)
          //  printf("xvars[%d]=%f\n", i, xvars[i]);
         //selectedrows++;
         for (i=0; i<fNx; i++){
            fXvar(selectedrows, i) = xvars[i];
            if (fXvar(selectedrows, i) < fMinmax(0, i))
               fMinmax(0, i)=fXvar(selectedrows, i);
            if (fXvar(selectedrows, i) > fMinmax(1, i))
               fMinmax(1, i)=fXvar(selectedrows, i);
         }
         for (i=0; i<fNy; i++){
            fYvar(selectedrows, i) = xvars[i+fNx];
            //printf("y_in_loop(%d, %d)=%f, xvars[%d]=%f\n", selectedrows, i, fYvar(selectedrows, i), i+fNx, xvars[i+fNx]);
            if (fYvar(selectedrows, i) < fMinmax(0, i+fNx))
               fMinmax(0, i+fNx) = fYvar(selectedrows, i);
            if (fYvar(selectedrows, i) > fMinmax(1, i+fNx))
               fMinmax(1, i+fNx) = fYvar(selectedrows, i);
            for (Int_t j=0; j<fNSpecies; j++){
               fYpdf(selectedrows, j*fNy + i)=xvars[j*fNy + i+fNx+fNy];
               if (fYpdf(selectedrows, j*fNy+i) < fMinmax(0, j*fNy+i+fNx+fNy))
                  fMinmax(0, j*fNy+i+fNx+fNy) = fYpdf(selectedrows, j*fNy+i);
               if (fYpdf(selectedrows, j*fNy+i) > fMinmax(1, j*fNy+i+fNx+fNy))
                  fMinmax(1, j*fNy+i+fNx+fNy) = fYpdf(selectedrows, j*fNy+i);
            }
         }
      selectedrows++;   
      }
   }
   fNevents=selectedrows;
  // for (i=0; i<fNevents; i++){
    //  printf("event#%d\n", i);
      //for (Int_t iy=0; iy<fNy; iy++)
        // printf("y[%d]=%f\n", iy, fYvar(i, iy));
      //for (Int_t ispecies=0; ispecies<fNSpecies; ispecies++){
      //   for (Int_t iy=0; iy<fNy; iy++)
        //    printf("ypdf[sp. %d, y %d]=%f\n", ispecies, iy, fYpdf(i, ispecies*fNy+iy));
     // }
   //}
   delete [] xvars;
   delete [] var;
}

//____________________________________________________________________
void Yields(Int_t &, Double_t *, Double_t &f, Double_t *x, Int_t /*iflag*/)
{
// FCN-function for Minuit

   Double_t lik;
   Int_t i, ispecies;

   TVirtualFitter *fitter = TVirtualFitter::GetFitter();
   TMatrixD *pdftot = (TMatrixD*)fitter->GetObjectFit();
   Int_t nev = pdftot->GetNrows();
   Int_t nes = pdftot->GetNcols();
   f=0;
   for (i=0; i<nev; i++){
      lik=0;
      for (ispecies=0; ispecies<nes; ispecies++)
         lik+=x[ispecies]*(*pdftot)(i, ispecies);
      if (lik<0) lik=1;
      f+=TMath::Log(lik);
   }
   //extended likelihood, equivalent to chi2
   Double_t ntot=0;
   for (i=0; i<nes; i++)
      ntot += x[i];
   f = -2*(f-ntot);
}

 TSPlot.cxx:1
 TSPlot.cxx:2
 TSPlot.cxx:3
 TSPlot.cxx:4
 TSPlot.cxx:5
 TSPlot.cxx:6
 TSPlot.cxx:7
 TSPlot.cxx:8
 TSPlot.cxx:9
 TSPlot.cxx:10
 TSPlot.cxx:11
 TSPlot.cxx:12
 TSPlot.cxx:13
 TSPlot.cxx:14
 TSPlot.cxx:15
 TSPlot.cxx:16
 TSPlot.cxx:17
 TSPlot.cxx:18
 TSPlot.cxx:19
 TSPlot.cxx:20
 TSPlot.cxx:21
 TSPlot.cxx:22
 TSPlot.cxx:23
 TSPlot.cxx:24
 TSPlot.cxx:25
 TSPlot.cxx:26
 TSPlot.cxx:27
 TSPlot.cxx:28
 TSPlot.cxx:29
 TSPlot.cxx:30
 TSPlot.cxx:31
 TSPlot.cxx:32
 TSPlot.cxx:33
 TSPlot.cxx:34
 TSPlot.cxx:35
 TSPlot.cxx:36
 TSPlot.cxx:37
 TSPlot.cxx:38
 TSPlot.cxx:39
 TSPlot.cxx:40
 TSPlot.cxx:41
 TSPlot.cxx:42
 TSPlot.cxx:43
 TSPlot.cxx:44
 TSPlot.cxx:45
 TSPlot.cxx:46
 TSPlot.cxx:47
 TSPlot.cxx:48
 TSPlot.cxx:49
 TSPlot.cxx:50
 TSPlot.cxx:51
 TSPlot.cxx:52
 TSPlot.cxx:53
 TSPlot.cxx:54
 TSPlot.cxx:55
 TSPlot.cxx:56
 TSPlot.cxx:57
 TSPlot.cxx:58
 TSPlot.cxx:59
 TSPlot.cxx:60
 TSPlot.cxx:61
 TSPlot.cxx:62
 TSPlot.cxx:63
 TSPlot.cxx:64
 TSPlot.cxx:65
 TSPlot.cxx:66
 TSPlot.cxx:67
 TSPlot.cxx:68
 TSPlot.cxx:69
 TSPlot.cxx:70
 TSPlot.cxx:71
 TSPlot.cxx:72
 TSPlot.cxx:73
 TSPlot.cxx:74
 TSPlot.cxx:75
 TSPlot.cxx:76
 TSPlot.cxx:77
 TSPlot.cxx:78
 TSPlot.cxx:79
 TSPlot.cxx:80
 TSPlot.cxx:81
 TSPlot.cxx:82
 TSPlot.cxx:83
 TSPlot.cxx:84
 TSPlot.cxx:85
 TSPlot.cxx:86
 TSPlot.cxx:87
 TSPlot.cxx:88
 TSPlot.cxx:89
 TSPlot.cxx:90
 TSPlot.cxx:91
 TSPlot.cxx:92
 TSPlot.cxx:93
 TSPlot.cxx:94
 TSPlot.cxx:95
 TSPlot.cxx:96
 TSPlot.cxx:97
 TSPlot.cxx:98
 TSPlot.cxx:99
 TSPlot.cxx:100
 TSPlot.cxx:101
 TSPlot.cxx:102
 TSPlot.cxx:103
 TSPlot.cxx:104
 TSPlot.cxx:105
 TSPlot.cxx:106
 TSPlot.cxx:107
 TSPlot.cxx:108
 TSPlot.cxx:109
 TSPlot.cxx:110
 TSPlot.cxx:111
 TSPlot.cxx:112
 TSPlot.cxx:113
 TSPlot.cxx:114
 TSPlot.cxx:115
 TSPlot.cxx:116
 TSPlot.cxx:117
 TSPlot.cxx:118
 TSPlot.cxx:119
 TSPlot.cxx:120
 TSPlot.cxx:121
 TSPlot.cxx:122
 TSPlot.cxx:123
 TSPlot.cxx:124
 TSPlot.cxx:125
 TSPlot.cxx:126
 TSPlot.cxx:127
 TSPlot.cxx:128
 TSPlot.cxx:129
 TSPlot.cxx:130
 TSPlot.cxx:131
 TSPlot.cxx:132
 TSPlot.cxx:133
 TSPlot.cxx:134
 TSPlot.cxx:135
 TSPlot.cxx:136
 TSPlot.cxx:137
 TSPlot.cxx:138
 TSPlot.cxx:139
 TSPlot.cxx:140
 TSPlot.cxx:141
 TSPlot.cxx:142
 TSPlot.cxx:143
 TSPlot.cxx:144
 TSPlot.cxx:145
 TSPlot.cxx:146
 TSPlot.cxx:147
 TSPlot.cxx:148
 TSPlot.cxx:149
 TSPlot.cxx:150
 TSPlot.cxx:151
 TSPlot.cxx:152
 TSPlot.cxx:153
 TSPlot.cxx:154
 TSPlot.cxx:155
 TSPlot.cxx:156
 TSPlot.cxx:157
 TSPlot.cxx:158
 TSPlot.cxx:159
 TSPlot.cxx:160
 TSPlot.cxx:161
 TSPlot.cxx:162
 TSPlot.cxx:163
 TSPlot.cxx:164
 TSPlot.cxx:165
 TSPlot.cxx:166
 TSPlot.cxx:167
 TSPlot.cxx:168
 TSPlot.cxx:169
 TSPlot.cxx:170
 TSPlot.cxx:171
 TSPlot.cxx:172
 TSPlot.cxx:173
 TSPlot.cxx:174
 TSPlot.cxx:175
 TSPlot.cxx:176
 TSPlot.cxx:177
 TSPlot.cxx:178
 TSPlot.cxx:179
 TSPlot.cxx:180
 TSPlot.cxx:181
 TSPlot.cxx:182
 TSPlot.cxx:183
 TSPlot.cxx:184
 TSPlot.cxx:185
 TSPlot.cxx:186
 TSPlot.cxx:187
 TSPlot.cxx:188
 TSPlot.cxx:189
 TSPlot.cxx:190
 TSPlot.cxx:191
 TSPlot.cxx:192
 TSPlot.cxx:193
 TSPlot.cxx:194
 TSPlot.cxx:195
 TSPlot.cxx:196
 TSPlot.cxx:197
 TSPlot.cxx:198
 TSPlot.cxx:199
 TSPlot.cxx:200
 TSPlot.cxx:201
 TSPlot.cxx:202
 TSPlot.cxx:203
 TSPlot.cxx:204
 TSPlot.cxx:205
 TSPlot.cxx:206
 TSPlot.cxx:207
 TSPlot.cxx:208
 TSPlot.cxx:209
 TSPlot.cxx:210
 TSPlot.cxx:211
 TSPlot.cxx:212
 TSPlot.cxx:213
 TSPlot.cxx:214
 TSPlot.cxx:215
 TSPlot.cxx:216
 TSPlot.cxx:217
 TSPlot.cxx:218
 TSPlot.cxx:219
 TSPlot.cxx:220
 TSPlot.cxx:221
 TSPlot.cxx:222
 TSPlot.cxx:223
 TSPlot.cxx:224
 TSPlot.cxx:225
 TSPlot.cxx:226
 TSPlot.cxx:227
 TSPlot.cxx:228
 TSPlot.cxx:229
 TSPlot.cxx:230
 TSPlot.cxx:231
 TSPlot.cxx:232
 TSPlot.cxx:233
 TSPlot.cxx:234
 TSPlot.cxx:235
 TSPlot.cxx:236
 TSPlot.cxx:237
 TSPlot.cxx:238
 TSPlot.cxx:239
 TSPlot.cxx:240
 TSPlot.cxx:241
 TSPlot.cxx:242
 TSPlot.cxx:243
 TSPlot.cxx:244
 TSPlot.cxx:245
 TSPlot.cxx:246
 TSPlot.cxx:247
 TSPlot.cxx:248
 TSPlot.cxx:249
 TSPlot.cxx:250
 TSPlot.cxx:251
 TSPlot.cxx:252
 TSPlot.cxx:253
 TSPlot.cxx:254
 TSPlot.cxx:255
 TSPlot.cxx:256
 TSPlot.cxx:257
 TSPlot.cxx:258
 TSPlot.cxx:259
 TSPlot.cxx:260
 TSPlot.cxx:261
 TSPlot.cxx:262
 TSPlot.cxx:263
 TSPlot.cxx:264
 TSPlot.cxx:265
 TSPlot.cxx:266
 TSPlot.cxx:267
 TSPlot.cxx:268
 TSPlot.cxx:269
 TSPlot.cxx:270
 TSPlot.cxx:271
 TSPlot.cxx:272
 TSPlot.cxx:273
 TSPlot.cxx:274
 TSPlot.cxx:275
 TSPlot.cxx:276
 TSPlot.cxx:277
 TSPlot.cxx:278
 TSPlot.cxx:279
 TSPlot.cxx:280
 TSPlot.cxx:281
 TSPlot.cxx:282
 TSPlot.cxx:283
 TSPlot.cxx:284
 TSPlot.cxx:285
 TSPlot.cxx:286
 TSPlot.cxx:287
 TSPlot.cxx:288
 TSPlot.cxx:289
 TSPlot.cxx:290
 TSPlot.cxx:291
 TSPlot.cxx:292
 TSPlot.cxx:293
 TSPlot.cxx:294
 TSPlot.cxx:295
 TSPlot.cxx:296
 TSPlot.cxx:297
 TSPlot.cxx:298
 TSPlot.cxx:299
 TSPlot.cxx:300
 TSPlot.cxx:301
 TSPlot.cxx:302
 TSPlot.cxx:303
 TSPlot.cxx:304
 TSPlot.cxx:305
 TSPlot.cxx:306
 TSPlot.cxx:307
 TSPlot.cxx:308
 TSPlot.cxx:309
 TSPlot.cxx:310
 TSPlot.cxx:311
 TSPlot.cxx:312
 TSPlot.cxx:313
 TSPlot.cxx:314
 TSPlot.cxx:315
 TSPlot.cxx:316
 TSPlot.cxx:317
 TSPlot.cxx:318
 TSPlot.cxx:319
 TSPlot.cxx:320
 TSPlot.cxx:321
 TSPlot.cxx:322
 TSPlot.cxx:323
 TSPlot.cxx:324
 TSPlot.cxx:325
 TSPlot.cxx:326
 TSPlot.cxx:327
 TSPlot.cxx:328
 TSPlot.cxx:329
 TSPlot.cxx:330
 TSPlot.cxx:331
 TSPlot.cxx:332
 TSPlot.cxx:333
 TSPlot.cxx:334
 TSPlot.cxx:335
 TSPlot.cxx:336
 TSPlot.cxx:337
 TSPlot.cxx:338
 TSPlot.cxx:339
 TSPlot.cxx:340
 TSPlot.cxx:341
 TSPlot.cxx:342
 TSPlot.cxx:343
 TSPlot.cxx:344
 TSPlot.cxx:345
 TSPlot.cxx:346
 TSPlot.cxx:347
 TSPlot.cxx:348
 TSPlot.cxx:349
 TSPlot.cxx:350
 TSPlot.cxx:351
 TSPlot.cxx:352
 TSPlot.cxx:353
 TSPlot.cxx:354
 TSPlot.cxx:355
 TSPlot.cxx:356
 TSPlot.cxx:357
 TSPlot.cxx:358
 TSPlot.cxx:359
 TSPlot.cxx:360
 TSPlot.cxx:361
 TSPlot.cxx:362
 TSPlot.cxx:363
 TSPlot.cxx:364
 TSPlot.cxx:365
 TSPlot.cxx:366
 TSPlot.cxx:367
 TSPlot.cxx:368
 TSPlot.cxx:369
 TSPlot.cxx:370
 TSPlot.cxx:371
 TSPlot.cxx:372
 TSPlot.cxx:373
 TSPlot.cxx:374
 TSPlot.cxx:375
 TSPlot.cxx:376
 TSPlot.cxx:377
 TSPlot.cxx:378
 TSPlot.cxx:379
 TSPlot.cxx:380
 TSPlot.cxx:381
 TSPlot.cxx:382
 TSPlot.cxx:383
 TSPlot.cxx:384
 TSPlot.cxx:385
 TSPlot.cxx:386
 TSPlot.cxx:387
 TSPlot.cxx:388
 TSPlot.cxx:389
 TSPlot.cxx:390
 TSPlot.cxx:391
 TSPlot.cxx:392
 TSPlot.cxx:393
 TSPlot.cxx:394
 TSPlot.cxx:395
 TSPlot.cxx:396
 TSPlot.cxx:397
 TSPlot.cxx:398
 TSPlot.cxx:399
 TSPlot.cxx:400
 TSPlot.cxx:401
 TSPlot.cxx:402
 TSPlot.cxx:403
 TSPlot.cxx:404
 TSPlot.cxx:405
 TSPlot.cxx:406
 TSPlot.cxx:407
 TSPlot.cxx:408
 TSPlot.cxx:409
 TSPlot.cxx:410
 TSPlot.cxx:411
 TSPlot.cxx:412
 TSPlot.cxx:413
 TSPlot.cxx:414
 TSPlot.cxx:415
 TSPlot.cxx:416
 TSPlot.cxx:417
 TSPlot.cxx:418
 TSPlot.cxx:419
 TSPlot.cxx:420
 TSPlot.cxx:421
 TSPlot.cxx:422
 TSPlot.cxx:423
 TSPlot.cxx:424
 TSPlot.cxx:425
 TSPlot.cxx:426
 TSPlot.cxx:427
 TSPlot.cxx:428
 TSPlot.cxx:429
 TSPlot.cxx:430
 TSPlot.cxx:431
 TSPlot.cxx:432
 TSPlot.cxx:433
 TSPlot.cxx:434
 TSPlot.cxx:435
 TSPlot.cxx:436
 TSPlot.cxx:437
 TSPlot.cxx:438
 TSPlot.cxx:439
 TSPlot.cxx:440
 TSPlot.cxx:441
 TSPlot.cxx:442
 TSPlot.cxx:443
 TSPlot.cxx:444
 TSPlot.cxx:445
 TSPlot.cxx:446
 TSPlot.cxx:447
 TSPlot.cxx:448
 TSPlot.cxx:449
 TSPlot.cxx:450
 TSPlot.cxx:451
 TSPlot.cxx:452
 TSPlot.cxx:453
 TSPlot.cxx:454
 TSPlot.cxx:455
 TSPlot.cxx:456
 TSPlot.cxx:457
 TSPlot.cxx:458
 TSPlot.cxx:459
 TSPlot.cxx:460
 TSPlot.cxx:461
 TSPlot.cxx:462
 TSPlot.cxx:463
 TSPlot.cxx:464
 TSPlot.cxx:465
 TSPlot.cxx:466
 TSPlot.cxx:467
 TSPlot.cxx:468
 TSPlot.cxx:469
 TSPlot.cxx:470
 TSPlot.cxx:471
 TSPlot.cxx:472
 TSPlot.cxx:473
 TSPlot.cxx:474
 TSPlot.cxx:475
 TSPlot.cxx:476
 TSPlot.cxx:477
 TSPlot.cxx:478
 TSPlot.cxx:479
 TSPlot.cxx:480
 TSPlot.cxx:481
 TSPlot.cxx:482
 TSPlot.cxx:483
 TSPlot.cxx:484
 TSPlot.cxx:485
 TSPlot.cxx:486
 TSPlot.cxx:487
 TSPlot.cxx:488
 TSPlot.cxx:489
 TSPlot.cxx:490
 TSPlot.cxx:491
 TSPlot.cxx:492
 TSPlot.cxx:493
 TSPlot.cxx:494
 TSPlot.cxx:495
 TSPlot.cxx:496
 TSPlot.cxx:497
 TSPlot.cxx:498
 TSPlot.cxx:499
 TSPlot.cxx:500
 TSPlot.cxx:501
 TSPlot.cxx:502
 TSPlot.cxx:503
 TSPlot.cxx:504
 TSPlot.cxx:505
 TSPlot.cxx:506
 TSPlot.cxx:507
 TSPlot.cxx:508
 TSPlot.cxx:509
 TSPlot.cxx:510
 TSPlot.cxx:511
 TSPlot.cxx:512
 TSPlot.cxx:513
 TSPlot.cxx:514
 TSPlot.cxx:515
 TSPlot.cxx:516
 TSPlot.cxx:517
 TSPlot.cxx:518
 TSPlot.cxx:519
 TSPlot.cxx:520
 TSPlot.cxx:521
 TSPlot.cxx:522
 TSPlot.cxx:523
 TSPlot.cxx:524
 TSPlot.cxx:525
 TSPlot.cxx:526
 TSPlot.cxx:527
 TSPlot.cxx:528
 TSPlot.cxx:529
 TSPlot.cxx:530
 TSPlot.cxx:531
 TSPlot.cxx:532
 TSPlot.cxx:533
 TSPlot.cxx:534
 TSPlot.cxx:535
 TSPlot.cxx:536
 TSPlot.cxx:537
 TSPlot.cxx:538
 TSPlot.cxx:539
 TSPlot.cxx:540
 TSPlot.cxx:541
 TSPlot.cxx:542
 TSPlot.cxx:543
 TSPlot.cxx:544
 TSPlot.cxx:545
 TSPlot.cxx:546
 TSPlot.cxx:547
 TSPlot.cxx:548
 TSPlot.cxx:549
 TSPlot.cxx:550
 TSPlot.cxx:551
 TSPlot.cxx:552
 TSPlot.cxx:553
 TSPlot.cxx:554
 TSPlot.cxx:555
 TSPlot.cxx:556
 TSPlot.cxx:557
 TSPlot.cxx:558
 TSPlot.cxx:559
 TSPlot.cxx:560
 TSPlot.cxx:561
 TSPlot.cxx:562
 TSPlot.cxx:563
 TSPlot.cxx:564
 TSPlot.cxx:565
 TSPlot.cxx:566
 TSPlot.cxx:567
 TSPlot.cxx:568
 TSPlot.cxx:569
 TSPlot.cxx:570
 TSPlot.cxx:571
 TSPlot.cxx:572
 TSPlot.cxx:573
 TSPlot.cxx:574
 TSPlot.cxx:575
 TSPlot.cxx:576
 TSPlot.cxx:577
 TSPlot.cxx:578
 TSPlot.cxx:579
 TSPlot.cxx:580
 TSPlot.cxx:581
 TSPlot.cxx:582
 TSPlot.cxx:583
 TSPlot.cxx:584
 TSPlot.cxx:585
 TSPlot.cxx:586
 TSPlot.cxx:587
 TSPlot.cxx:588
 TSPlot.cxx:589
 TSPlot.cxx:590
 TSPlot.cxx:591
 TSPlot.cxx:592
 TSPlot.cxx:593
 TSPlot.cxx:594
 TSPlot.cxx:595
 TSPlot.cxx:596
 TSPlot.cxx:597
 TSPlot.cxx:598
 TSPlot.cxx:599
 TSPlot.cxx:600
 TSPlot.cxx:601
 TSPlot.cxx:602
 TSPlot.cxx:603
 TSPlot.cxx:604
 TSPlot.cxx:605
 TSPlot.cxx:606
 TSPlot.cxx:607
 TSPlot.cxx:608
 TSPlot.cxx:609
 TSPlot.cxx:610
 TSPlot.cxx:611
 TSPlot.cxx:612
 TSPlot.cxx:613
 TSPlot.cxx:614
 TSPlot.cxx:615
 TSPlot.cxx:616
 TSPlot.cxx:617
 TSPlot.cxx:618
 TSPlot.cxx:619
 TSPlot.cxx:620
 TSPlot.cxx:621
 TSPlot.cxx:622
 TSPlot.cxx:623
 TSPlot.cxx:624
 TSPlot.cxx:625
 TSPlot.cxx:626
 TSPlot.cxx:627
 TSPlot.cxx:628
 TSPlot.cxx:629
 TSPlot.cxx:630
 TSPlot.cxx:631
 TSPlot.cxx:632
 TSPlot.cxx:633
 TSPlot.cxx:634
 TSPlot.cxx:635
 TSPlot.cxx:636
 TSPlot.cxx:637
 TSPlot.cxx:638
 TSPlot.cxx:639
 TSPlot.cxx:640
 TSPlot.cxx:641
 TSPlot.cxx:642
 TSPlot.cxx:643
 TSPlot.cxx:644
 TSPlot.cxx:645
 TSPlot.cxx:646
 TSPlot.cxx:647
 TSPlot.cxx:648
 TSPlot.cxx:649
 TSPlot.cxx:650
 TSPlot.cxx:651
 TSPlot.cxx:652
 TSPlot.cxx:653
 TSPlot.cxx:654
 TSPlot.cxx:655
 TSPlot.cxx:656
 TSPlot.cxx:657
 TSPlot.cxx:658
 TSPlot.cxx:659
 TSPlot.cxx:660
 TSPlot.cxx:661
 TSPlot.cxx:662
 TSPlot.cxx:663
 TSPlot.cxx:664
 TSPlot.cxx:665
 TSPlot.cxx:666
 TSPlot.cxx:667
 TSPlot.cxx:668
 TSPlot.cxx:669
 TSPlot.cxx:670
 TSPlot.cxx:671
 TSPlot.cxx:672
 TSPlot.cxx:673
 TSPlot.cxx:674
 TSPlot.cxx:675
 TSPlot.cxx:676
 TSPlot.cxx:677
 TSPlot.cxx:678
 TSPlot.cxx:679
 TSPlot.cxx:680
 TSPlot.cxx:681
 TSPlot.cxx:682
 TSPlot.cxx:683
 TSPlot.cxx:684
 TSPlot.cxx:685
 TSPlot.cxx:686
 TSPlot.cxx:687
 TSPlot.cxx:688
 TSPlot.cxx:689
 TSPlot.cxx:690
 TSPlot.cxx:691
 TSPlot.cxx:692
 TSPlot.cxx:693
 TSPlot.cxx:694
 TSPlot.cxx:695
 TSPlot.cxx:696
 TSPlot.cxx:697
 TSPlot.cxx:698
 TSPlot.cxx:699
 TSPlot.cxx:700
 TSPlot.cxx:701
 TSPlot.cxx:702
 TSPlot.cxx:703
 TSPlot.cxx:704
 TSPlot.cxx:705
 TSPlot.cxx:706
 TSPlot.cxx:707
 TSPlot.cxx:708
 TSPlot.cxx:709
 TSPlot.cxx:710
 TSPlot.cxx:711
 TSPlot.cxx:712
 TSPlot.cxx:713
 TSPlot.cxx:714
 TSPlot.cxx:715
 TSPlot.cxx:716
 TSPlot.cxx:717
 TSPlot.cxx:718
 TSPlot.cxx:719
 TSPlot.cxx:720
 TSPlot.cxx:721
 TSPlot.cxx:722
 TSPlot.cxx:723
 TSPlot.cxx:724
 TSPlot.cxx:725
 TSPlot.cxx:726
 TSPlot.cxx:727
 TSPlot.cxx:728
 TSPlot.cxx:729
 TSPlot.cxx:730
 TSPlot.cxx:731
 TSPlot.cxx:732
 TSPlot.cxx:733
 TSPlot.cxx:734
 TSPlot.cxx:735
 TSPlot.cxx:736
 TSPlot.cxx:737
 TSPlot.cxx:738
 TSPlot.cxx:739
 TSPlot.cxx:740
 TSPlot.cxx:741
 TSPlot.cxx:742
 TSPlot.cxx:743
 TSPlot.cxx:744
 TSPlot.cxx:745
 TSPlot.cxx:746
 TSPlot.cxx:747
 TSPlot.cxx:748
 TSPlot.cxx:749
 TSPlot.cxx:750
 TSPlot.cxx:751
 TSPlot.cxx:752
 TSPlot.cxx:753
 TSPlot.cxx:754
 TSPlot.cxx:755
 TSPlot.cxx:756
 TSPlot.cxx:757
 TSPlot.cxx:758
 TSPlot.cxx:759
 TSPlot.cxx:760
 TSPlot.cxx:761
 TSPlot.cxx:762
 TSPlot.cxx:763
 TSPlot.cxx:764
 TSPlot.cxx:765
 TSPlot.cxx:766
 TSPlot.cxx:767
 TSPlot.cxx:768
 TSPlot.cxx:769
 TSPlot.cxx:770
 TSPlot.cxx:771
 TSPlot.cxx:772
 TSPlot.cxx:773
 TSPlot.cxx:774
 TSPlot.cxx:775
 TSPlot.cxx:776
 TSPlot.cxx:777
 TSPlot.cxx:778
 TSPlot.cxx:779
 TSPlot.cxx:780
 TSPlot.cxx:781
 TSPlot.cxx:782
 TSPlot.cxx:783
 TSPlot.cxx:784
 TSPlot.cxx:785
 TSPlot.cxx:786
 TSPlot.cxx:787
 TSPlot.cxx:788
 TSPlot.cxx:789
 TSPlot.cxx:790
 TSPlot.cxx:791
 TSPlot.cxx:792
 TSPlot.cxx:793
 TSPlot.cxx:794
 TSPlot.cxx:795
 TSPlot.cxx:796
 TSPlot.cxx:797
 TSPlot.cxx:798
 TSPlot.cxx:799
 TSPlot.cxx:800
 TSPlot.cxx:801
 TSPlot.cxx:802
 TSPlot.cxx:803
 TSPlot.cxx:804
 TSPlot.cxx:805
 TSPlot.cxx:806
 TSPlot.cxx:807
 TSPlot.cxx:808
 TSPlot.cxx:809
 TSPlot.cxx:810
 TSPlot.cxx:811
 TSPlot.cxx:812
 TSPlot.cxx:813
 TSPlot.cxx:814
 TSPlot.cxx:815
 TSPlot.cxx:816
 TSPlot.cxx:817
 TSPlot.cxx:818
 TSPlot.cxx:819
 TSPlot.cxx:820
 TSPlot.cxx:821
 TSPlot.cxx:822
 TSPlot.cxx:823
 TSPlot.cxx:824
 TSPlot.cxx:825
 TSPlot.cxx:826
 TSPlot.cxx:827
 TSPlot.cxx:828
 TSPlot.cxx:829
 TSPlot.cxx:830
 TSPlot.cxx:831
 TSPlot.cxx:832
 TSPlot.cxx:833
 TSPlot.cxx:834
 TSPlot.cxx:835
 TSPlot.cxx:836
 TSPlot.cxx:837
 TSPlot.cxx:838
 TSPlot.cxx:839
 TSPlot.cxx:840
 TSPlot.cxx:841
 TSPlot.cxx:842
 TSPlot.cxx:843
 TSPlot.cxx:844
 TSPlot.cxx:845
 TSPlot.cxx:846
 TSPlot.cxx:847
 TSPlot.cxx:848
 TSPlot.cxx:849
 TSPlot.cxx:850
 TSPlot.cxx:851
 TSPlot.cxx:852
 TSPlot.cxx:853
 TSPlot.cxx:854
 TSPlot.cxx:855
 TSPlot.cxx:856
 TSPlot.cxx:857
 TSPlot.cxx:858
 TSPlot.cxx:859
 TSPlot.cxx:860
 TSPlot.cxx:861
 TSPlot.cxx:862
 TSPlot.cxx:863
 TSPlot.cxx:864
 TSPlot.cxx:865
 TSPlot.cxx:866
 TSPlot.cxx:867
 TSPlot.cxx:868
 TSPlot.cxx:869
 TSPlot.cxx:870
 TSPlot.cxx:871
 TSPlot.cxx:872
 TSPlot.cxx:873
 TSPlot.cxx:874
 TSPlot.cxx:875
 TSPlot.cxx:876
 TSPlot.cxx:877
 TSPlot.cxx:878
 TSPlot.cxx:879
 TSPlot.cxx:880
 TSPlot.cxx:881
 TSPlot.cxx:882
 TSPlot.cxx:883
 TSPlot.cxx:884
 TSPlot.cxx:885
 TSPlot.cxx:886
 TSPlot.cxx:887
 TSPlot.cxx:888
 TSPlot.cxx:889
 TSPlot.cxx:890
 TSPlot.cxx:891
 TSPlot.cxx:892
 TSPlot.cxx:893
 TSPlot.cxx:894
 TSPlot.cxx:895
 TSPlot.cxx:896
 TSPlot.cxx:897
 TSPlot.cxx:898
 TSPlot.cxx:899
 TSPlot.cxx:900
 TSPlot.cxx:901
 TSPlot.cxx:902
 TSPlot.cxx:903
 TSPlot.cxx:904
 TSPlot.cxx:905
 TSPlot.cxx:906
 TSPlot.cxx:907
 TSPlot.cxx:908
 TSPlot.cxx:909
 TSPlot.cxx:910
 TSPlot.cxx:911
 TSPlot.cxx:912
 TSPlot.cxx:913
 TSPlot.cxx:914
 TSPlot.cxx:915
 TSPlot.cxx:916
 TSPlot.cxx:917
 TSPlot.cxx:918
 TSPlot.cxx:919
 TSPlot.cxx:920
 TSPlot.cxx:921
 TSPlot.cxx:922
 TSPlot.cxx:923
 TSPlot.cxx:924
 TSPlot.cxx:925
 TSPlot.cxx:926
 TSPlot.cxx:927
 TSPlot.cxx:928
 TSPlot.cxx:929
 TSPlot.cxx:930
 TSPlot.cxx:931
 TSPlot.cxx:932
 TSPlot.cxx:933
 TSPlot.cxx:934
 TSPlot.cxx:935
 TSPlot.cxx:936
 TSPlot.cxx:937
 TSPlot.cxx:938
 TSPlot.cxx:939
 TSPlot.cxx:940
 TSPlot.cxx:941
 TSPlot.cxx:942
 TSPlot.cxx:943
 TSPlot.cxx:944
 TSPlot.cxx:945
 TSPlot.cxx:946
 TSPlot.cxx:947
 TSPlot.cxx:948
 TSPlot.cxx:949
 TSPlot.cxx:950
 TSPlot.cxx:951
 TSPlot.cxx:952
 TSPlot.cxx:953
 TSPlot.cxx:954
 TSPlot.cxx:955
 TSPlot.cxx:956
 TSPlot.cxx:957
 TSPlot.cxx:958
 TSPlot.cxx:959
 TSPlot.cxx:960
 TSPlot.cxx:961
 TSPlot.cxx:962
 TSPlot.cxx:963
 TSPlot.cxx:964
 TSPlot.cxx:965
 TSPlot.cxx:966
 TSPlot.cxx:967
 TSPlot.cxx:968
 TSPlot.cxx:969
 TSPlot.cxx:970
 TSPlot.cxx:971
 TSPlot.cxx:972
 TSPlot.cxx:973
 TSPlot.cxx:974
 TSPlot.cxx:975
 TSPlot.cxx:976
 TSPlot.cxx:977
 TSPlot.cxx:978
 TSPlot.cxx:979
 TSPlot.cxx:980
 TSPlot.cxx:981
 TSPlot.cxx:982
 TSPlot.cxx:983
 TSPlot.cxx:984
 TSPlot.cxx:985
 TSPlot.cxx:986
 TSPlot.cxx:987
 TSPlot.cxx:988
 TSPlot.cxx:989
 TSPlot.cxx:990
 TSPlot.cxx:991
 TSPlot.cxx:992
 TSPlot.cxx:993
 TSPlot.cxx:994
 TSPlot.cxx:995
 TSPlot.cxx:996
 TSPlot.cxx:997
 TSPlot.cxx:998
 TSPlot.cxx:999
 TSPlot.cxx:1000
 TSPlot.cxx:1001
 TSPlot.cxx:1002
 TSPlot.cxx:1003
 TSPlot.cxx:1004
 TSPlot.cxx:1005
 TSPlot.cxx:1006
 TSPlot.cxx:1007
 TSPlot.cxx:1008
 TSPlot.cxx:1009
 TSPlot.cxx:1010
 TSPlot.cxx:1011
 TSPlot.cxx:1012
 TSPlot.cxx:1013
 TSPlot.cxx:1014
 TSPlot.cxx:1015
 TSPlot.cxx:1016
 TSPlot.cxx:1017
 TSPlot.cxx:1018
 TSPlot.cxx:1019
 TSPlot.cxx:1020
 TSPlot.cxx:1021
 TSPlot.cxx:1022
 TSPlot.cxx:1023
 TSPlot.cxx:1024
 TSPlot.cxx:1025
 TSPlot.cxx:1026
 TSPlot.cxx:1027
 TSPlot.cxx:1028
 TSPlot.cxx:1029
 TSPlot.cxx:1030
 TSPlot.cxx:1031
 TSPlot.cxx:1032
 TSPlot.cxx:1033
 TSPlot.cxx:1034
 TSPlot.cxx:1035
 TSPlot.cxx:1036
 TSPlot.cxx:1037
 TSPlot.cxx:1038
 TSPlot.cxx:1039
 TSPlot.cxx:1040
 TSPlot.cxx:1041
 TSPlot.cxx:1042
 TSPlot.cxx:1043
 TSPlot.cxx:1044
 TSPlot.cxx:1045
 TSPlot.cxx:1046
 TSPlot.cxx:1047
 TSPlot.cxx:1048
 TSPlot.cxx:1049
 TSPlot.cxx:1050
 TSPlot.cxx:1051
 TSPlot.cxx:1052
 TSPlot.cxx:1053
 TSPlot.cxx:1054
 TSPlot.cxx:1055
 TSPlot.cxx:1056
 TSPlot.cxx:1057
 TSPlot.cxx:1058
 TSPlot.cxx:1059
 TSPlot.cxx:1060
 TSPlot.cxx:1061
 TSPlot.cxx:1062
 TSPlot.cxx:1063
 TSPlot.cxx:1064
 TSPlot.cxx:1065
 TSPlot.cxx:1066
 TSPlot.cxx:1067
 TSPlot.cxx:1068
 TSPlot.cxx:1069
 TSPlot.cxx:1070
 TSPlot.cxx:1071
 TSPlot.cxx:1072
 TSPlot.cxx:1073
 TSPlot.cxx:1074
 TSPlot.cxx:1075
 TSPlot.cxx:1076
 TSPlot.cxx:1077
 TSPlot.cxx:1078
 TSPlot.cxx:1079
 TSPlot.cxx:1080
 TSPlot.cxx:1081
 TSPlot.cxx:1082
 TSPlot.cxx:1083
 TSPlot.cxx:1084
 TSPlot.cxx:1085
 TSPlot.cxx:1086
 TSPlot.cxx:1087
 TSPlot.cxx:1088
 TSPlot.cxx:1089
 TSPlot.cxx:1090
 TSPlot.cxx:1091
 TSPlot.cxx:1092
 TSPlot.cxx:1093
 TSPlot.cxx:1094
 TSPlot.cxx:1095
 TSPlot.cxx:1096
 TSPlot.cxx:1097
 TSPlot.cxx:1098
 TSPlot.cxx:1099
 TSPlot.cxx:1100
 TSPlot.cxx:1101
 TSPlot.cxx:1102
 TSPlot.cxx:1103
 TSPlot.cxx:1104
 TSPlot.cxx:1105
 TSPlot.cxx:1106
 TSPlot.cxx:1107
 TSPlot.cxx:1108
 TSPlot.cxx:1109
 TSPlot.cxx:1110
 TSPlot.cxx:1111
 TSPlot.cxx:1112
 TSPlot.cxx:1113
 TSPlot.cxx:1114
 TSPlot.cxx:1115
 TSPlot.cxx:1116
 TSPlot.cxx:1117
 TSPlot.cxx:1118
 TSPlot.cxx:1119
 TSPlot.cxx:1120
 TSPlot.cxx:1121
 TSPlot.cxx:1122
 TSPlot.cxx:1123
 TSPlot.cxx:1124
 TSPlot.cxx:1125
 TSPlot.cxx:1126
 TSPlot.cxx:1127
 TSPlot.cxx:1128
 TSPlot.cxx:1129
 TSPlot.cxx:1130
 TSPlot.cxx:1131
 TSPlot.cxx:1132
 TSPlot.cxx:1133
 TSPlot.cxx:1134
 TSPlot.cxx:1135
 TSPlot.cxx:1136
 TSPlot.cxx:1137
 TSPlot.cxx:1138
 TSPlot.cxx:1139
 TSPlot.cxx:1140
 TSPlot.cxx:1141
 TSPlot.cxx:1142
 TSPlot.cxx:1143
 TSPlot.cxx:1144
 TSPlot.cxx:1145
 TSPlot.cxx:1146
 TSPlot.cxx:1147
 TSPlot.cxx:1148
 TSPlot.cxx:1149
 TSPlot.cxx:1150
 TSPlot.cxx:1151
 TSPlot.cxx:1152
 TSPlot.cxx:1153
 TSPlot.cxx:1154
 TSPlot.cxx:1155
 TSPlot.cxx:1156
 TSPlot.cxx:1157
 TSPlot.cxx:1158
 TSPlot.cxx:1159
 TSPlot.cxx:1160
 TSPlot.cxx:1161
 TSPlot.cxx:1162
 TSPlot.cxx:1163
 TSPlot.cxx:1164
 TSPlot.cxx:1165
 TSPlot.cxx:1166
 TSPlot.cxx:1167
 TSPlot.cxx:1168
 TSPlot.cxx:1169
 TSPlot.cxx:1170
 TSPlot.cxx:1171