// @(#)root/roostats:$Id$
// Author: Kyle Cranmer, Sven Kreiss   23/05/10
/*************************************************************************
 * 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.             *
 *************************************************************************/

/**
   Performs hypothesis tests using aysmptotic formula for the profile likelihood and 
   Asimov data set
*/


#include "RooStats/AsymptoticCalculator.h"
#include "RooStats/ToyMCSampler.h"
#include "RooStats/ModelConfig.h"
#include "RooStats/ProfileLikelihoodTestStat.h"
#include "RooStats/RooStatsUtils.h"

#include "RooArgSet.h"
#include "RooArgList.h"
#include "RooProdPdf.h"
#include "RooSimultaneous.h"
#include "RooDataSet.h"
#include "RooCategory.h"
#include "RooRealVar.h"
#include "RooMinimizer.h"
#include "RooFitResult.h"
#include "RooNLLVar.h"
#include "Math/MinimizerOptions.h"
#include "RooPoisson.h"
#include "RooUniform.h"
#include "RooGamma.h"
#include "RooGaussian.h"
#include "RooBifurGauss.h"
#include "RooLognormal.h"
#include "RooDataHist.h"
#include <cmath>
#include <typeinfo>

#include "Math/BrentRootFinder.h"
#include "Math/WrappedFunction.h"

#include "TStopwatch.h"

using namespace RooStats;
using namespace std;


ClassImp(RooStats::AsymptoticCalculator);

int AsymptoticCalculator::fgPrintLevel = 1;


void AsymptoticCalculator::SetPrintLevel(int level) { 
   // set print level (static function)
   // 0 minimal, 1 normal,  2 debug
   fgPrintLevel = level;
}


AsymptoticCalculator::AsymptoticCalculator(
   RooAbsData &data,
   const ModelConfig &altModel,
   const ModelConfig &nullModel, bool nominalAsimov) :
      HypoTestCalculatorGeneric(data, altModel, nullModel, 0), 
      fOneSided(false), fOneSidedDiscovery(false), fNominalAsimov(nominalAsimov),
      fUseQTilde(-1), 
      fNLLObs(0), fNLLAsimov(0), 
      fAsimovData(0)   
{
   // constructor for asymptotic calculator from Data set  and ModelConfig
   if (!Initialize()) return; 

   int verbose = fgPrintLevel; 
   // try to guess default configuration
   // (this part should be only in constructor because the null snapshot might change during HypoTestInversion
   const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
   assert(nullSnapshot);
   RooRealVar * muNull  = dynamic_cast<RooRealVar*>( nullSnapshot->first() );
   assert(muNull);
   if (muNull->getVal() == muNull->getMin()) { 
      fOneSidedDiscovery = true; 
      if (verbose > 0) 
         oocoutI((TObject*)0,InputArguments) << "AsymptotiCalculator: Minimum of POI is " << muNull->getMin() << " corresponds to null  snapshot   - default configuration is  one-sided discovery formulae  " << std::endl;
   }

}

   
bool AsymptoticCalculator::Initialize() const {
   // Initialize the calculator 
   // The initialization will perform a global fit of the model to the data 
   // and build an Asimov data set. 
   // It will then also fit the model to the Asimov data set to find the likelihood value  
   // of the Asimov data set
   // nominalAsimov is an option for using Asimov data set obtained using nominal nuisance parameter values 
   // By default the nuisance parameters are fitted to the data  
   // NOTE: If a fit has been done before, one for speeding up could set all the initial prameters 
   // to the fit value and in addition set the null snapshot to the best fit

   int verbose = fgPrintLevel; 
   if (verbose >= 0)
      oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize...." << std::endl;


   RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
   if (!nullPdf) {
      oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has not a pdf defined" << std::endl;
      return false;       
   }
   RooAbsData * obsData = const_cast<RooAbsData *>(GetData() );
   if (!obsData ) {
      oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - data set has not been defined" << std::endl;
      return false;       
   }
   RooAbsData & data = *obsData;

   

   const RooArgSet * poi = GetNullModel()->GetParametersOfInterest(); 
   if (!poi || poi->getSize() == 0) { 
      oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize -  ModelConfig has not POI defined." << endl;
      return false;
   }
   if (poi->getSize() > 1) { 
      oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has more than one POI defined \n\t" 
                                          << "The asymptotic calculator works for only one POI - consider as POI only the first parameter" 
                                          << std::endl;
   }
 

   // This will set the poi value to the null snapshot value in the ModelConfig
   const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
   if(nullSnapshot == NULL || nullSnapshot->getSize() == 0) {
      oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::Initialize - Null model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
      return false;
   }
   
   // GetNullModel()->Print();
   // printf("ASymptotic calc: null snapshot\n");
   // nullSnapshot->Print("v");
   // printf("PDF  variables " );
   // nullPdf->getVariables()->Print("v");
   
   // keep snapshot for the initial parameter values (need for nominal Asimov)
   RooArgSet nominalParams; 
   RooArgSet * allParams = nullPdf->getParameters(data);
   RemoveConstantParameters(allParams);
   if (fNominalAsimov) { 
      allParams->snapshot(nominalParams);
   }
   fBestFitPoi.removeAll(); 
   fBestFitParams.removeAll();
   fAsimovGlobObs.removeAll();
      
   // evaluate the unconditional nll for the full model on the  observed data 
   if (verbose >= 0)
      oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize - Find  best unconditional NLL on observed data" << endl;
   fNLLObs = EvaluateNLL( *nullPdf, data, GetNullModel()->GetConditionalObservables());
   // fill also snapshot of best poi
   poi->snapshot(fBestFitPoi);
   RooRealVar * muBest = dynamic_cast<RooRealVar*>(fBestFitPoi.first());
   assert(muBest);
   if (verbose >= 0)  
      oocoutP((TObject*)0,Eval) << "Best fitted POI value = " << muBest->getVal() << " +/- " << muBest->getError() << std::endl;   
   // keep snapshot of all best fit parameters
   allParams->snapshot(fBestFitParams);
   delete allParams;
   
   // compute Asimov data set for the background (alt poi ) value
   const RooArgSet * altSnapshot = GetAlternateModel()->GetSnapshot();
   if(altSnapshot == NULL || altSnapshot->getSize() == 0) {
      oocoutE((TObject*)0,InputArguments) << "Alt (Background)  model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
      return false;
   }

   RooArgSet poiAlt(*altSnapshot);  // this is the poi snapshot of B (i.e. for mu=0)

   oocoutP((TObject*)0,Eval) << "AsymptoticCalculator: Building Asimov data Set" << endl;

   // check that in case of binned models th ennumber of bins of the observables are consistent 
   // with the number of bins  in the observed data 
   // This number will be used for making the Asimov data set so it will be more consistent with the 
   // observed data
   int prevBins = 0; 
   RooRealVar * xobs = 0;
   if (GetNullModel()->GetObservables() && GetNullModel()->GetObservables()->getSize() == 1 ) {
      xobs = (RooRealVar*) (GetNullModel()->GetObservables())->first();
      if (data.IsA() == RooDataHist::Class() ) { 
         if (data.numEntries() != xobs->getBins() ) { 
            prevBins = xobs->getBins();
            oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator: number of bins in " << xobs->GetName() << " are different than data bins " 
                                                << " set the same data bins " << data.numEntries() << " in range " 
                                                << " [ " << xobs->getMin() << " , " << xobs->getMax() << " ]" << std::endl;
            xobs->setBins(data.numEntries());
         }
      }
   }

   if (!fNominalAsimov) {
      if (verbose >= 0) 
         oocoutI((TObject*)0,InputArguments) << "AsymptoticCalculator: Asimov data will be generated using fitted nuisance parameter values" << endl;
      RooArgSet * tmp = (RooArgSet*) poiAlt.snapshot(); 
      fAsimovData = MakeAsimovData( data, *GetNullModel(), poiAlt, fAsimovGlobObs,tmp);
   }

   else {
      // assume use current value of nuisance as nominal ones
      if (verbose >= 0) 
         oocoutI((TObject*)0,InputArguments) << "AsymptoticCalculator: Asimovdata set will be generated using nominal (current) nuisance parameter values" << endl;
      nominalParams = poiAlt; // set poi to alt value but keep nuisance at the nominal one
      fAsimovData = MakeAsimovData( *GetNullModel(), nominalParams, fAsimovGlobObs);
   }

   if (!fAsimovData) { 
      oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator: Error : Asimov data set could not be generated " << endl;
      return false;
   }

   // set global observables to their Asimov values 
   RooArgSet globObs;
   RooArgSet globObsSnapshot;
   if (GetNullModel()->GetGlobalObservables()  ) {
      globObs.add(*GetNullModel()->GetGlobalObservables());
      assert(globObs.getSize() == fAsimovGlobObs.getSize() );
      // store previous snapshot value
      globObs.snapshot(globObsSnapshot);
      globObs = fAsimovGlobObs; 
   }


   // evaluate  the likelihood. Since we use on Asimov data , conditional and unconditional values should be the same
   // do conditional fit since is faster

   RooRealVar * muAlt = (RooRealVar*) poiAlt.first();
   assert(muAlt);
   if (verbose>=0)
      oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::Initialize Find  best conditional NLL on ASIMOV data set for given alt POI ( " << 
         muAlt->GetName() << " ) = " << muAlt->getVal() << std::endl;

   fNLLAsimov =  EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), &poiAlt );
   // for unconditional fit 
   //fNLLAsimov =  EvaluateNLL( *nullPdf, *fAsimovData);
   //poi->Print("v");

   // restore previous value 
   globObs = globObsSnapshot;

   // restore number of bins 
   if (prevBins > 0 && xobs) xobs->setBins(prevBins);

   fIsInitialized = true; 
   return true; 
}

