// @(#)root/roostats:$Id$
// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
/*************************************************************************
 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/


/**
   HypoTestInverterResult class: holds the array of hypothesis test results and compute a confidence interval.
   Based on the RatioFinder code available in the RooStatsCms package developed by Gregory Schott and Danilo Piparo
   Ported and adapted to RooStats by Gregory Schott
   Some contributions to this class have been written by Matthias Wolf (error estimation)
**/


// include header file of this class 
#include "RooStats/HypoTestInverterResult.h"

#include "RooStats/HybridResult.h"
#include "RooStats/SamplingDistribution.h"
#include "RooStats/AsymptoticCalculator.h"
#include "RooMsgService.h"
#include "RooGlobalFunc.h"

#include "TF1.h"
#include "TGraphErrors.h"
#include <cmath>
#include "Math/BrentRootFinder.h"
#include "Math/WrappedFunction.h"

#include "TCanvas.h"
#include "RooStats/SamplingDistPlot.h"

#include <algorithm>

ClassImp(RooStats::HypoTestInverterResult)

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


// initialize static value 
double HypoTestInverterResult::fgAsymptoticMaxSigma      = 5; 


HypoTestInverterResult::HypoTestInverterResult(const char * name ) :
   SimpleInterval(name),
   fUseCLs(false),
   fIsTwoSided(false),
   fInterpolateLowerLimit(true),
   fInterpolateUpperLimit(true),
   fFittedLowerLimit(false),
   fFittedUpperLimit(false),
   fInterpolOption(kLinear),
   fLowerLimitError(-1),
   fUpperLimitError(-1),
   fCLsCleanupThreshold(0.005)
{
  // default constructor
   fLowerLimit = TMath::QuietNaN();
   fUpperLimit = TMath::QuietNaN();
}


HypoTestInverterResult::HypoTestInverterResult( const HypoTestInverterResult& other, const char * name ) :
   SimpleInterval(other,name),
   fUseCLs(other.fUseCLs),
   fIsTwoSided(other.fIsTwoSided),
   fInterpolateLowerLimit(other.fInterpolateLowerLimit),
   fInterpolateUpperLimit(other.fInterpolateUpperLimit),
   fFittedLowerLimit(other.fFittedLowerLimit),
   fFittedUpperLimit(other.fFittedUpperLimit),
   fInterpolOption(other.fInterpolOption),
   fLowerLimitError(other.fLowerLimitError),
   fUpperLimitError(other.fUpperLimitError),
   fCLsCleanupThreshold(other.fCLsCleanupThreshold)
{
   // copy constructor
   fLowerLimit = TMath::QuietNaN();
   fUpperLimit = TMath::QuietNaN();
   int nOther = other.ArraySize();

   fXValues = other.fXValues;
   for (int i = 0; i < nOther; ++i) 
     fYObjects.Add( other.fYObjects.At(i)->Clone() );
   for (int i = 0; i <  fExpPValues.GetSize() ; ++i)
     fExpPValues.Add( other.fExpPValues.At(i)->Clone() );   
}


HypoTestInverterResult&
HypoTestInverterResult::operator=(const HypoTestInverterResult& other)
{
  if (&other==this) {
    return *this ;
  } 

  SimpleInterval::operator = (other);
  fLowerLimit = other.fLowerLimit;
  fUpperLimit = other.fUpperLimit;
  fUseCLs = other.fUseCLs;
  fIsTwoSided = other.fIsTwoSided;
  fInterpolateLowerLimit = other.fInterpolateLowerLimit;
  fInterpolateUpperLimit = other.fInterpolateUpperLimit;
  fFittedLowerLimit = other.fFittedLowerLimit;
  fFittedUpperLimit = other.fFittedUpperLimit;
  fInterpolOption = other.fInterpolOption;
  fLowerLimitError = other.fLowerLimitError;
  fUpperLimitError = other.fUpperLimitError;
  fCLsCleanupThreshold = other.fCLsCleanupThreshold;

  int nOther = other.ArraySize();
  fXValues = other.fXValues;

  fYObjects.RemoveAll();
  for (int i=0; i < nOther; ++i) {
    fYObjects.Add( other.fYObjects.At(i)->Clone() );
  }
  fExpPValues.RemoveAll();
  for (int i=0; i <  fExpPValues.GetSize() ; ++i) {
    fExpPValues.Add( other.fExpPValues.At(i)->Clone() );   
  }
  
  return *this;
}


HypoTestInverterResult::HypoTestInverterResult( const char* name,
						const RooRealVar& scannedVariable,
						double cl ) :
   SimpleInterval(name,scannedVariable,TMath::QuietNaN(),TMath::QuietNaN(),cl), 
   fUseCLs(false),
   fIsTwoSided(false),
   fInterpolateLowerLimit(true),
   fInterpolateUpperLimit(true),
   fFittedLowerLimit(false),
   fFittedUpperLimit(false),
   fInterpolOption(kLinear),
   fLowerLimitError(-1),
   fUpperLimitError(-1),
   fCLsCleanupThreshold(0.005)
{
   // constructor 
   fYObjects.SetOwner();
   // put a cloned copy of scanned variable to set in the interval
   // to avoid I/O problem of the Result class - 
   // make the set owning the cloned copy (use clone istead of Clone to not copying all links)
   fParameters.removeAll();
   fParameters.takeOwnership();
   fParameters.addOwned(*((RooRealVar *) scannedVariable.clone(scannedVariable.GetName()) ));
}