//_________________________________________________________________
Double_t AsymptoticCalculator::EvaluateNLL(RooAbsPdf & pdf, RooAbsData& data,   const RooArgSet * condObs, const RooArgSet *poiSet) {

    int verbose = fgPrintLevel;
      
    
    RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
    if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);


    RooArgSet* allParams = pdf.getParameters(data);
    RooStats::RemoveConstantParameters(allParams);
    // add constraint terms for all non-constant parameters

    RooArgSet conditionalObs;
    if (condObs) conditionalObs.add(*condObs);

    // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
    RooAbsReal* nll = pdf.createNLL(data, RooFit::CloneData(kFALSE),RooFit::Constrain(*allParams),RooFit::ConditionalObservables(conditionalObs), RooFit::Offset(RooStats::IsNLLOffset()));

    RooArgSet* attachedSet = nll->getVariables();

    // if poi are specified - do a conditional fit 
    RooArgSet paramsSetConstant;
    // support now only one POI 
    if (poiSet && poiSet->getSize() > 0) { 
       RooRealVar * muTest = (RooRealVar*) (poiSet->first());
       RooRealVar * poiVar = dynamic_cast<RooRealVar*>( attachedSet->find( muTest->GetName() ) );
       if (poiVar && !poiVar->isConstant() ) {
          poiVar->setVal(  muTest->getVal() );
          poiVar->setConstant(); 
          paramsSetConstant.add(*poiVar);
       }
       if (poiSet->getSize() > 1) 
          std::cout << "Model with more than one POI are not supported - ignore extra parameters, consider only first one" << std::endl;

 

       // This for more than one POI (not yet supported)
       //
       // RooLinkedListIter it = poiSet->iterator();
       // RooRealVar* tmpPar = NULL, *tmpParA=NULL;
       // while((tmpPar = (RooRealVar*)it.Next())){
       //    tmpParA =  ((RooRealVar*)attachedSet->find(tmpPar->GetName()));
       //    tmpParA->setVal( tmpPar->getVal() );
       //    if (!tmpParA->isConstant() ) { 
       //       tmpParA->setConstant();
       //       paramsSetConstant.add(*tmpParA);
       //    }
       // }

       // check if there are non-const parameters so it is worth to do the minimization

    }

    TStopwatch tw; 
    tw.Start();
    double val =  -1;

    //check if needed to skip the fit 
    RooArgSet nllParams(*attachedSet); 
    RooStats::RemoveConstantParameters(&nllParams);
    delete attachedSet;
    bool skipFit = (nllParams.getSize() == 0);

    if (skipFit) 
       val = nll->getVal(); // just evaluate nll in conditional fits with model without nuisance params
    else {

       
       int minimPrintLevel = verbose;
       
       RooMinimizer minim(*nll);
       int strategy = ROOT::Math::MinimizerOptions::DefaultStrategy();
       minim.setStrategy( strategy);
       // use tolerance - but never smaller than 1 (default in RooMinimizer)
       double tol =  ROOT::Math::MinimizerOptions::DefaultTolerance();
       tol = std::max(tol,1.0); // 1.0 is the minimum value used in RooMinimizer
       minim.setEps( tol );
       //LM: RooMinimizer.setPrintLevel has +1 offset - so subtruct  here -1
       minim.setPrintLevel(minimPrintLevel-1);
       int status = -1;
       minim.optimizeConst(2);
       TString minimizer = ROOT::Math::MinimizerOptions::DefaultMinimizerType(); 
       TString algorithm = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo(); 

       if (verbose > 0 )
          std::cout << "AsymptoticCalculator::EvaluateNLL  ........ using " << minimizer << " / " << algorithm 
                    << " with strategy  " << strategy << " and tolerance " << tol << std::endl;


       for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
          //	 status = minim.minimize(fMinimizer, ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str());
          status = minim.minimize(minimizer, algorithm);  
          if (status%1000 == 0) {  // ignore erros from Improve 
             break;
          } else { 
             if (tries == 1) {
                printf("    ----> Doing a re-scan first\n");
                minim.minimize(minimizer,"Scan");
             }
             if (tries == 2) {
                if (ROOT::Math::MinimizerOptions::DefaultStrategy() == 0 ) { 
                   printf("    ----> trying with strategy = 1\n");
                   minim.setStrategy(1);
                }
                else 
                   tries++; // skip this trial if stratehy is already 1 
             }
             if (tries == 3) {
                printf("    ----> trying with improve\n");
                minimizer = "Minuit";
                algorithm = "migradimproved";
             }
          }
       }
       
       RooFitResult * result = 0; 


       if (status%100 == 0) { // ignore errors in Hesse or in Improve
          result = minim.save();
       }
       if (result){
          if (!RooStats::IsNLLOffset() ) 
             val = result->minNll();
          else {
             bool previous = RooAbsReal::hideOffset();
             RooAbsReal::setHideOffset(kTRUE) ;
             val = nll->getVal();
             if (!previous)  RooAbsReal::setHideOffset(kFALSE) ;
          }
             
       }
       else { 
          oocoutE((TObject*)0,Fitting) << "FIT FAILED !- return a NaN NLL " << std::endl;
          val =  TMath::QuietNaN();       
       }

       minim.optimizeConst(false);
       if (result) delete result;


    }

    double muTest = 0; 
    if (verbose > 0) { 
       std::cout << "AsymptoticCalculator::EvaluateNLL -  value = " << val;
       if (poiSet) { 
          muTest = ( (RooRealVar*) poiSet->first() )->getVal();
          std::cout << " for poi fixed at = " << muTest; 
       }
       if (!skipFit) {
          std::cout << "\tfit time : ";  
          tw.Print();
       }
       else 
          std::cout << std::endl;
    }

    // reset the parameter free which where set as constant
    if (poiSet && paramsSetConstant.getSize() > 0) SetAllConstant(paramsSetConstant,false);


    if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);

    delete allParams;
    delete nll;

    return val;
}

//____________________________________________________
HypoTestResult* AsymptoticCalculator::GetHypoTest() const {
   // It performs an hypothesis tests using the likelihood function
   // and computes the p values for the null and the alternate using the asymptotic 
   // formulae for the profile likelihood ratio.
   // See G. Cowan, K. Cranmer, E. Gross and O. Vitells.
   // Asymptotic formulae for likelihood- based tests of new physics. Eur. Phys. J., C71:1–19, 2011.
   // The formulae are valid only for one POI. If more than one POI exists consider as POI only the 
   // first one

   int verbose = fgPrintLevel;

   // re-initialized the calculator in case it is needed (pdf or data modifided)
   if (!fIsInitialized) {
      if (!Initialize() ) {
         oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest - Error initializing Asymptotic calculator - return NULL result " << endl;
         return 0;
      }
   }

   if (!fAsimovData) { 
       oocoutE((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest - Asimov data set has not been generated - return NULL result " << endl;
       return 0;
   }

   assert(GetNullModel() );
   assert(GetData() );

   RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
   assert(nullPdf); 

   // make conditional fit on null snapshot of poi

   const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
   assert(nullSnapshot && nullSnapshot->getSize() > 0);
   
   // use as POI the nullSnapshot
   // if more than one POI exists, consider only the first one
   RooArgSet poiTest(*nullSnapshot);

   if (poiTest.getSize() > 1)  { 
      oocoutW((TObject*)0,InputArguments) << "AsymptoticCalculator::GetHypoTest: snapshot has more than one POI - assume as POI first parameter " << std::endl;         
   }

   RooArgSet * allParams = nullPdf->getParameters(*GetData() );
   *allParams = fBestFitParams;
   delete allParams;

   // set the one-side condition
   // (this works when we have only one params of interest 
   RooRealVar * muHat =  dynamic_cast<RooRealVar*> (  fBestFitPoi.first() );
   assert(muHat && "no best fit parameter defined"); 
   RooRealVar * muTest = dynamic_cast<RooRealVar*> ( nullSnapshot->find(muHat->GetName() ) );
   assert(muTest && "poi snapshot is not existing"); 



   if (verbose> 0) {
      std::cout << std::endl;
      oocoutI((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest: - perform  an hypothesis test for  POI ( " << muTest->GetName() << " ) = " << muTest->getVal() << std::endl;
      oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest -  Find  best conditional NLL on OBSERVED data set ..... " << std::endl;
   }

   // evaluate the conditional NLL on the observed data for the snapshot value
   double condNLL = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()), GetNullModel()->GetConditionalObservables(), &poiTest);

   double qmu = 2.*(condNLL - fNLLObs); 
   
   

   if (verbose > 0) 
      oocoutP((TObject*)0,Eval) << "\t OBSERVED DATA :  qmu   = " << qmu << " condNLL = " << condNLL << " uncond " << fNLLObs << std::endl;


   // this tolerance is used to avoid having negative qmu due to numerical errors
   double tol = 1.E-4 * ROOT::Math::MinimizerOptions::DefaultTolerance();
   if (qmu < -tol || TMath::IsNaN(fNLLObs) ) {

      if (qmu < 0) 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  Found a negative value of the qmu - retry to do the unconditional fit " 
                                           << std::endl;         
      else 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  unconditional fit failed before - retry to do it now " 
                                           << std::endl;         
      
      
      double nll = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()),GetNullModel()->GetConditionalObservables());
      
      if (nll < fNLLObs || (TMath::IsNaN(fNLLObs) && !TMath::IsNaN(nll) ) ) { 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  Found a better unconditional minimum "
                                           << " old NLL = " << fNLLObs << " old muHat " << muHat->getVal() << std::endl;            

         // update values
         fNLLObs = nll; 
         const RooArgSet * poi = GetNullModel()->GetParametersOfInterest(); 
         assert(poi);
         fBestFitPoi.removeAll();
         poi->snapshot(fBestFitPoi);
         // restore also muHad since previous pointr has been deleted
         muHat =  dynamic_cast<RooRealVar*> (  fBestFitPoi.first() );
         assert(muHat);

        oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  New minimum  found for                       "
                                          << "    NLL = " << fNLLObs << "    muHat  " << muHat->getVal() << std::endl;            


        qmu = 2.*(condNLL - fNLLObs); 

        if (verbose > 0) 
           oocoutP((TObject*)0,Eval) << "After unconditional refit,  new qmu value is " << qmu << std::endl;

      }
   }

   if (qmu < -tol ) {       
      oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator:  qmu is still < 0  for mu = " 
                                        <<  muTest->getVal() << " return a dummy result "  
                                        << std::endl;         
      return new HypoTestResult();
   }
   if (TMath::IsNaN(qmu) ) {       
      oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator:  failure in fitting for qmu or qmuA " 
                                        <<  muTest->getVal() << " return a dummy result "  
                                        << std::endl;         
      return new HypoTestResult();
   }





   // compute conditional ML on Asimov data set
   // (need to const cast because it uses fitTo which is a non const method
   // RooArgSet asimovGlobObs;
   // RooAbsData * asimovData = (const_cast<AsymptoticCalculator*>(this))->MakeAsimovData( poi, asimovGlobObs);
   // set global observables to their Asimov values 
   RooArgSet globObs;
   RooArgSet globObsSnapshot;
   if (GetNullModel()->GetGlobalObservables()  ) {
      globObs.add(*GetNullModel()->GetGlobalObservables());
      // store previous snapshot value
      globObs.snapshot(globObsSnapshot);
      globObs = fAsimovGlobObs; 
   }


   if (verbose > 0) oocoutP((TObject*)0,Eval) << "AsymptoticCalculator::GetHypoTest -- Find  best conditional NLL on ASIMOV data set .... " << std::endl;

   double condNLL_A = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), &poiTest);


   double qmu_A = 2.*(condNLL_A - fNLLAsimov  );

   if (verbose > 0) 
      oocoutP((TObject*)0,Eval) << "\t ASIMOV data qmu_A = " << qmu_A << " condNLL = " << condNLL_A << " uncond " << fNLLAsimov << std::endl;

   if (qmu_A < -tol || TMath::IsNaN(fNLLAsimov) ) {

      if (qmu_A < 0) 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  Found a negative value of the qmu Asimov- retry to do the unconditional fit " 
                                           << std::endl;         
      else 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  Fit failed for  unconditional the qmu Asimov- retry  unconditional fit " 
                                           << std::endl;         
      
      
      double nll = EvaluateNLL( *nullPdf, *fAsimovData,  GetNullModel()->GetConditionalObservables() );
      
      if (nll < fNLLAsimov || (TMath::IsNaN(fNLLAsimov) && !TMath::IsNaN(nll) )) { 
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  Found a better unconditional minimum for Asimov data set"
                                           << " old NLL = " << fNLLAsimov << std::endl;            

         // update values
         fNLLAsimov = nll; 
         
         oocoutW((TObject*)0,Minimization) << "AsymptoticCalculator:  New minimum  found for                       "
                                           << "    NLL = " << fNLLAsimov << std::endl;            
         qmu_A = 2.*(condNLL_A - fNLLAsimov); 

        if (verbose > 0) 
           oocoutP((TObject*)0,Eval) << "After unconditional Asimov refit,  new qmu_A value is " << qmu_A << std::endl;

      }
   }

   if (qmu_A < - tol) {       
      oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator:  qmu_A is still < 0  for mu = " 
                                        <<  muTest->getVal() << " return a dummy result "  
                                        << std::endl;         
      return new HypoTestResult();
   }
   if (TMath::IsNaN(qmu) ) {       
      oocoutE((TObject*)0,Minimization) << "AsymptoticCalculator:  failure in fitting for qmu or qmuA " 
                                        <<  muTest->getVal() << " return a dummy result "  
                                        << std::endl;         
      return new HypoTestResult();
   }


   // restore previous value of global observables
   globObs = globObsSnapshot;

   // now we compute p-values using the asumptotic formulae 
   // described in the paper 
   //  Cowan et al, Eur.Phys.J. C (2011) 71:1554

   // first try to guess autoatically if needed to use qtilde (or ttilde in case of two sided) 
   // if explicitly fUseQTilde this was not set
   // qtilde is in this case used if poi is bounded at the value of the alt hypothesis
   //  for Qtilde (need to distinguish case when qmu > qmuA = mu^2/ sigma^2) 
   // (see Cowan et al, Eur.Phys.J. C(2011) 71:1554 paper equations 64 and 65
   // (remember qmu_A = mu^2/sigma^2 )
   bool useQTilde = false; 
   // default case (check if poi is limited or not to a zero value)
   if (!fOneSidedDiscovery) { // qtilde is not a discovery test 
      if (fUseQTilde == -1 && !fOneSidedDiscovery) { 
         // alternate snapshot is value for which background is zero (for limits)
         RooRealVar * muAlt = dynamic_cast<RooRealVar*>( GetAlternateModel()->GetSnapshot()->first() );
         // null snapshot is value for which background is zero (for discovery)
         //RooRealVar * muNull = dynamic_cast<RooRealVar*>( GetNullModel()->GetSnapshot()->first() );
         assert(muAlt != 0 );
         if (muTest->getMin() == muAlt->getVal()   ) { 
            fUseQTilde = 1;  
            oocoutI((TObject*)0,InputArguments) << "Minimum of POI is " << muTest->getMin() << " corresponds to alt  snapshot   - using qtilde asymptotic formulae  " << std::endl;
         } else {
            fUseQTilde = 0;  
            oocoutI((TObject*)0,InputArguments) << "Minimum of POI is " << muTest->getMin() << " is different to alt snapshot " << muAlt->getVal() 
                                                << " - using standard q asymptotic formulae  " << std::endl;
         }         
      }
      useQTilde = fUseQTilde;
   }


   //check for one side condition (remember this is valid only for one poi)
   if (fOneSided ) { 
      if ( muHat->getVal() > muTest->getVal() ) { 
         oocoutI((TObject*)0,Eval) << "Using one-sided qmu - setting qmu to zero  muHat = " << muHat->getVal() 
                                   << " muTest = " << muTest->getVal() << std::endl;
         qmu = 0;
      }
   }
   if (fOneSidedDiscovery ) { 
      if ( muHat->getVal() < muTest->getVal() ) { 
         oocoutI((TObject*)0,Eval) << "Using one-sided discovery qmu - setting qmu to zero  muHat = " << muHat->getVal() 
                                   << " muTest = " << muTest->getVal() << std::endl;
         qmu = 0;
      }
   }

   // fix for negative qmu values due to numerical errors
   if (qmu < 0 && qmu > -tol) qmu = 0; 
   if (qmu_A < 0 && qmu_A > -tol) qmu_A = 0; 

   // asymptotic formula for pnull and from  paper Eur.Phys.J C 2011  71:1554
   // we have 4 different cases: 
   //          t(mu), t_tilde(mu) for the 2-sided 
   //          q(mu) and q_tilde(mu) for the one -sided test statistics

   double pnull = -1, palt = -1;

   // asymptotic formula for pnull (for only one POI) 
   // From fact that qmu is a chi2 with ndf=1

   double sqrtqmu = (qmu > 0) ? std::sqrt(qmu) : 0; 
   double sqrtqmu_A = (qmu_A > 0) ? std::sqrt(qmu_A) : 0; 

   
   if (fOneSided || fOneSidedDiscovery) {
      // for one-sided PL (q_mu : equations 56,57)
      if (verbose>2) {
         if (fOneSided) 
            oocoutI((TObject*)0,Eval) << "Using one-sided limit asymptotic formula (qmu)" << endl;
         else
            oocoutI((TObject*)0,Eval) << "Using one-sided discovery asymptotic formula (q0)" << endl;
      }
      pnull = ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
      palt = ROOT::Math::normal_cdf( sqrtqmu_A - sqrtqmu, 1.);
   }
   else  {
      // for 2-sided PL (t_mu : equations 35,36 in asymptotic paper) 
      if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using two-sided asimptotic  formula (tmu)" << endl;
      pnull = 2.*ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
      palt = ROOT::Math::normal_cdf_c( sqrtqmu + sqrtqmu_A, 1.) + 
         ROOT::Math::normal_cdf_c( sqrtqmu - sqrtqmu_A, 1.); 
         
   }

   if (useQTilde ) { 
      if (fOneSided) { 
         // for bounded one-sided (q_mu_tilde: equations 64,65)
         if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) { // to avoid case 0/0
            if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using qmu_tilde (qmu is greater than qmu_A)" << endl;
            pnull = ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
            palt = ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
         }
      }
      else {  
         // for 2 sided bounded test statistic  (N.B there is no one sided discovery qtilde)
         // t_mu_tilde: equations 43,44 in asymptotic paper
         if ( qmu >  qmu_A  && (qmu_A > 0 || qmu > tol)  ) { 
            if (verbose > 2) oocoutI((TObject*)0,Eval) << "Using tmu_tilde (qmu is greater than qmu_A)" << endl;
            pnull = ROOT::Math::normal_cdf_c(sqrtqmu,1.) + 
                    ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
            palt = ROOT::Math::normal_cdf_c( sqrtqmu_A + sqrtqmu, 1.) + 
                   ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
         }
      }
   }



   // create an HypoTest result but where the sampling distributions are set to zero
   string resultname = "HypoTestAsymptotic_result";
   HypoTestResult* res = new HypoTestResult(resultname.c_str(), pnull, palt);

   if (verbose > 0) 
      //std::cout 
      oocoutP((TObject*)0,Eval) 
         << "poi = " << muTest->getVal() << " qmu = " << qmu << " qmu_A = " << qmu_A 
         << " sigma = " << muTest->getVal()/sqrtqmu_A
         << "  CLsplusb = " << pnull << " CLb = " << palt << " CLs = " <<  res->CLs() << std::endl; 

   return res; 

}