int
HypoTestInverterResult::ExclusionCleanup()
{
  const int nEntries  = ArraySize();

  // initialization
  double nsig1(1.0);
  double nsig2(2.0);
  double p[5]; 
  double q[5];
  std::vector<double> qv;
  qv.resize(11,-1.0);

  p[0] = ROOT::Math::normal_cdf(-nsig2);
  p[1] = ROOT::Math::normal_cdf(-nsig1);
  p[2] = 0.5;
  p[3] = ROOT::Math::normal_cdf(nsig1);
  p[4] = ROOT::Math::normal_cdf(nsig2);

  bool resultIsAsymptotic(false);
  if (nEntries>=1) {
    HypoTestResult* r = dynamic_cast<HypoTestResult *> ( GetResult(0) );
    assert(r!=0);
    if ( !r->GetNullDistribution() && !r->GetAltDistribution() ) { 
      resultIsAsymptotic = true;
    }
  }

  int nPointsRemoved(0);

  double CLsobsprev(1.0);
  std::vector<double>::iterator itr = fXValues.begin();

  for (; itr!=fXValues.end();) {

    double x = (*itr);
    int i = FindIndex(x);
    //HypoTestResult* oneresult = GetResult(i);

    SamplingDistribution * s = GetExpectedPValueDist(i);
    if (!s) break; 

    /////////////////////////////////////////////////////////////////////////////////////////

    const std::vector<double> & values = s->GetSamplingDistribution();
    
    /// expected p-values
    // special case for asymptotic results (cannot use TMath::quantile in that case)
    if (resultIsAsymptotic) { 
      double maxSigma = 5; // == HypoTestInverterResult::fgAsymptoticMaxSigma; // MB: HACK
      double dsig = 2.*maxSigma / (values.size() -1) ;         
      int  i0 = (int) TMath::Floor ( ( -nsig2 + maxSigma )/dsig + 0.5 );
      int  i1 = (int) TMath::Floor ( ( -nsig1 + maxSigma )/dsig + 0.5 );
      int  i2 = (int) TMath::Floor ( ( maxSigma )/dsig + 0.5 );
      int  i3 = (int) TMath::Floor ( ( nsig1 + maxSigma )/dsig + 0.5 );
      int  i4 = (int) TMath::Floor ( ( nsig2 + maxSigma )/dsig + 0.5 );
      //
      q[0] = values[i0];
      q[1] = values[i1];
      q[2] = values[i2];
      q[3] = values[i3];
      q[4] = values[i4];
    } else { 
      double * z = const_cast<double *>( &values[0] ); // need to change TMath::Quantiles
      TMath::Quantiles(values.size(), 5, z, q, p, false);
    }
    
    delete s;

    /// store useful quantities for reuse later ...
    /// http://root.cern.ch/root/html532/src/RooStats__HypoTestInverterPlot.cxx.html#197
    for (int j=0; j<5; ++j) { qv[j]=q[j]; }
    qv[5]  = CLs(i) ; //
    qv[6]  = CLsError(i) ; //
    qv[7]  = CLb(i) ; //
    qv[8]  = CLbError(i) ; //
    qv[9]  = CLsplusb(i) ; //
    qv[10] = CLsplusbError(i) ; //
    double CLsobs = qv[5];

    /////////////////////////////////////////////////////////////////////////////////////////

    bool removeThisPoint(false);

    // 1. CLs should drop, else skip this point
    if (!removeThisPoint && resultIsAsymptotic && i>=1 && CLsobs>CLsobsprev) {
      //StatToolsLogger << kERROR << "Asymptotic. CLs not dropping: " << CLsobs << ". Remove this point." << GEndl;
      removeThisPoint = true; 
    } else { CLsobsprev = CLsobs; }
      
    // 2. CLs should not spike, else skip this point
    if (!removeThisPoint && i>=1 && CLsobs>=0.9999) {
      //StatToolsLogger << kERROR << "CLs spiking at 1.0: " << CLsobs << ". Remove this point." << GEndl;      
      removeThisPoint = true; 
    }
    // 3. Not interested in CLs values that become too low.
    if (!removeThisPoint && i>=1 && qv[4]<fCLsCleanupThreshold) { removeThisPoint = true; }

    // to remove or not to remove
    if (removeThisPoint) {
      itr = fXValues.erase(itr); // returned itr has been updated.     
      fYObjects.RemoveAt(i);
      fExpPValues.RemoveAt(i);
      nPointsRemoved++;
      continue;
    } else { // keep
      CLsobsprev = CLsobs;
      itr++;
    }
  }

  // after cleanup, reset existing limits
  fFittedUpperLimit = false; 
  fFittedLowerLimit = false;
  FindInterpolatedLimit(1-ConfidenceLevel(),true);

  return nPointsRemoved;
}


HypoTestInverterResult::~HypoTestInverterResult()
{
   // destructor
   // no need to delete explictly the objects in the TList since the TList owns the objects
}


bool HypoTestInverterResult::Add( const HypoTestInverterResult& otherResult   )
{
   // Merge this HypoTestInverterResult with another
   // HypoTestInverterResult passed as argument
   // The merge is done by combining the HypoTestResult when the same point value exist in both results.
   // If results exist at different points these are added in the new result
   // NOTE: Merging of the expected p-values obtained with pseudo-data.
   //  When expected p-values exist in the result (i.e. when rebuild option is used when getting the expected 
   // limit distribuiton in the HYpoTestInverter) then the expected p-values are also merged. This is equivalent 
   // at merging the pseudo-data. However there can be an inconsistency if the expected p-values have been 
   // obtained with different toys. In this case the merge is done but a warning message is printed.


   int nThis = ArraySize();
   int nOther = otherResult.ArraySize();
   if (nOther == 0) return true;
   if (nOther != otherResult.fYObjects.GetSize() ) return false; 
   if (nThis != fYObjects.GetSize() ) return false; 

   // cannot merge in case of inconsistent memebrs
   if (fExpPValues.GetSize() > 0 && fExpPValues.GetSize() != nThis ) return false;
   if (otherResult.fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() != nOther ) return false;

   oocoutI(this,Eval) << "HypoTestInverterResult::Add - merging result from " << otherResult.GetName()  
                                << " in " << GetName() << std::endl;

   bool addExpPValues = (fExpPValues.GetSize() == 0 && otherResult.fExpPValues.GetSize() > 0);
   bool mergeExpPValues = (fExpPValues.GetSize() > 0 && otherResult.fExpPValues.GetSize() > 0);

   if (addExpPValues || mergeExpPValues)
      oocoutI(this,Eval) << "HypoTestInverterResult::Add - merging also the expected p-values from pseudo-data" << std::endl;


   // case current result is empty
   // just make a simple copy of the other result
   if (nThis == 0) {
      fXValues = otherResult.fXValues;
      for (int i = 0; i < nOther; ++i) 
         fYObjects.Add( otherResult.fYObjects.At(i)->Clone() );
      for (int i = 0; i <  fExpPValues.GetSize() ; ++i)
         fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
   }
   // now do teh real mergemerge combining point with same value or adding extra ones 
   else { 
      for (int i = 0; i < nOther; ++i) {
         double otherVal = otherResult.fXValues[i];
         HypoTestResult * otherHTR = (HypoTestResult*) otherResult.fYObjects.At(i);
         if (otherHTR == 0) continue;
         bool sameXFound = false;
         for (int j = 0; j < nThis; ++j) {
            double thisVal = fXValues[j];
            
            // if same value merge the result 
            if ( (std::abs(otherVal) < 1  && TMath::AreEqualAbs(otherVal, thisVal,1.E-12) ) || 
                 (std::abs(otherVal) >= 1 && TMath::AreEqualRel(otherVal, thisVal,1.E-12) ) ) {
               HypoTestResult * thisHTR = (HypoTestResult*) fYObjects.At(j);               
               thisHTR->Append(otherHTR);
               sameXFound = true;
               if (mergeExpPValues) { 
                  ((SamplingDistribution*) fExpPValues.At(j))->Add( (SamplingDistribution*)otherResult.fExpPValues.At(i) );
                  // check if same toys have been used for the test statistic distribution
                  int thisNToys = (thisHTR->GetNullDistribution() ) ? thisHTR->GetNullDistribution()->GetSize() : 0;
                  int otherNToys = (otherHTR->GetNullDistribution() ) ? otherHTR->GetNullDistribution()->GetSize() : 0;
                  if (thisNToys != otherNToys ) 
                     oocoutW(this,Eval) << "HypoTestInverterResult::Add expexted p values have been generated with different toys " << thisNToys << " , " << otherNToys << std::endl;                  
               }
               break;
            }
         }
         if (!sameXFound) { 
            // add the new result 
            fYObjects.Add(otherHTR->Clone() );
            fXValues.push_back( otherVal);
         }
         // add in any case also when same x found
         if (addExpPValues)  
            fExpPValues.Add( otherResult.fExpPValues.At(i)->Clone() );
         
         
      }
   }
   
   if (ArraySize() > nThis) 
      oocoutI(this,Eval) << "HypoTestInverterResult::Add  - new number of points is " << fXValues.size()
                         << std::endl;
   else 
      oocoutI(this,Eval) << "HypoTestInverterResult::Add  - new toys/point is " 
                         <<  ((HypoTestResult*) fYObjects.At(0))->GetNullDistribution()->GetSize() 
                         << std::endl;
      
   return true;
}

bool HypoTestInverterResult::Add (Double_t x, const HypoTestResult & res) 
{
   // Add a single point result (an HypoTestResult)
   int i= FindIndex(x);
   if (i<0) {
      fXValues.push_back(x);
      fYObjects.Add(res.Clone());
   } else {
      HypoTestResult* r= GetResult(i);
      if (!r) return false;
      r->Append(&res);
   }
   return true;
}


 
double HypoTestInverterResult::GetXValue( int index ) const
{
 // function to return the value of the parameter of interest for the i^th entry in the results
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }

  return fXValues[index];
}

double HypoTestInverterResult::GetYValue( int index ) const
{
// function to return the value of the confidence level for the i^th entry in the results
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }

  if (fUseCLs) 
    return ((HypoTestResult*)fYObjects.At(index))->CLs();
  else 
     return ((HypoTestResult*)fYObjects.At(index))->CLsplusb();  // CLs+b
}

double HypoTestInverterResult::GetYError( int index ) const
{
// function to return the estimated error on the value of the confidence level for the i^th entry in the results
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }

  if (fUseCLs) 
    return ((HypoTestResult*)fYObjects.At(index))->CLsError();
  else 
    return ((HypoTestResult*)fYObjects.At(index))->CLsplusbError();
}

double HypoTestInverterResult::CLb( int index ) const
{
  // function to return the observed CLb value  for the i-th entry
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }
  return ((HypoTestResult*)fYObjects.At(index))->CLb();  
}

double HypoTestInverterResult::CLsplusb( int index ) const
{
  // function to return the observed CLs+b value  for the i-th entry
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }
  return ((HypoTestResult*)fYObjects.At(index))->CLsplusb();  
}
double HypoTestInverterResult::CLs( int index ) const
{
  // function to return the observed CLs value  for the i-th entry
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }
  return ((HypoTestResult*)fYObjects.At(index))->CLs();  
}

double HypoTestInverterResult::CLbError( int index ) const
{
  // function to return the error on the observed CLb value  for the i-th entry
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }
  return ((HypoTestResult*)fYObjects.At(index))->CLbError();  
}

double HypoTestInverterResult::CLsplusbError( int index ) const
{
  // function to return the error on the observed CLs+b value  for the i-th entry
  if ( index >= ArraySize() || index<0 ) {
    oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
    return -999;
  }
  return ((HypoTestResult*)fYObjects.At(index))->CLsplusbError();  
}
double HypoTestInverterResult::CLsError( int index ) const
{
   // function to return the error on the observed CLs value  for the i-th entry
   if ( index >= ArraySize() || index<0 ) {
      oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
      return -999;
   }
   return ((HypoTestResult*)fYObjects.At(index))->CLsError();  
}


HypoTestResult* HypoTestInverterResult::GetResult( int index ) const
{
   // get the HypoTestResult object at the given index point
   if ( index >= ArraySize() || index<0 ) {
      oocoutE(this,InputArguments) << "Problem: You are asking for an impossible array index value\n";
      return 0;
   }
   
   return ((HypoTestResult*) fYObjects.At(index));
}