struct PaltFunction { 
   PaltFunction( double offset, double pval, int icase) : 
      fOffset(offset), fPval(pval), fCase(icase) {}
   double operator() (double x) const { 
      return ROOT::Math::normal_cdf_c(x + fOffset) + ROOT::Math::normal_cdf_c(fCase*(x - fOffset)) - fPval;
   }
   double fOffset;
   double fPval;
   int fCase;
};

double AsymptoticCalculator::GetExpectedPValues(double pnull, double palt, double nsigma, bool useCls, bool oneSided ) { 
   // function given the null and the alt p value - return the expected one given the N - sigma value
   if (oneSided) { 
      double sqrtqmu =  ROOT::Math::normal_quantile_c( pnull,1.);
      double sqrtqmu_A =  ROOT::Math::normal_quantile( palt,1.) + sqrtqmu;
      double clsplusb = ROOT::Math::normal_cdf_c( sqrtqmu_A - nsigma, 1.);
      if (!useCls) return clsplusb; 
      double clb = ROOT::Math::normal_cdf( nsigma, 1.);
      return (clb == 0) ? -1 : clsplusb / clb;  
   }

   // case of 2 sided test statistic
   // need to compute numerically
   double sqrttmu =  ROOT::Math::normal_quantile_c( 0.5*pnull,1.);
   if (sqrttmu == 0) { 
      // here cannot invert the function - skip the point
      return -1; 
   }      
   // invert formula for palt to get sqrttmu_A
   PaltFunction f( sqrttmu, palt, -1);       
   ROOT::Math::BrentRootFinder brf;
   ROOT::Math::WrappedFunction<PaltFunction> wf(f);
   brf.SetFunction( wf, 0, 20);
   bool ret = brf.Solve();
   if (!ret) { 
      oocoutE((TObject*)0,Eval)  << "Error finding expected p-values - return -1" << std::endl;
      return -1;
   }
   double sqrttmu_A = brf.Root();

   // now invert for expected value
   PaltFunction f2( sqrttmu_A,  ROOT::Math::normal_cdf( nsigma, 1.), 1);
   ROOT::Math::WrappedFunction<PaltFunction> wf2(f2);
   brf.SetFunction(wf2,0,20);
   ret = brf.Solve();
   if (!ret) { 
      oocoutE((TObject*)0,Eval)  << "Error finding expected p-values - return -1" << std::endl;
      return -1;
   }   
   return  2*ROOT::Math::normal_cdf_c( brf.Root(),1.);
}