int HypoTestInverterResult::FindIndex(double xvalue) const
{
   // find the index corresponding at the poi value xvalue
   // If no points is found return -1
   // Note that a tolerance is used of 10^-12 to find the closest point
  const double tol = 1.E-12;
  for (int i=0; i<ArraySize(); i++) {
     double xpoint = fXValues[i];
     if ( (std::abs(xvalue) > 1 && TMath::AreEqualRel( xvalue, xpoint, tol) ) ||
          (std::abs(xvalue) < 1 && TMath::AreEqualAbs( xvalue, xpoint, tol) ) )
        return i; 
  }
  return -1;
}


struct InterpolatedGraph { 
   InterpolatedGraph( const TGraph & g, double target, const char * interpOpt) : 
      fGraph(g), fTarget(target), fInterpOpt(interpOpt) {}

   // return interpolated value for x - target
   double operator() (double x) const { 
      return fGraph.Eval(x, (TSpline*) 0, fInterpOpt) - fTarget;
   }
   const TGraph & fGraph;
   double fTarget;
   TString fInterpOpt;
}; 

double HypoTestInverterResult::GetGraphX(const TGraph & graph, double y0, bool lowSearch, double &axmin, double &axmax) const  {
   // return the X value of the given graph for the target value y0
   // the graph is evaluated using linear interpolation by default. 
   // if option = "S" a TSpline3 is used 



   TString opt = "";
   if (fInterpolOption == kSpline)  opt = "S";

   InterpolatedGraph f(graph,y0,opt);
   ROOT::Math::BrentRootFinder brf;
   ROOT::Math::WrappedFunction<InterpolatedGraph> wf(f);

   // find reasanable xmin and xmax if not given
   const double * y = graph.GetY(); 
   int n = graph.GetN();
   if (n < 2) {
      ooccoutE(this,Eval) << "HypoTestInverterResult::GetGraphX - need at least 2 points for interpolation (n=" << n << ")\n";
      return (n>0) ?  y[0] : 0;
   } 

   double varmin = - TMath::Infinity();
   double varmax = TMath::Infinity();
   const RooRealVar* var = dynamic_cast<RooRealVar*>( fParameters.first() );
   if (var) { 
       varmin = var->getMin();
       varmax = var->getMax();
   }

   double xmin = axmin; 
   double xmax = axmax;
   // case no range is given check if need to extrapolate to lower/upper values 
   if (axmin >= axmax ) { 
      xmin = graph.GetX()[0];
      xmax = graph.GetX()[n-1];

      // case we have lower/upper limits

      // find ymin and ymax  and corresponding values 
      double ymin = TMath::MinElement(n,y);
      double ymax = TMath::MaxElement(n,y);
      // do lower extrapolation 
      if ( (ymax < y0 && !lowSearch) || ( ymin > y0 && lowSearch) ) { 
         xmin = varmin;
      }
      // do upper extrapolation  
      if ( (ymax < y0 && lowSearch) || ( ymin > y0 && !lowSearch) ) { 
         xmax = varmax;
      }
   }
   else { 

#ifdef ISREALLYNEEDED //??
      // in case of range given, check if all points are above or below the confidence level line
      bool isCross = false;
      bool first = true; 
      double prod = 1;
      for (int i = 0; i< n; ++i) { 
         double xp, yp;
         graph.GetPoint(i,xp,yp);
         if (xp >= xmin && xp <= xmax) { 
            prod *= TMath::Sign(1., (yp - y0) );
            if (prod < 0 && !first) { 
               isCross = true; 
               break;
            }
            first = false;
         }
      }
      if (!isCross) { 
         return (lowSearch) ? varmin : varmax;
      }
#endif
   }
   


   brf.SetFunction(wf,xmin,xmax);
   brf.SetNpx(20);
   bool ret = brf.Solve(100, 1.E-6, 1.E-6);
   if (!ret) { 
      ooccoutE(this,Eval) << "HypoTestInverterResult - interpolation failed - return inf" << std::endl;
         return TMath::Infinity();
   }
   double limit =  brf.Root();

//#define DO_DEBUG
#ifdef DO_DEBUG
   if (lowSearch) std::cout << "lower limit search : ";
   else std::cout << "Upper limit search :  ";
   std::cout << "do interpolation between " << xmin << "  and " << xmax << " limit is " << limit << std::endl;
#endif

   // look in case if a new interseption exists
   // only when boundaries are not given
   if (axmin >= axmax) { 
      int index = TMath::BinarySearch(n, graph.GetX(), limit);
#ifdef DO_DEBUG
   std::cout << "do new interpolation dividing from " << index << "  and " << y[index] << std::endl;
#endif

      if (lowSearch && index >= 1 && (y[0] - y0) * ( y[index]- y0) < 0) { 
         //search if  another interseption exists at a lower value
         limit =  GetGraphX(graph, y0, lowSearch, graph.GetX()[0], graph.GetX()[index] );
      }
      else if (!lowSearch && index < n-2 && (y[n-1] - y0) * ( y[index+1]- y0) < 0) { 
         // another interseption exists at an higher value
         limit =  GetGraphX(graph, y0, lowSearch, graph.GetX()[index+1], graph.GetX()[n-1] );
      }
   }
   // return also xmin, xmax values
   axmin = xmin;
   axmax = xmax;

   return limit;
}

double HypoTestInverterResult::FindInterpolatedLimit(double target, bool lowSearch, double xmin, double xmax )
{
   // interpolate to find a limit value
   // Use a linear or a spline interpolation depending on the interpolation option

   ooccoutD(this,Eval) << "HypoTestInverterResult - " 
                       << "Interpolate the upper limit between the 2 results closest to the target confidence level" 
                       << std::endl;

   // variable minimum and maximum
   double varmin = - TMath::Infinity();
   double varmax = TMath::Infinity();
   const RooRealVar* var = dynamic_cast<RooRealVar*>( fParameters.first() );
   if (var) { 
       varmin = var->getMin();
       varmax = var->getMax();
   }

   if (ArraySize()<2) {
      double val =  (lowSearch) ? xmin : xmax;
      oocoutW(this,Eval) << "HypoTestInverterResult::FindInterpolatedLimit" 
                         << " - not enough points to get the inverted interval - return " 
                         <<  val << std::endl;
      fLowerLimit = varmin;
      fUpperLimit = varmax; 
      return (lowSearch) ? fLowerLimit : fUpperLimit;  
   }

   // sort the values in x 
   int n = ArraySize();
   std::vector<unsigned int> index(n );
   TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false); 
   // make a graph with the sorted point
   TGraph graph(n);
   for (int i = 0; i < n; ++i) 
      graph.SetPoint(i, GetXValue(index[i]), GetYValue(index[i] ) );


   //std::cout << " search for " << lowSearch << std::endl;
   

   // search first for min/max in the given range
   if (xmin >= xmax) {
      

      // search for maximum between the point
      double * itrmax = std::max_element(graph.GetY() , graph.GetY() +n);
      double ymax = *itrmax; 
      int iymax = itrmax - graph.GetY(); 
      double xwithymax = graph.GetX()[iymax];

      //std::cout << " max of y " << iymax << "  " << xwithymax << "  " << ymax << " target is " << target << std::endl;

      // look if maximum is above/belove target
      if (ymax > target) {
         if (lowSearch)  {
            if ( iymax > 0) { 
                  // low search (minimmum is first point or minimum range)
               xmin = ( graph.GetY()[0] <= target ) ? graph.GetX()[0] : varmin;  
               xmax = xwithymax;
               } 
            else { 
               // no room for lower limit
               fLowerLimit = varmin; 
               lowSearch = false; // search now for upper limits
            }
         }
         if (!lowSearch ) {
            // up search 
            if ( iymax < n-1 ) { 
               xmin = xwithymax; 
               xmax = ( graph.GetY()[n-1] <= target ) ? graph.GetX()[n-1] : varmax;
            }
            else { 
               // no room for upper limit
               fUpperLimit = varmax; 
               lowSearch = true; // search now for lower limits
               xmin = varmin; 
               xmax = xwithymax;
               }
         }
      }
      else { 
         // in case is below the target
         // find out if is a lower or upper search
         if (iymax <= (n-1)/2 ) { 
            lowSearch = false; 
            fLowerLimit = varmin; 
         }
         else { 
            lowSearch = true; 
            fUpperLimit = varmax;
         }
      }

#ifdef DO_DEBUG   
      std::cout << " found xmin, xmax  = " << xmin << "  " << xmax << " for search " << lowSearch << std::endl;
#endif

      // now come here if I have already found a lower/upper limit 
      // i.e. I am calling routine for the second time
#ifdef ISNEEDED
      // should not really come here
      if (lowSearch &&  fUpperLimit < varmax) {
         xmin = fXValues[ index.front() ];
         // find xmax (is first point before upper limit)
         int upI = FindClosestPointIndex(target, 2, fUpperLimit);
         if (upI < 1) return xmin; 
         xmax = GetXValue(upI);
      }
      else if (!lowSearch && fLowerLimit > varmin ) {
         // find xmin (is first point after lower limit)
         int lowI = FindClosestPointIndex(target, 3, fLowerLimit);
         if (lowI >= n-1) return xmax;
         xmin = GetXValue(lowI);
         xmax = fXValues[ index.back() ];
      }
#endif
   }

#ifdef DO_DEBUG   
   std::cout << "finding " << lowSearch << " limit betweem " << xmin << "  " << xmax << endl;
#endif

   double limit =  GetGraphX(graph, target, lowSearch, xmin, xmax);
   if (lowSearch) fLowerLimit = limit; 
   else fUpperLimit = limit;
   CalculateEstimatedError( target, lowSearch, xmin, xmax);

#ifdef DO_DEBUG
   std::cout << "limit is " << limit << std::endl;
#endif

   if (lowSearch && !TMath::IsNaN(fUpperLimit)) return fLowerLimit;
   if (!lowSearch && !TMath::IsNaN(fLowerLimit)) return fUpperLimit;

   // now perform the opposite search on the complement interval
   if (lowSearch) { 
      xmin = xmax;
      xmax = varmax;         
   } else { 
      xmax = xmin;         
      xmin = varmin;
   }
   double limit2 =  GetGraphX(graph, target, !lowSearch, xmin, xmax);
   if (!lowSearch) fLowerLimit = limit2; 
   else fUpperLimit = limit2;

   CalculateEstimatedError( target, !lowSearch, xmin, xmax);

#ifdef DO_DEBUG
   std::cout << "other limit is " << limit2 << std::endl;
#endif
 
   return (lowSearch) ? fLowerLimit : fUpperLimit;

}


int HypoTestInverterResult::FindClosestPointIndex(double target, int mode, double xtarget)
{
   // if mode = 0 
   // find closest point to target in Y, the object closest to the target which is 3 sigma from the target 
   // and has smaller error
   // if mode = 1
   // find 2 closest point to target in X and between these two take the one closer to the target
   // if mode = 2  as in mode = 1 but return the lower point not the closest one
   // if mode = 3  as in mode = 1 but return the upper point not the closest one

  int bestIndex = -1;
  int closestIndex = -1;
  if (mode == 0) { 
     double smallestError = 2; // error must be < 1
     double bestValue = 2;
     for (int i=0; i<ArraySize(); i++) {
        double dist = fabs(GetYValue(i)-target);
        if ( dist <3 *GetYError(i) ) { // less than 1 sigma from target CL
           if (GetYError(i) < smallestError ) { 
              smallestError = GetYError(i);
              bestIndex = i; 
           }
        }
        if ( dist < bestValue) { 
           bestValue = dist;  
           closestIndex = i; 
        }
     }
     if (bestIndex >=0) return bestIndex;
     // if no points found just return the closest one to the target
     return closestIndex;
  }
  // else mode = 1,2,3
  // find the two closest points to limit value
  // sort the array first 
  int n = fXValues.size();
  std::vector<unsigned int> indx(n);
  TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false); 
  std::vector<double> xsorted( n);
  for (int i = 0; i < n; ++i) xsorted[i] = fXValues[indx[i] ];
  int index1 = TMath::BinarySearch( n, &xsorted[0], xtarget);

#ifdef DO_DEBUG
  std::cout << "finding closest point to " << xtarget << " is " << index1 << "  " << indx[index1] << std::endl; 