// void GetExpectedLimit(double nsigma, double alpha, double &clsblimit, double &clslimit) { 
//    // get expected limit 
//    double 
// }


void AsymptoticCalculator::FillBins(const RooAbsPdf & pdf, const RooArgList &obs, RooAbsData & data, int &index,  double &binVolume, int &ibin) { 
   /// fill bins by looping recursivly on observables 

   bool debug = (fgPrintLevel >= 2);  

   RooRealVar * v = dynamic_cast<RooRealVar*>( &(obs[index]) );
   if (!v) return;

   RooArgSet obstmp(obs);
   double expectedEvents = pdf.expectedEvents(obstmp);
   // if (debug)  { 
   //    std::cout << "expected events = " << expectedEvents << std::endl;
   // }

   if (debug) cout << "looping on observable " << v->GetName() << endl;
   for (int i = 0; i < v->getBins(); ++i) {
      v->setBin(i);
      if (index < obs.getSize() -1) {
         index++;  // increase index
         double prevBinVolume = binVolume; 
         binVolume *= v->getBinWidth(i); // increase bin volume
         FillBins(pdf, obs, data, index,  binVolume, ibin);
         index--; // decrease index
         binVolume = prevBinVolume; // decrease also bin volume
      }
      else {

         // this is now a new bin - compute the pdf in this bin 
         double totBinVolume = binVolume * v->getBinWidth(i);
         double fval = pdf.getVal(&obstmp)*totBinVolume;

         //if (debug) std::cout << "pdf value in the bin " << fval << " bin volume = " << totBinVolume << "   " << fval*expectedEvents << std::endl;
         if (fval*expectedEvents <= 0)
         {
            if (fval*expectedEvents < 0)
               cout << "WARNING::Detected a bin with negative expected events! Please check your inputs." << endl;
            else
               cout << "WARNING::Detected a bin with zero expected events- skip it" << endl;
         }
         // have a cut off for overflows ??
         else 
            data.add(obs, fval*expectedEvents);

         if (debug) { 
            cout << "bin " << ibin << "\t";
            for (int j=0; j < obs.getSize(); ++j) { cout << "  " <<  ((RooRealVar&) obs[j]).getVal(); }
            cout << " w = " << fval*expectedEvents;
            cout << endl;
         }
         // RooArgSet xxx(obs);
         // h3->Fill(((RooRealVar&) obs[0]).getVal(), ((RooRealVar&) obs[1]).getVal(), ((RooRealVar&) obs[2]).getVal() ,
         //          pdf->getVal(&xxx) );
         ibin++;
      }
   }
   //reset bin values
   if (debug) 
      cout << "ending loop on .. " << v->GetName() << endl;

   v->setBin(0);
   
}


bool AsymptoticCalculator::SetObsToExpected(RooProdPdf &prod, const RooArgSet &obs) 
{
   // iterate a Prod pdf to find all the Poisson or Gaussian part to set the observed value to expected one
    RooLinkedListIter  iter(prod.pdfList().iterator());
    bool ret = false;
    for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
        if (!a->dependsOn(obs)) continue;
        RooPoisson *pois = 0;
        RooGaussian * gaus = 0;
        if ((pois = dynamic_cast<RooPoisson *>(a)) != 0) {
            SetObsToExpected(*pois, obs);
            pois->setNoRounding(true);  //needed since expecteed value is not an integer
        } else if ((gaus = dynamic_cast<RooGaussian *>(a)) != 0) {
            SetObsToExpected(*gaus, obs);
        } else {
           // should try to add also lognormal case ? 
            RooProdPdf *subprod = dynamic_cast<RooProdPdf *>(a);
            if (subprod) 
               return SetObsToExpected(*subprod, obs);
            else {
               oocoutE((TObject*)0,InputArguments) << "Illegal term in counting model: depends on observables, but not Poisson or Gaussian or Product" 
                                                   << endl;
               return false;
            }
        }
        ret = (pois != 0 || gaus != 0 ); 
    }
    return ret;
}