#endif

   // case xtarget is outside the range (bbefore or afterwards)
  if (index1 < 0) return indx[0]; 
  if (index1 >= n-1) return indx[n-1];  
  int index2 = index1 +1;
  
  if (mode == 2) return (GetXValue(indx[index1]) < GetXValue(indx[index2])) ? indx[index1] : indx[index2];
  if (mode == 3) return (GetXValue(indx[index1]) > GetXValue(indx[index2])) ? indx[index1] : indx[index2];
  // get smaller point of the two (mode == 1)
  if (fabs(GetYValue(indx[index1])-target) <= fabs(GetYValue(indx[index2])-target) )
     return indx[index1];
  return indx[index2];

}
        

Double_t HypoTestInverterResult::LowerLimit()
{
  if (fFittedLowerLimit) return fLowerLimit;
  //std::cout << "finding point with cl = " << 1-(1-ConfidenceLevel())/2 << endl;
  if ( fInterpolateLowerLimit ) { 
     // find both lower/upper limit
     if (TMath::IsNaN(fLowerLimit) )  FindInterpolatedLimit(1-ConfidenceLevel(),true);
  } else {
     //LM: I think this is never called
     fLowerLimit = GetXValue( FindClosestPointIndex((1-ConfidenceLevel())) );
  }
  return fLowerLimit;
}

Double_t HypoTestInverterResult::UpperLimit()
{
   //std::cout << "finding point with cl = " << (1-ConfidenceLevel())/2 << endl;
  if (fFittedUpperLimit) return fUpperLimit;
  if ( fInterpolateUpperLimit ) {
     if (TMath::IsNaN(fUpperLimit) )  FindInterpolatedLimit(1-ConfidenceLevel(),false);
  } else {
     //LM: I think this is never called
     fUpperLimit = GetXValue( FindClosestPointIndex((1-ConfidenceLevel())) );
  }
  return fUpperLimit;
}

Double_t HypoTestInverterResult::CalculateEstimatedError(double target, bool lower, double xmin, double xmax)
{
  // Return an error estimate on the upper(lower) limit.  This is the error on
  // either CLs or CLsplusb divided by an estimate of the slope at this
  // point.
   
  if (ArraySize()==0) {
     oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimateError" 
                        << "Empty result \n";
    return 0;
  }

  if (ArraySize()<2) {
     oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimateError" 
                        << " only  points - return its error\n"; 
     return GetYError(0);
  }

  // it does not make sense in case of asymptotic which do not have point errors 
  if (!GetNullTestStatDist(0) ) return 0; 
 
  TString type = (!lower) ? "upper" : "lower";

#ifdef DO_DEBUG
  std::cout << "calculate estimate error " << type << " between " << xmin << " and " << xmax << std::endl;
  std::cout << "computed limit is " << ( (lower) ? fLowerLimit : fUpperLimit ) << std::endl;
#endif

    // make a TGraph Errors with the sorted points
  std::vector<unsigned int> indx(fXValues.size());
  TMath::SortItr(fXValues.begin(), fXValues.end(), indx.begin(), false); 
  // make a graph with the sorted point
  TGraphErrors graph;
  int ip = 0, np = 0;
  for (int i = 0; i < ArraySize(); ++i) {
     if ( (xmin < xmax) && ( GetXValue(indx[i]) >= xmin && GetXValue(indx[i]) <= xmax) ) {
        np++;
        // exclude points with zero or very small errors 
        if (GetYError(indx[i] ) > 1.E-6) {  
           graph.SetPoint(ip, GetXValue(indx[i]), GetYValue(indx[i] ) );
           graph.SetPointError(ip, 0.,  GetYError(indx[i]) );
           ip++;
        }
     }
  }
  if (graph.GetN() < 2) {
     if (np >= 2) oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimatedError - no valid points - cannot estimate  the " << type << " limit error " << std::endl;
     return 0; 
  }

  double minX = xmin;
  double maxX = xmax;
  if (xmin >= xmax) {
     minX = fXValues[ indx.front() ];
     maxX = fXValues[ indx.back() ];
  }



  TF1 fct("fct", "exp([0] * (x - [2] ) + [1] * (x-[2])**2)", minX, maxX);
  double scale = maxX-minX; 
  if (lower) { 
     fct.SetParameters( 2./scale, 0.1/scale, graph.GetX()[0] ); 
     fct.SetParLimits(0,0,100./scale); 
     fct.SetParLimits(1,0, 10./scale); }
  else  { 
     fct.SetParameters( -2./scale, -0.1/scale ); 
     fct.SetParLimits(0,-100./scale, 0); 
     fct.SetParLimits(1,-100./scale, 0); }

  if (graph.GetN() < 3) fct.FixParameter(1,0.);

  // find the point closest to the limit
  double limit = (!lower) ? fUpperLimit : fLowerLimit;
  if (TMath::IsNaN(limit)) return 0; // cannot do if limit not computed  


#ifdef DO_DEBUG
  TCanvas * c1 = new TCanvas(); 
  std::cout << "fitting for limit " << type << "between " << minX << " , " << maxX << " points considered " << graph.GetN() <<  std::endl;
  int fitstat = graph.Fit(&fct," EX0");
  graph.SetMarkerStyle(20);
  graph.Draw("AP");
  graph.Print(); 
  c1->SaveAs(TString::Format("graphFit_%s.pdf",type.Data()) );
  delete c1; 
#else
  int fitstat = graph.Fit(&fct,"Q EX0");
#endif 

  int index = FindClosestPointIndex(target, 1, limit);
  double theError = 0;
  if (fitstat == 0) { 
     double errY = GetYError(index);
     if (errY >  0) {
        double m = fct.Derivative( GetXValue(index) );
        theError = std::min(fabs( GetYError(index) / m), maxX-minX);     
     }
  }
  else { 
     oocoutW(this,Eval) << "HypoTestInverterResult::CalculateEstimatedError - cannot estimate  the " << type << " limit error " << std::endl;
     theError = 0; 
  }
  if (lower) 
     fLowerLimitError = theError;
  else
     fUpperLimitError = theError;

#ifdef DO_DEBUG
  std::cout << "closes point to the limit is " << index << "  " << GetXValue(index) << " and has error " << GetYError(index) << std::endl;
#endif

  return theError;
}