bool AsymptoticCalculator::SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs) 
{
   // set observed value to the expected one 
   // works for Gaussian, Poisson or LogNormal
   // assumes mean parameter value is the argument not constant and not depoending on observables
   // (if more than two arguments are not constant will use first one but printr a warning !)
   // need to iterate on the components of the POisson to get n and nu (nu can be a RooAbsReal)
   // (code from G. Petrucciani and extended by L.M.)
   RooRealVar *myobs = 0;
   RooAbsReal *myexp = 0;
   const char * pdfName = pdf.IsA()->GetName();
   RooFIter iter(pdf.serverMIterator());
   for (RooAbsArg *a =  iter.next(); a != 0; a = iter.next()) {
      if (obs.contains(*a)) {
         if (myobs != 0) { 
            oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two observables ?? " << endl;
            return false;
         }
         myobs = dynamic_cast<RooRealVar *>(a);
         if (myobs == 0) { 
            oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Observable is not a RooRealVar??" << endl;
            return false; 
         }
      } else {
         if (!a->isConstant() ) {
            if (myexp != 0) { 
               oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Has two non-const arguments  " << endl;
               return false;
            }
            myexp = dynamic_cast<RooAbsReal *>(a);
            if (myexp == 0) { 
               oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : Expected is not a RooAbsReal??" << endl;
               return false; 
            }
         }
      }
   }
   if (myobs == 0)  { 
      oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
      return false;
   }
   if (myexp == 0) {    
      oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::SetObsExpected( " << pdfName << " ) : No observable?" << endl;
      return false;
   }

   myobs->setVal(myexp->getVal());

   if (fgPrintLevel > 2) { 
      std::cout << "SetObsToExpected : setting " << myobs->GetName() << " to expected value " << myexp->getVal() << " of " << myexp->GetName() << std::endl;
   }

   return true;
}


RooAbsData * AsymptoticCalculator::GenerateCountingAsimovData(RooAbsPdf & pdf, const RooArgSet & observables,  const RooRealVar & , RooCategory * channelCat) { 
   // generate counting Asimov data for the case when the pdf cannot be extended
   // assume pdf is a RooPoisson or can be decomposed in a product of RooPoisson, 
   // otherwise we cannot know how to make the Asimov data sets in the other cases
    RooArgSet obs(observables);
    RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
    RooPoisson *pois = 0;
    RooGaussian *gaus = 0;

    if (fgPrintLevel > 1) 
       std::cout << "generate counting Asimov data for pdf of type " << pdf.IsA()->GetName() << std::endl;

    bool r = false;
    if (prod != 0) {
        r = SetObsToExpected(*prod, observables);
    } else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != 0) {
        r = SetObsToExpected(*pois, observables);
        // we need in this case to set Poisson to real values
        pois->setNoRounding(true);
    } else if ((gaus = dynamic_cast<RooGaussian *>(&pdf)) != 0) {
        r = SetObsToExpected(*gaus, observables);
    } else {
       oocoutE((TObject*)0,InputArguments) << "A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << endl;
    }
    if (!r) return 0;
    int icat = 0;
    if (channelCat) {
       icat = channelCat->getIndex(); 
    }

    RooDataSet *ret = new RooDataSet(TString::Format("CountingAsimovData%d",icat),TString::Format("CountingAsimovData%d",icat), obs);
    ret->add(obs);
    return ret;
}

RooAbsData * AsymptoticCalculator::GenerateAsimovDataSinglePdf(const RooAbsPdf & pdf, const RooArgSet & allobs,  const RooRealVar & weightVar, RooCategory * channelCat) { 
   // compute the asimov data set for an observable of a pdf 
   // use the number of bins sets in the observables 
   // to do :  (possibility to change number of bins)
   // implement integration over bin content

   int printLevel = fgPrintLevel;

   // Get observables defined by the pdf associated with this state
   std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );


   // if pdf cannot be extended assume is then a counting experiment
   if (!pdf.canBeExtended() ) return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);

   RooArgSet obsAndWeight(*obs); 
   obsAndWeight.add(weightVar);

   RooDataSet* asimovData = 0; 
   if (channelCat) {
      int icat = channelCat->getIndex(); 
      asimovData = new RooDataSet(TString::Format("AsimovData%d",icat),TString::Format("combAsimovData%d",icat),
                                  RooArgSet(obsAndWeight,*channelCat),RooFit::WeightVar(weightVar));
   }
   else 
      asimovData = new RooDataSet("AsimovData","AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));

    // This works ony for 1D observables 
    //RooRealVar* thisObs = ((RooRealVar*)obstmp->first());

    RooArgList obsList(*obs);

    // loop on observables and on the bins 
    if (printLevel >= 2) { 
       cout << "Generating Asimov data for pdf " << pdf.GetName() << endl;
       cout << "list of observables  " << endl;
       obsList.Print();
    }

    int obsIndex = 0; 
    double binVolume = 1; 
    int nbins = 0; 
    FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
    if (printLevel >= 2) 
       cout << "filled from " << pdf.GetName() << "   " << nbins << " nbins " << " volume is " << binVolume << endl;

    // for (int iobs = 0; iobs < obsList.getSize(); ++iobs) { 
    //    RooRealVar * thisObs = dynamic_cast<RooRealVar*> &obsList[i];
    //    if (thisObs == 0) continue; 
    //    // loop on the bin contents
    //    for(int  ibin=0; ibin<thisObs->numBins(); ++ibin){
    //       thisObs->setBin(ibin);

    //   thisNorm=pdftmp->getVal(obstmp)*thisObs->getBinWidth(jj);
    //   if (thisNorm*expectedEvents <= 0)
    //   {
    //     cout << "WARNING::Detected bin with zero expected events! Please check your inputs." << endl;
    //   }
    //   // have a cut off for overflows ??
    //   obsDataUnbinned->add(*mc->GetObservables(), thisNorm*expectedEvents);
    // }
    
    if (printLevel >= 1)
    {
      asimovData->Print();
      //cout <<"sum entries "<< asimovData->sumEntries()<<endl;
    }
    if( TMath::IsNaN(asimovData->sumEntries()) ){
      cout << "sum entries is nan"<<endl;
      assert(0);
      delete asimovData;
      asimovData = 0;
    }

    return asimovData;

}

RooAbsData * AsymptoticCalculator::GenerateAsimovData(const RooAbsPdf & pdf, const RooArgSet & observables  )  { 
   // generate the asimov data for the observables (not the global ones) 
   // need to deal with the case of a sim pdf 

   int printLevel = fgPrintLevel;

   RooRealVar * weightVar = new RooRealVar("binWeightAsimov", "binWeightAsimov", 1, 0, 1.E30 );

   if (printLevel > 1) cout <<" Generate Asimov data for observables"<<endl;
  //RooDataSet* simData=NULL;
   const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf);
   if (!simPdf) { 
      // generate data for non sim pdf
      return GenerateAsimovDataSinglePdf( pdf, observables, *weightVar, 0);
   }

   std::map<std::string, RooDataSet*> asimovDataMap;
    
  //look at category of simpdf 
  RooCategory& channelCat = (RooCategory&)simPdf->indexCat();
  int nrIndices = channelCat.numTypes();
  if( nrIndices == 0 ) {
    oocoutW((TObject*)0,Generation) << "Simultaneous pdf does not contain any categories." << endl;
  }
  for (int i=0;i<nrIndices;i++){
    channelCat.setIndex(i);
    //iFrame++;
    // Get pdf associated with state from simpdf
    RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getLabel()) ;
    assert(pdftmp != 0);
	
    if (printLevel > 1)
    {
      cout << "on type " << channelCat.getLabel() << " " << channelCat.getIndex() << endl;
    }

    RooAbsData * dataSinglePdf = GenerateAsimovDataSinglePdf( *pdftmp, observables, *weightVar, &channelCat);
    //((RooRealVar*)obstmp->first())->Print();
    //cout << "expected events " << pdftmp->expectedEvents(*obstmp) << endl;
    if (!dataSinglePdf) { 
       oocoutE((TObject*)0,Generation) << "Error generating an Asimov data set for pdf " << pdftmp->GetName() << endl;
       return 0;
    }
     

    asimovDataMap[string(channelCat.getLabel())] = (RooDataSet*) dataSinglePdf;

    if (printLevel > 1)
    {
      cout << "channel: " << channelCat.getLabel() << ", data: ";
      dataSinglePdf->Print();
      cout << endl;
    }
  }

  RooArgSet obsAndWeight(observables); 
  obsAndWeight.add(*weightVar);


  RooDataSet* asimovData = new RooDataSet("asimovDataFullModel","asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
                                          RooFit::Index(channelCat),RooFit::Import(asimovDataMap),RooFit::WeightVar(*weightVar));

  delete weightVar; 
  return asimovData;

}