Double_t HypoTestInverterResult::LowerLimitEstimatedError()
{
   // need to have compute first lower limit
   if (TMath::IsNaN(fLowerLimit) ) LowerLimit();
   if (fLowerLimitError >= 0) return fLowerLimitError;    
   // try to recompute the error
   return CalculateEstimatedError(1-ConfidenceLevel(), true);
}


Double_t HypoTestInverterResult::UpperLimitEstimatedError()
{
   if (TMath::IsNaN(fUpperLimit) ) UpperLimit();
   if (fUpperLimitError >= 0) return fUpperLimitError;
   // try to recompute the error
   return CalculateEstimatedError(1-ConfidenceLevel(), false);
}

SamplingDistribution *  HypoTestInverterResult::GetBackgroundTestStatDist(int index ) const { 
   // get the background test statistic distribution

   HypoTestResult * firstResult = (HypoTestResult*) fYObjects.At(index);
   if (!firstResult) return 0;
   return firstResult->GetBackGroundIsAlt() ? firstResult->GetAltDistribution() : firstResult->GetNullDistribution();
}

SamplingDistribution *  HypoTestInverterResult::GetSignalAndBackgroundTestStatDist(int index) const { 
   // get the signal and background test statistic distribution
   HypoTestResult * result = (HypoTestResult*) fYObjects.At(index);
   if (!result) return 0;
   return !result->GetBackGroundIsAlt() ? result->GetAltDistribution() : result->GetNullDistribution();       
}

SamplingDistribution *  HypoTestInverterResult::GetExpectedPValueDist(int index) const { 
   // get the expected p-value distribution at the scanned point index

   if (index < 0 || index >=  ArraySize() ) return 0; 

   if (fExpPValues.GetSize() == ArraySize()  ) {
      return (SamplingDistribution*)  fExpPValues.At(index)->Clone();
   }
   
   static bool useFirstB = false;
   // get the S+B distribution 
   int bIndex = (useFirstB) ? 0 : index;
 
   SamplingDistribution * bDistribution = GetBackgroundTestStatDist(bIndex);
   SamplingDistribution * sbDistribution = GetSignalAndBackgroundTestStatDist(index);

   HypoTestResult * result = (HypoTestResult*) fYObjects.At(index);

   if (bDistribution && sbDistribution) { 

      // create a new HypoTestResult
      HypoTestResult tempResult; 
      tempResult.SetPValueIsRightTail( result->GetPValueIsRightTail() );
      tempResult.SetBackgroundAsAlt( true);
      // ownership of SamplingDistribution is in HypoTestResult class now
      tempResult.SetNullDistribution( new SamplingDistribution(*sbDistribution) );
      tempResult.SetAltDistribution( new SamplingDistribution(*bDistribution ) );
      
      std::vector<double> values(bDistribution->GetSize()); 
      for (int i = 0; i < bDistribution->GetSize(); ++i) { 
         tempResult.SetTestStatisticData( bDistribution->GetSamplingDistribution()[i] );
         values[i] = (fUseCLs) ? tempResult.CLs() : tempResult.CLsplusb();
      }
      return new SamplingDistribution("expected values","expected values",values);
   }  
   // in case b abs sbDistribution are null assume is coming from the asymptotic calculator 
   // hard -coded this value (no really needed to be used by user)
   fgAsymptoticMaxSigma = 5; 
   const int npoints = 11;
   const double smax = fgAsymptoticMaxSigma;
   const double dsig = 2* smax/ (npoints-1) ;
   std::vector<double> values(npoints);
   for (int i = 0; i < npoints; ++i) { 
      double nsig = -smax + dsig*i; 
      double pval = AsymptoticCalculator::GetExpectedPValues( result->NullPValue(), result->AlternatePValue(), nsig, fUseCLs, !fIsTwoSided);
      if (pval < 0) { return 0;}
      values[i] = pval;
   }
   return new SamplingDistribution("Asymptotic expected values","Asymptotic expected values",values);
   
}



SamplingDistribution *  HypoTestInverterResult::GetLimitDistribution(bool lower ) const { 
   // get the limit distribution (lower/upper depending on the flag)
   // by interpolating  the expected p values for each point

   
   if (ArraySize()<2) {
      oocoutE(this,Eval) << "HypoTestInverterResult::GetLimitDistribution" 
                         << " not  enought points -  return 0 " << std::endl; 
      return 0; 
   }
   
   ooccoutD(this,Eval) << "HypoTestInverterResult - computing  limit distribution...." << std::endl;

      
   // find optimal size by looking at the PValue distribution obtained
   int npoints = ArraySize();
   std::vector<SamplingDistribution*> distVec( npoints );
   double sum = 0;
   for (unsigned int i = 0; i < distVec.size(); ++i) {
      distVec[i] =  GetExpectedPValueDist(i);
      // sort the distributions
      // hack (by calling InverseCDF(0) will sort the sampling distribution
      if (distVec[i] ) { 
         distVec[i]->InverseCDF(0);
         sum += distVec[i]->GetSize();     
      }
   }
   int  size =  int( sum/ npoints);
  
   if (size < 10) { 
      ooccoutW(this,InputArguments) << "HypoTestInverterResult - set a minimum size of 10 for limit distribution" <<   std::endl;
      size = 10;
   }


  double target = 1-fConfidenceLevel;

  // vector with the quantiles of the p-values for each scanned poi point
  std::vector< std::vector<double>  > quantVec(npoints );
  for (int i = 0; i <  npoints; ++i) {

     if (!distVec[i]) continue;

     // make quantiles from the sampling distributions of the expected p values 
     std::vector<double> pvalues = distVec[i]->GetSamplingDistribution();
     delete distVec[i];  distVec[i] = 0;
     std::sort(pvalues.begin(), pvalues.end());
     // find the quantiles of the distribution 
     double p[1] = {0}; 
     double q[1] = {0};

     quantVec[i] = std::vector<double>(size);
     for (int ibin = 0; ibin < size; ++ibin) { 
        // exclude for a bug in TMath::Quantiles last item 
        p[0] = std::min( (ibin+1) * 1./double(size), 1.0);
        // use the type 1 which give the point value
        TMath::Quantiles(pvalues.size(), 1, &pvalues[0], q, p, true, 0, 1 );
        (quantVec[i])[ibin] = q[0]; 
     } 
  
  }

  // sort the values in x 
  std::vector<unsigned int> index( npoints );
  TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false); 

  // SamplingDistribution * dist0 = distVec[index[0]];
  // SamplingDistribution * dist1 = distVec[index[1]];

  
  std::vector<double> limits(size);
  // loop on the p values and find the limit for each expected point in the quantiles vector  
  for (int j = 0; j < size; ++j ) {

     TGraph g;
     for (int k = 0; k < npoints ; ++k) { 
        if (quantVec[index[k]].size()  > 0 )
           g.SetPoint(g.GetN(), GetXValue(index[k]), (quantVec[index[k]])[j] );
     }

     limits[j] = GetGraphX(g, target, lower);

  }


  if (lower) 
     return new SamplingDistribution("Expected lower Limit","Expected lower limits",limits);
  else
     return new SamplingDistribution("Expected upper Limit","Expected upper limits",limits);
  
}