//______________________________________________________________________________
RooAbsData * AsymptoticCalculator::MakeAsimovData(RooAbsData & realData, const ModelConfig & model, const  RooArgSet & paramValues, RooArgSet & asimovGlobObs, const RooArgSet * genPoiValues )  {
   // static function to the an Asimov data set
   // given an observed dat set, a model and a snapshot of poi. 
   // Return the asimov data set + global observables set to values satisfying the constraints


   int verbose = fgPrintLevel;


   RooArgSet  poi(*model.GetParametersOfInterest());
   poi = paramValues; 
   //RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
   // set poi constant for conditional MLE 
   // need to fit nuisance parameters at their conditional MLE value
   RooLinkedListIter it = poi.iterator();
   RooRealVar*  tmpPar = NULL;
   RooArgSet paramsSetConstant; 
   while((tmpPar = (RooRealVar*)it.Next())){
      tmpPar->setConstant();
      if (verbose>0)  
         std::cout << "MakeAsimov: Setting poi " << tmpPar->GetName() << " to a constant value = " << tmpPar->getVal() << std::endl;
      paramsSetConstant.add(*tmpPar); 
   }

   // find conditional value of the nuisance parameters
   bool hasFloatParams = false;
   RooArgSet  constrainParams;
   if (model.GetNuisanceParameters()) {
      constrainParams.add(*model.GetNuisanceParameters());
      RooStats::RemoveConstantParameters(&constrainParams);
      if (constrainParams.getSize() > 0) hasFloatParams = true; 

   } else {
      // Do we have free parameters anyway that need fitting?
      std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
      RooLinkedListIter iter(params->iterator());
      for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
         RooRealVar *rrv = dynamic_cast<RooRealVar *>(a);
         if ( rrv != 0 && rrv->isConstant() == false ) { hasFloatParams = true; break; }
      } 
   }
   if (hasFloatParams) { 
      // models need to be fitted to find best nuisance parameter values

      TStopwatch tw2; tw2.Start(); 
      int minimPrintLevel = ROOT::Math::MinimizerOptions::DefaultPrintLevel();
      if (verbose>0) { 
         std::cout << "MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
         minimPrintLevel = verbose;
         if (verbose>1) {
            std::cout << "POI values:\n"; poi.Print("v");
            if (verbose > 2) { 
               std::cout << "Nuis param values:\n"; 
               constrainParams.Print("v");
            }
         }
      }         
      RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
      if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);

      RooArgSet conditionalObs;
      if (model.GetConditionalObservables()) conditionalObs.add(*model.GetConditionalObservables());

      std::string minimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
      std::string minimizerAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo();
      model.GetPdf()->fitTo(realData, RooFit::Minimizer(minimizerType.c_str(),minimizerAlgo.c_str()), 
                 RooFit::Strategy(ROOT::Math::MinimizerOptions::DefaultStrategy()),
                 RooFit::PrintLevel(minimPrintLevel-1), RooFit::Hesse(false),
                            RooFit::Constrain(constrainParams),RooFit::ConditionalObservables(conditionalObs), RooFit::Offset(RooStats::IsNLLOffset()));
      if (verbose>0) { std::cout << "fit time "; tw2.Print();}
      if (verbose > 1) { 
         // after the fit the nuisance parameters will have their best fit value
         if (model.GetNuisanceParameters() ) {
            std::cout << "Nuisance parameters after fit for asimov dataset: " << std::endl;
            model.GetNuisanceParameters()->Print("V");
         }
      }

      if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);

   }

   // restore the parameters which were set constant
   SetAllConstant(paramsSetConstant, false);


   RooArgSet *  allParams = model.GetPdf()->getParameters(realData);
   RooStats::RemoveConstantParameters( allParams );

   // if a RooArgSet of poi is passed , different poi will be used for generating the Asimov data set 
   if (genPoiValues) *allParams = *genPoiValues;

   // now do the actual generation of the AsimovData Set
   // no need to pass parameters values since we have set them before
   RooAbsData * asimovData =  MakeAsimovData(model, *allParams, asimovGlobObs);

   delete allParams;

   return asimovData;

}

//______________________________________________________________________________
RooAbsData * AsymptoticCalculator::MakeAsimovData(const ModelConfig & model, const  RooArgSet & allParamValues, RooArgSet & asimovGlobObs)  {
   // static function to the an Asimov data set
   // given the model and the values of all parameters including the nuisance 
   // Return the asimov data set + global observables set to values satisfying the constraints


   int verbose = fgPrintLevel;
 
   TStopwatch tw; 
   tw.Start();

   // set the parameter values (do I need the poi to be constant ? )
   // the nuisance parameter values could be set at their fitted value (the MLE)
   if (allParamValues.getSize() > 0) { 
      RooArgSet *  allVars = model.GetPdf()->getVariables();
      *allVars = allParamValues;
      delete allVars;
   }


   // generate the Asimov data set for the observables 
   RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
   
   if (verbose>0) {
      std::cout << "Generated Asimov data for observables "; (model.GetObservables() )->Print(); 
      if (verbose > 1) { 
         if (asimov->numEntries() == 1 ) {
            std::cout << "--- Asimov data values \n";
            asimov->get()->Print("v");
         }
         else { 
            std::cout << "--- Asimov data numEntries = " << asimov->numEntries() << " sumOfEntries = " << asimov->sumEntries() << std::endl;
         }
         std::cout << "\ttime for generating : ";  tw.Print(); 
      }
   }


    // Now need to have in ASIMOV the data sets also the global observables
   // Their values must be the one satisfying the constraint. 
   // to do it make a nuisance pdf with all product of constraints and then 
   // assign to each constraint a glob observable value = to the current fitted nuisance parameter value
   // IN general  one should solve in general the system of equations f( gobs| nuispar ) = 0 where f are the 
   //  derivatives of the constraint with respect the nuisance parameter and they are evaluated at the best fit nuisance
   // parameter points
   // As simple solution assume that constrain has a direct dependence on the nuisance parameter, i.e. 
   // Constraint (gobs, func( nuispar) ) and the condifion is satisfied for 
   // gobs = func( nuispar) where nunispar is at the MLE value


   if (model.GetGlobalObservables() && model.GetGlobalObservables()->getSize() > 0) {

      if (verbose>1) {
         std::cout << "Generating Asimov data for global observables " << std::endl; 
      }

      RooArgSet gobs(*model.GetGlobalObservables());

      // snapshot data global observables
      RooArgSet snapGlobalObsData;
      SetAllConstant(gobs, true);
      gobs.snapshot(snapGlobalObsData);


      RooArgSet nuis; 
      if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
      if (nuis.getSize() == 0) { 
            oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
                                            << " set global observales to model values " << endl;
            asimovGlobObs = gobs;
            return asimov;
      }

      // part 1: create the nuisance pdf
      std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,"TempNuisPdf") );
      if (nuispdf.get() == 0) { 
         oocoutF((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and global obs but no nuisance pdf "
                                         << std::endl;
      }
      // unfold the nuisance pdf if it is a prod pdf
      RooArgList pdfList;  
      RooProdPdf *prod = dynamic_cast<RooProdPdf *>(nuispdf.get());
      if (prod ) { 
         pdfList.add(prod->pdfList());
      }
      else
         // nothing to unfold - just use the pdf
         pdfList.add(*nuispdf.get());

      RooLinkedListIter iter(pdfList.iterator());
      for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
         RooAbsPdf *cterm = dynamic_cast<RooAbsPdf *>(a); 
         assert(cterm && "AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
         if (!cterm->dependsOn(nuis)) continue; // dummy constraints
         // skip also the case of uniform components
         if (typeid(*cterm) == typeid(RooUniform)) continue;

         std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
         std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
         if (cgobs->getSize() > 1) {
            oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term  " <<  cterm->GetName() 
                                            << " has multiple global observables -cannot generate - skip it" << std::endl;
            continue; 
         }
         else if (cgobs->getSize() == 0) { 
            oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term  " <<  cterm->GetName() 
                                            << " has no global observables - skip it" << std::endl;
            continue; 
         }
         // the variable representing the global observable
         RooRealVar &rrv = dynamic_cast<RooRealVar &>(*cgobs->first());

         // remove the constant parameters in cpars 
         RooStats::RemoveConstantParameters(cpars.get());
         if (cpars->getSize() != 1) {
            oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term " 
                                            << cterm->GetName() << " has multiple floating params - cannot generate - skip it " << std::endl;
            continue;
         }

         bool foundServer = false;
         // note : this will work only for this type of constraints
         // expressed as RooPoisson, RooGaussian, RooLognormal, RooGamma
         TClass * cClass = cterm->IsA();
         if (verbose > 2) std::cout << "Constraint " << cterm->GetName() << " of type " << cClass->GetName() << std::endl;
         if ( cClass != RooGaussian::Class() && cClass != RooPoisson::Class() &&
              cClass != RooGamma::Class() && cClass != RooLognormal::Class() &&
              cClass != RooBifurGauss::Class()  ) {          
            TString className =  (cClass) ?  cClass->GetName() : "undefined"; 
            oocoutW((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term " 
                                            << cterm->GetName() << " of type " << className 
                                            << " is a non-supported type - result might be not correct " << std::endl;
         }

         // in case of a Poisson constraint make sure the rounding is not set 
         if (cClass == RooPoisson::Class() ) { 
            RooPoisson * pois = dynamic_cast<RooPoisson*>(cterm); 
            assert(pois); 
            pois->setNoRounding(true); 
         }

         // look at server of the constraint term and check if the global observable is part of the server
         RooAbsArg * arg = cterm->findServer(rrv); 
         if (!arg) {
            // special case is for the Gamma where one might define the global observable n and you have a Gamma(b, n+1, ...._
            // in this case n+1 is the server and we don;t have a direct dependency, but we want to set n to the b value
            // so in case of the Gamma ignore this test
            if ( cClass != RooGamma::Class() ) {
               oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term " 
                                               << cterm->GetName() << " has no direct dependence on global observable- cannot generate it " << std::endl;
               continue;
            }
         }

         // loop on the server of the constraint term
         // neeed to treat the Gamma as a special case
         // the mode of the Gamma is (k-1)*theta where theta is the inverse of the rate parameter.
         // we assume that the global observable is defined as ngobs = k-1 and the theta parameter has the name theta otherwise we use other procedure which might be wrong
         RooAbsReal * thetaGamma = 0;
         if ( cClass == RooGamma::Class() ) {
            RooFIter itc(cterm->serverMIterator() );
            for (RooAbsArg *a2 = itc.next(); a2 != 0; a2 = itc.next()) {
               if (TString(a2->GetName()).Contains("theta") ) {
                  thetaGamma = dynamic_cast<RooAbsReal*>(a2);
                  break;
               }
            }
            if (thetaGamma == 0) {
               oocoutI((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term " 
                                               << cterm->GetName() << " is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is  1 " << std::endl;
            }
            else {
               if (verbose>2)
                  std::cout << "Gamma constraint has a scale " << thetaGamma->GetName() << "  = " << thetaGamma->getVal() << std::endl;
            }
         }         
         RooFIter iter2(cterm->serverMIterator() );
         for (RooAbsArg *a2 = iter2.next(); a2 != 0; a2 = iter2.next()) {
            RooAbsReal * rrv2 = dynamic_cast<RooAbsReal *>(a2); 
            if (verbose > 2) std::cout << "Loop on constraint server term  " << a2->GetName() << std::endl;
            if (rrv2 && rrv2->dependsOn(nuis) ) { 


               // found server depending on nuisance               
               if (foundServer) { 
                  oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term " 
                                            << cterm->GetName() << " constraint term has more server depending on nuisance- cannot generate it " <<
                     std::endl;
                  foundServer = false;
                  break;
               }
               if (thetaGamma && thetaGamma->getVal() > 0)
                  rrv.setVal( rrv2->getVal() / thetaGamma->getVal() ); 
               else 
                  rrv.setVal( rrv2->getVal() );
               foundServer = true;

               if (verbose>2) 
                  std::cout << "setting global observable "<< rrv.GetName() << " to value " << rrv.getVal() 
                         << " which comes from " << rrv2->GetName() << std::endl;
            }
         }

         if (!foundServer) {  
            oocoutE((TObject*)0,Generation) << "AsymptoticCalculator::MakeAsimovData - can't find nuisance for constraint term - global observales will not be set to Asimov value " << cterm->GetName() << std::endl;
            std::cerr << "Parameters: " << std::endl;
            cpars->Print("V");
            std::cerr << "Observables: " << std::endl;
            cgobs->Print("V");
         }
//         rrv.setVal(match->getVal());
      }

      // make a snapshot of global observables 
      // needed this ?? (LM) 

      asimovGlobObs.removeAll();
      SetAllConstant(gobs, true);
      gobs.snapshot(asimovGlobObs);

      // revert global observables to the data value
      gobs = snapGlobalObsData;

      if (verbose>0) {
         std::cout << "Generated Asimov data for global observables ";
         if (verbose == 1) gobs.Print();  
      }

      if (verbose > 1) {
         std::cout << "\nGlobal observables for data: " << std::endl;
         gobs.Print("V");
         std::cout << "\nGlobal observables for asimov: " << std::endl;
         asimovGlobObs.Print("V");
      }


   }

   return asimov;

}






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