double  HypoTestInverterResult::GetExpectedLowerLimit(double nsig, const char * opt ) const {
   // Get the expected lower limit  
   // nsig is used to specify which expected value of the UpperLimitDistribution
   // For example 
   // nsig = 0 (default value) returns the expected value 
   // nsig = -1 returns the lower band value at -1 sigma 
   // nsig + 1 return the upper value
   // opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
   //                         and then find the quantiles in the limit distribution 
   // ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
   // interpolates them
   
   return GetExpectedLimit(nsig, true, opt);
} 

double  HypoTestInverterResult::GetExpectedUpperLimit(double nsig, const char * opt ) const {
   // Get the expected upper limit  
   // nsig is used to specify which expected value of the UpperLimitDistribution
   // For example 
   // nsig = 0 (default value) returns the expected value 
   // nsig = -1 returns the lower band value at -1 sigma 
   // nsig + 1 return the upper value
   // opt is an option specifying the type of method used for computing the upper limit
   // opt = "" (default) : compute limit by interpolating all the p values, find the corresponding limit distribution
   //                         and then find the quantiles in the limit distribution 
   // ioption = "P" is the method used for plotting. One Finds the corresponding nsig quantile in the p values and then
   // interpolates them
   //              
   
   return GetExpectedLimit(nsig, false, opt);
} 


   
double  HypoTestInverterResult::GetExpectedLimit(double nsig, bool lower, const char * opt ) const {
   // get expected limit (lower/upper) depending on the flag 
   // for asymptotic is a special case (the distribution is generated an step in sigma values)
   // distringuish asymptotic looking at the hypotest results
   // if option = "P" get expected limit using directly quantiles of p value distribution 
   // else (default) find expected limit by obtaining first a full limit distributions
   // The last one is in general more correct

   const int nEntries = ArraySize();
   if (nEntries <= 0)  return (lower) ? 1 : 0;  // return 1 for lower, 0 for upper
   
   HypoTestResult * r = dynamic_cast<HypoTestResult *> (fYObjects.First() );
   assert(r != 0);
   if (!r->GetNullDistribution() && !r->GetAltDistribution() ) { 
      // we are in the asymptotic case 
      // get the limits obtained at the different sigma values 
      SamplingDistribution * limitDist = GetLimitDistribution(lower); 
      if (!limitDist) return 0;
      const std::vector<double> & values = limitDist->GetSamplingDistribution();
      if (values.size() <= 1) return 0; 
      double dsig = 2* fgAsymptoticMaxSigma/ (values.size() -1) ;
      int  i = (int) TMath::Floor ( (nsig +  fgAsymptoticMaxSigma)/dsig + 0.5);
      return values[i];
   }

   double p[1] = {0};
   double q[1] = {0};
   p[0] = ROOT::Math::normal_cdf(nsig,1);

   // for CLs+b can get the quantiles of p-value distribution and 
   // interpolate them
   // In the case of CLs (since it is not a real p-value anymore but a ratio) 
   // then it is needed to get first all limit distribution values and then the quantiles 

   TString option(opt);
   option.ToUpper();
   if (option.Contains("P")) {

      TGraph  g;   

      // sort the arrays based on the x values
      std::vector<unsigned int> index(nEntries);
      TMath::SortItr(fXValues.begin(), fXValues.end(), index.begin(), false);

      for (int j=0; j<nEntries; ++j) {
         int i = index[j]; // i is the order index 
         SamplingDistribution * s = GetExpectedPValueDist(i);
         if (!s) { 
            ooccoutI(this,Eval) << "HypoTestInverterResult - cannot compute expected p value distribution for point, x = " 
                                << GetXValue(i)  << " skip it " << std::endl;
            continue;
         } 
         const std::vector<double> & values = s->GetSamplingDistribution();
         double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
         TMath::Quantiles(values.size(), 1, x,q,p,false);
         g.SetPoint(g.GetN(), fXValues[i], q[0] );
         delete s;
      }
      if (g.GetN() < 2) { 
         ooccoutE(this,Eval) << "HypoTestInverterResult - cannot compute limits , not enough points, n =  " << g.GetN() << std::endl;
         return 0;
      }

      // interpolate the graph to obtain the limit
      double target = 1-fConfidenceLevel;
      return GetGraphX(g, target, lower);

   }
   // here need to use the limit distribution 
   SamplingDistribution * limitDist = GetLimitDistribution(lower); 
   if (!limitDist) return 0;
   const std::vector<double> & values = limitDist->GetSamplingDistribution();
   double * x = const_cast<double *>(&values[0]); // need to change TMath::Quantiles
   TMath::Quantiles(values.size(), 1, x,q,p,false);
   return q[0];
   
}

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