From $ROOTSYS/tutorials/roostats/StandardHypoTestDemo.C

// Standard tutorial macro for hypothesis test (for computing the discovery significance) using all
// RooStats hypotheiss tests calculators and test statistics
//
//
//Author:  L. Moneta
//
// Usage:
//
// root>.L StandardHypoTestDemo.C
// root> StandardHypoTestDemo("fileName","workspace name","S+B modelconfig name","B model name","data set name",calculator type, test statistic type, //                             number of toys)
//
//
// type = 0 Freq calculator
// type = 1 Hybrid calculator
// type = 2 Asymptotic calculator
// type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
//
// testStatType = 0 LEP
//              = 1 Tevatron
//              = 2 Profile Likelihood
//              = 3 Profile Likelihood one sided (i.e. = 0 if mu_hat < 0)
//


#include "TFile.h"
#include "RooWorkspace.h"
#include "RooAbsPdf.h"
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooStats/ModelConfig.h"
#include "RooRandom.h"
#include "TGraphErrors.h"
#include "TGraphAsymmErrors.h"
#include "TCanvas.h"
#include "TLine.h"
#include "TSystem.h"
#include "TROOT.h"

#include "RooStats/AsymptoticCalculator.h"
#include "RooStats/HybridCalculator.h"
#include "RooStats/FrequentistCalculator.h"
#include "RooStats/ToyMCSampler.h"
#include "RooStats/HypoTestPlot.h"

#include "RooStats/NumEventsTestStat.h"
#include "RooStats/ProfileLikelihoodTestStat.h"
#include "RooStats/SimpleLikelihoodRatioTestStat.h"
#include "RooStats/RatioOfProfiledLikelihoodsTestStat.h"
#include "RooStats/MaxLikelihoodEstimateTestStat.h"

#include "RooStats/HypoTestInverter.h"
#include "RooStats/HypoTestInverterResult.h"
#include "RooStats/HypoTestInverterPlot.h"

using namespace RooFit;
using namespace RooStats;


bool noSystematics = false;              // force all systematics to be off (i.e. set all nuisance parameters as constat
double nToysRatio = 4;                   // ratio Ntoys Null/ntoys ALT
double poiValue = -1;                    // change poi snapshot value for S+B model (needed for expected p0 values)
int  printLevel=0;
bool generateBinned = false;             // for binned generation

void StandardHypoTestDemo(const char* infile = "",
                          const char* workspaceName = "combined",
                          const char* modelSBName = "ModelConfig",
                          const char* modelBName = "",
                          const char* dataName = "obsData",
                          int calcType = 0, // 0 freq 1 hybrid, 2 asymptotic
                          int testStatType = 3,   // 0 LEP, 1 TeV, 2 LHC, 3 LHC - one sided
                          int ntoys = 5000,
                          bool useNC = false,
                          const char * nuisPriorName = 0)
{

/*

  Other Parameter to pass in tutorial
  apart from standard for filename, ws, modelconfig and data

  type = 0 Freq calculator
  type = 1 Hybrid calculator
  type = 2 Asymptotic calculator

  testStatType = 0 LEP
  = 1 Tevatron
  = 2 Profile Likelihood
  = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)

  ntoys:         number of toys to use

  useNumberCounting:  set to true when using number counting events

  nuisPriorName:   name of prior for the nnuisance. This is often expressed as constraint term in the global model
  It is needed only when using the HybridCalculator (type=1)
  If not given by default the prior pdf from ModelConfig is used.

  extra options are available as global paramwters of the macro. They major ones are:

  generateBinned       generate binned data sets for toys (default is false) - be careful not to activate with
  a too large (>=3) number of observables
  nToyRatio            ratio of S+B/B toys (default is 2)
  printLevel

*/

   // disable - can cause some problems
   //ToyMCSampler::SetAlwaysUseMultiGen(true);

   SimpleLikelihoodRatioTestStat::SetAlwaysReuseNLL(true);
   ProfileLikelihoodTestStat::SetAlwaysReuseNLL(true);
   RatioOfProfiledLikelihoodsTestStat::SetAlwaysReuseNLL(true);

   //RooRandom::randomGenerator()->SetSeed(0);

   // to change minimizers
   // ROOT::Math::MinimizerOptions::SetDefaultStrategy(0);
   // ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2");
   // ROOT::Math::MinimizerOptions::SetDefaultTolerance(1);

  /////////////////////////////////////////////////////////////
  // First part is just to access a user-defined file
  // or create the standard example file if it doesn't exist
  ////////////////////////////////////////////////////////////
   const char* filename = "";
   if (!strcmp(infile,"")) {
      filename = "results/example_combined_GaussExample_model.root";
      bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
      // if file does not exists generate with histfactory
      if (!fileExist) {
#ifdef _WIN32
         cout << "HistFactory file cannot be generated on Windows - exit" << endl;
         return;
#endif
         // Normally this would be run on the command line
         cout <<"will run standard hist2workspace example"<<endl;
         gROOT->ProcessLine(".! prepareHistFactory .");
         gROOT->ProcessLine(".! hist2workspace config/example.xml");
         cout <<"\n\n---------------------"<<endl;
         cout <<"Done creating example input"<<endl;
         cout <<"---------------------\n\n"<<endl;
      }

   }
   else
      filename = infile;

   // Try to open the file
   TFile *file = TFile::Open(filename);

   // if input file was specified byt not found, quit
   if(!file ){
      cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
      return;
   }


  /////////////////////////////////////////////////////////////
  // Tutorial starts here
  ////////////////////////////////////////////////////////////

  // get the workspace out of the file
  RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
  if(!w){
    cout <<"workspace not found" << endl;
    return;
  }
  w->Print();

  // get the modelConfig out of the file
  ModelConfig* sbModel = (ModelConfig*) w->obj(modelSBName);


  // get the modelConfig out of the file
  RooAbsData* data = w->data(dataName);

  // make sure ingredients are found
  if(!data || !sbModel){
    w->Print();
    cout << "data or ModelConfig was not found" <<endl;
    return;
  }
  // make b model
  ModelConfig* bModel = (ModelConfig*) w->obj(modelBName);


   // case of no systematics
   // remove nuisance parameters from model
   if (noSystematics) {
      const RooArgSet * nuisPar = sbModel->GetNuisanceParameters();
      if (nuisPar && nuisPar->getSize() > 0) {
         std::cout << "StandardHypoTestInvDemo" << "  -  Switch off all systematics by setting them constant to their initial values" << std::endl;
         RooStats::SetAllConstant(*nuisPar);
      }
      if (bModel) {
         const RooArgSet * bnuisPar = bModel->GetNuisanceParameters();
         if (bnuisPar)
            RooStats::SetAllConstant(*bnuisPar);
      }
   }


  if (!bModel ) {
      Info("StandardHypoTestInvDemo","The background model %s does not exist",modelBName);
      Info("StandardHypoTestInvDemo","Copy it from ModelConfig %s and set POI to zero",modelSBName);
      bModel = (ModelConfig*) sbModel->Clone();
      bModel->SetName(TString(modelSBName)+TString("B_only"));
      RooRealVar * var = dynamic_cast<RooRealVar*>(bModel->GetParametersOfInterest()->first());
      if (!var) return;
      double oldval = var->getVal();
      var->setVal(0);
      //bModel->SetSnapshot( RooArgSet(*var, *w->var("lumi"))  );
      bModel->SetSnapshot( RooArgSet(*var)  );
      var->setVal(oldval);
  }

   if (!sbModel->GetSnapshot() || poiValue > 0) {
      Info("StandardHypoTestDemo","Model %s has no snapshot  - make one using model poi",modelSBName);
      RooRealVar * var = dynamic_cast<RooRealVar*>(sbModel->GetParametersOfInterest()->first());
      if (!var) return;
      double oldval = var->getVal();
      if (poiValue > 0)  var->setVal(poiValue);
      //sbModel->SetSnapshot( RooArgSet(*var, *w->var("lumi") ) );
      sbModel->SetSnapshot( RooArgSet(*var) );
      if (poiValue > 0) var->setVal(oldval);
      //sbModel->SetSnapshot( *sbModel->GetParametersOfInterest() );
   }





   // part 1, hypothesis testing
   SimpleLikelihoodRatioTestStat * slrts = new SimpleLikelihoodRatioTestStat(*bModel->GetPdf(), *sbModel->GetPdf());
   // null parameters must includes snapshot of poi plus the nuisance values
   RooArgSet nullParams(*bModel->GetSnapshot());
   if (bModel->GetNuisanceParameters()) nullParams.add(*bModel->GetNuisanceParameters());

   slrts->SetNullParameters(nullParams);
   RooArgSet altParams(*sbModel->GetSnapshot());
   if (sbModel->GetNuisanceParameters()) altParams.add(*sbModel->GetNuisanceParameters());
   slrts->SetAltParameters(altParams);


   ProfileLikelihoodTestStat * profll = new ProfileLikelihoodTestStat(*bModel->GetPdf());


   RatioOfProfiledLikelihoodsTestStat *
      ropl = new RatioOfProfiledLikelihoodsTestStat(*bModel->GetPdf(), *sbModel->GetPdf(), sbModel->GetSnapshot());
   ropl->SetSubtractMLE(false);

   if (testStatType == 3) profll->SetOneSidedDiscovery(1);
   profll->SetPrintLevel(printLevel);

   // profll.SetReuseNLL(mOptimize);
   // slrts.SetReuseNLL(mOptimize);
   // ropl.SetReuseNLL(mOptimize);

   AsymptoticCalculator::SetPrintLevel(printLevel);

   HypoTestCalculatorGeneric *  hypoCalc = 0;
   // note here Null is B and Alt is S+B
   if (calcType == 0) hypoCalc = new  FrequentistCalculator(*data, *sbModel, *bModel);
   else if (calcType == 1) hypoCalc= new  HybridCalculator(*data, *sbModel, *bModel);
   else if (calcType == 2) hypoCalc= new  AsymptoticCalculator(*data, *sbModel, *bModel);

   if (calcType == 0)
       ((FrequentistCalculator*)hypoCalc)->SetToys(ntoys, ntoys/nToysRatio);
   if (calcType == 1)
       ((HybridCalculator*)hypoCalc)->SetToys(ntoys, ntoys/nToysRatio);
   if (calcType == 2 ) {
      if (testStatType == 3) ((AsymptoticCalculator*) hypoCalc)->SetOneSidedDiscovery(true);
      if (testStatType != 2 && testStatType != 3)
         Warning("StandardHypoTestDemo","Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");


   }


   // check for nuisance prior pdf in case of nuisance parameters
   if (calcType == 1 && (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters() )) {
         RooAbsPdf * nuisPdf = 0;
         if (nuisPriorName) nuisPdf = w->pdf(nuisPriorName);
         // use prior defined first in bModel (then in SbModel)
         if (!nuisPdf)  {
            Info("StandardHypoTestDemo","No nuisance pdf given for the HybridCalculator - try to deduce  pdf from the   model");
            if (bModel->GetPdf() && bModel->GetObservables() )
               nuisPdf = RooStats::MakeNuisancePdf(*bModel,"nuisancePdf_bmodel");
            else
               nuisPdf = RooStats::MakeNuisancePdf(*sbModel,"nuisancePdf_sbmodel");
         }
         if (!nuisPdf ) {
            if (bModel->GetPriorPdf())  {
               nuisPdf = bModel->GetPriorPdf();
               Info("StandardHypoTestDemo","No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",nuisPdf->GetName());
            }
            else {
               Error("StandardHypoTestDemo","Cannnot run Hybrid calculator because no prior on the nuisance parameter is specified or can be derived");
               return;
            }
         }
         assert(nuisPdf);
         Info("StandardHypoTestDemo","Using as nuisance Pdf ... " );
         nuisPdf->Print();

         const RooArgSet * nuisParams = (bModel->GetNuisanceParameters() ) ? bModel->GetNuisanceParameters() : sbModel->GetNuisanceParameters();
         RooArgSet * np = nuisPdf->getObservables(*nuisParams);
         if (np->getSize() == 0) {
            Warning("StandardHypoTestDemo","Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
         }
         delete np;

         ((HybridCalculator*)hypoCalc)->ForcePriorNuisanceAlt(*nuisPdf);
         ((HybridCalculator*)hypoCalc)->ForcePriorNuisanceNull(*nuisPdf);
   }

   // hypoCalc->ForcePriorNuisanceAlt(*sbModel->GetPriorPdf());
   // hypoCalc->ForcePriorNuisanceNull(*bModel->GetPriorPdf());

   ToyMCSampler * sampler = (ToyMCSampler *)hypoCalc->GetTestStatSampler();

   if (sampler && (calcType == 0 || calcType == 1) ) {

      // look if pdf is number counting or extended
      if (sbModel->GetPdf()->canBeExtended() ) {
         if (useNC)   Warning("StandardHypoTestDemo","Pdf is extended: but number counting flag is set: ignore it ");
      }
      else {
         // for not extended pdf
         if (!useNC)  {
            int nEvents = data->numEntries();
            Info("StandardHypoTestDemo","Pdf is not extended: number of events to generate taken  from observed data set is %d",nEvents);
            sampler->SetNEventsPerToy(nEvents);
         }
         else {
            Info("StandardHypoTestDemo","using a number counting pdf");
            sampler->SetNEventsPerToy(1);
         }
      }

      if (data->isWeighted() && !generateBinned) {
         Info("StandardHypoTestDemo","Data set is weighted, nentries = %d and sum of weights = %8.1f but toy generation is unbinned - it would be faster to set generateBinned to true\n",data->numEntries(), data->sumEntries());
      }
      if (generateBinned)  sampler->SetGenerateBinned(generateBinned);


      // set the test statistic
      if (testStatType == 0) sampler->SetTestStatistic(slrts);
      if (testStatType == 1) sampler->SetTestStatistic(ropl);
      if (testStatType == 2 || testStatType == 3) sampler->SetTestStatistic(profll);

   }

   HypoTestResult *  htr = hypoCalc->GetHypoTest();
   htr->SetPValueIsRightTail(true);
   htr->SetBackgroundAsAlt(false);
   htr->Print(); // how to get meaningfull CLs at this point?

   delete sampler;
   delete slrts;
   delete ropl;
   delete profll;

   if (calcType != 2) {
      HypoTestPlot * plot = new HypoTestPlot(*htr,100);
      plot->SetLogYaxis(true);
      plot->Draw();
   }
   else {
      std::cout << "Asymptotic results " << std::endl;

   }

   // look at expected significances
   // found median of S+B distribution
   if (calcType != 2) {

      SamplingDistribution * altDist = htr->GetAltDistribution();
      HypoTestResult htExp("Expected Result");
      htExp.Append(htr);
      // find quantiles in alt (S+B) distribution
      double p[5];
      double q[5];
      for (int i = 0; i < 5; ++i) {
         double sig = -2  + i;
         p[i] = ROOT::Math::normal_cdf(sig,1);
      }
      std::vector<double> values = altDist->GetSamplingDistribution();
      TMath::Quantiles( values.size(), 5, &values[0], q, p, false);

      for (int i = 0; i < 5; ++i) {
         htExp.SetTestStatisticData( q[i] );
         double sig = -2  + i;
         std::cout << " Expected p -value and significance at " << sig << " sigma = "
                   << htExp.NullPValue() << " significance " << htExp.Significance() << " sigma " << std::endl;

      }
   }
   else {
      // case of asymptotic calculator
      for (int i = 0; i < 5; ++i) {
         double sig = -2  + i;
         // sigma is inverted here
         double pval = AsymptoticCalculator::GetExpectedPValues( htr->NullPValue(), htr->AlternatePValue(), -sig, false);
         std::cout << " Expected p -value and significance at " << sig << " sigma = "
                   << pval << " significance " << ROOT::Math::normal_quantile_c(pval,1) << " sigma " << std::endl;

      }
   }

}

 StandardHypoTestDemo.C:1
 StandardHypoTestDemo.C:2
 StandardHypoTestDemo.C:3
 StandardHypoTestDemo.C:4
 StandardHypoTestDemo.C:5
 StandardHypoTestDemo.C:6
 StandardHypoTestDemo.C:7
 StandardHypoTestDemo.C:8
 StandardHypoTestDemo.C:9
 StandardHypoTestDemo.C:10
 StandardHypoTestDemo.C:11
 StandardHypoTestDemo.C:12
 StandardHypoTestDemo.C:13
 StandardHypoTestDemo.C:14
 StandardHypoTestDemo.C:15
 StandardHypoTestDemo.C:16
 StandardHypoTestDemo.C:17
 StandardHypoTestDemo.C:18
 StandardHypoTestDemo.C:19
 StandardHypoTestDemo.C:20
 StandardHypoTestDemo.C:21
 StandardHypoTestDemo.C:22
 StandardHypoTestDemo.C:23
 StandardHypoTestDemo.C:24
 StandardHypoTestDemo.C:25
 StandardHypoTestDemo.C:26
 StandardHypoTestDemo.C:27
 StandardHypoTestDemo.C:28
 StandardHypoTestDemo.C:29
 StandardHypoTestDemo.C:30
 StandardHypoTestDemo.C:31
 StandardHypoTestDemo.C:32
 StandardHypoTestDemo.C:33
 StandardHypoTestDemo.C:34
 StandardHypoTestDemo.C:35
 StandardHypoTestDemo.C:36
 StandardHypoTestDemo.C:37
 StandardHypoTestDemo.C:38
 StandardHypoTestDemo.C:39
 StandardHypoTestDemo.C:40
 StandardHypoTestDemo.C:41
 StandardHypoTestDemo.C:42
 StandardHypoTestDemo.C:43
 StandardHypoTestDemo.C:44
 StandardHypoTestDemo.C:45
 StandardHypoTestDemo.C:46
 StandardHypoTestDemo.C:47
 StandardHypoTestDemo.C:48
 StandardHypoTestDemo.C:49
 StandardHypoTestDemo.C:50
 StandardHypoTestDemo.C:51
 StandardHypoTestDemo.C:52
 StandardHypoTestDemo.C:53
 StandardHypoTestDemo.C:54
 StandardHypoTestDemo.C:55
 StandardHypoTestDemo.C:56
 StandardHypoTestDemo.C:57
 StandardHypoTestDemo.C:58
 StandardHypoTestDemo.C:59
 StandardHypoTestDemo.C:60
 StandardHypoTestDemo.C:61
 StandardHypoTestDemo.C:62
 StandardHypoTestDemo.C:63
 StandardHypoTestDemo.C:64
 StandardHypoTestDemo.C:65
 StandardHypoTestDemo.C:66
 StandardHypoTestDemo.C:67
 StandardHypoTestDemo.C:68
 StandardHypoTestDemo.C:69
 StandardHypoTestDemo.C:70
 StandardHypoTestDemo.C:71
 StandardHypoTestDemo.C:72
 StandardHypoTestDemo.C:73
 StandardHypoTestDemo.C:74
 StandardHypoTestDemo.C:75
 StandardHypoTestDemo.C:76
 StandardHypoTestDemo.C:77
 StandardHypoTestDemo.C:78
 StandardHypoTestDemo.C:79
 StandardHypoTestDemo.C:80
 StandardHypoTestDemo.C:81
 StandardHypoTestDemo.C:82
 StandardHypoTestDemo.C:83
 StandardHypoTestDemo.C:84
 StandardHypoTestDemo.C:85
 StandardHypoTestDemo.C:86
 StandardHypoTestDemo.C:87
 StandardHypoTestDemo.C:88
 StandardHypoTestDemo.C:89
 StandardHypoTestDemo.C:90
 StandardHypoTestDemo.C:91
 StandardHypoTestDemo.C:92
 StandardHypoTestDemo.C:93
 StandardHypoTestDemo.C:94
 StandardHypoTestDemo.C:95
 StandardHypoTestDemo.C:96
 StandardHypoTestDemo.C:97
 StandardHypoTestDemo.C:98
 StandardHypoTestDemo.C:99
 StandardHypoTestDemo.C:100
 StandardHypoTestDemo.C:101
 StandardHypoTestDemo.C:102
 StandardHypoTestDemo.C:103
 StandardHypoTestDemo.C:104
 StandardHypoTestDemo.C:105
 StandardHypoTestDemo.C:106
 StandardHypoTestDemo.C:107
 StandardHypoTestDemo.C:108
 StandardHypoTestDemo.C:109
 StandardHypoTestDemo.C:110
 StandardHypoTestDemo.C:111
 StandardHypoTestDemo.C:112
 StandardHypoTestDemo.C:113
 StandardHypoTestDemo.C:114
 StandardHypoTestDemo.C:115
 StandardHypoTestDemo.C:116
 StandardHypoTestDemo.C:117
 StandardHypoTestDemo.C:118
 StandardHypoTestDemo.C:119
 StandardHypoTestDemo.C:120
 StandardHypoTestDemo.C:121
 StandardHypoTestDemo.C:122
 StandardHypoTestDemo.C:123
 StandardHypoTestDemo.C:124
 StandardHypoTestDemo.C:125
 StandardHypoTestDemo.C:126
 StandardHypoTestDemo.C:127
 StandardHypoTestDemo.C:128
 StandardHypoTestDemo.C:129
 StandardHypoTestDemo.C:130
 StandardHypoTestDemo.C:131
 StandardHypoTestDemo.C:132
 StandardHypoTestDemo.C:133
 StandardHypoTestDemo.C:134
 StandardHypoTestDemo.C:135
 StandardHypoTestDemo.C:136
 StandardHypoTestDemo.C:137
 StandardHypoTestDemo.C:138
 StandardHypoTestDemo.C:139
 StandardHypoTestDemo.C:140
 StandardHypoTestDemo.C:141
 StandardHypoTestDemo.C:142
 StandardHypoTestDemo.C:143
 StandardHypoTestDemo.C:144
 StandardHypoTestDemo.C:145
 StandardHypoTestDemo.C:146
 StandardHypoTestDemo.C:147
 StandardHypoTestDemo.C:148
 StandardHypoTestDemo.C:149
 StandardHypoTestDemo.C:150
 StandardHypoTestDemo.C:151
 StandardHypoTestDemo.C:152
 StandardHypoTestDemo.C:153
 StandardHypoTestDemo.C:154
 StandardHypoTestDemo.C:155
 StandardHypoTestDemo.C:156
 StandardHypoTestDemo.C:157
 StandardHypoTestDemo.C:158
 StandardHypoTestDemo.C:159
 StandardHypoTestDemo.C:160
 StandardHypoTestDemo.C:161
 StandardHypoTestDemo.C:162
 StandardHypoTestDemo.C:163
 StandardHypoTestDemo.C:164
 StandardHypoTestDemo.C:165
 StandardHypoTestDemo.C:166
 StandardHypoTestDemo.C:167
 StandardHypoTestDemo.C:168
 StandardHypoTestDemo.C:169
 StandardHypoTestDemo.C:170
 StandardHypoTestDemo.C:171
 StandardHypoTestDemo.C:172
 StandardHypoTestDemo.C:173
 StandardHypoTestDemo.C:174
 StandardHypoTestDemo.C:175
 StandardHypoTestDemo.C:176
 StandardHypoTestDemo.C:177
 StandardHypoTestDemo.C:178
 StandardHypoTestDemo.C:179
 StandardHypoTestDemo.C:180
 StandardHypoTestDemo.C:181
 StandardHypoTestDemo.C:182
 StandardHypoTestDemo.C:183
 StandardHypoTestDemo.C:184
 StandardHypoTestDemo.C:185
 StandardHypoTestDemo.C:186
 StandardHypoTestDemo.C:187
 StandardHypoTestDemo.C:188
 StandardHypoTestDemo.C:189
 StandardHypoTestDemo.C:190
 StandardHypoTestDemo.C:191
 StandardHypoTestDemo.C:192
 StandardHypoTestDemo.C:193
 StandardHypoTestDemo.C:194
 StandardHypoTestDemo.C:195
 StandardHypoTestDemo.C:196
 StandardHypoTestDemo.C:197
 StandardHypoTestDemo.C:198
 StandardHypoTestDemo.C:199
 StandardHypoTestDemo.C:200
 StandardHypoTestDemo.C:201
 StandardHypoTestDemo.C:202
 StandardHypoTestDemo.C:203
 StandardHypoTestDemo.C:204
 StandardHypoTestDemo.C:205
 StandardHypoTestDemo.C:206
 StandardHypoTestDemo.C:207
 StandardHypoTestDemo.C:208
 StandardHypoTestDemo.C:209
 StandardHypoTestDemo.C:210
 StandardHypoTestDemo.C:211
 StandardHypoTestDemo.C:212
 StandardHypoTestDemo.C:213
 StandardHypoTestDemo.C:214
 StandardHypoTestDemo.C:215
 StandardHypoTestDemo.C:216
 StandardHypoTestDemo.C:217
 StandardHypoTestDemo.C:218
 StandardHypoTestDemo.C:219
 StandardHypoTestDemo.C:220
 StandardHypoTestDemo.C:221
 StandardHypoTestDemo.C:222
 StandardHypoTestDemo.C:223
 StandardHypoTestDemo.C:224
 StandardHypoTestDemo.C:225
 StandardHypoTestDemo.C:226
 StandardHypoTestDemo.C:227
 StandardHypoTestDemo.C:228
 StandardHypoTestDemo.C:229
 StandardHypoTestDemo.C:230
 StandardHypoTestDemo.C:231
 StandardHypoTestDemo.C:232
 StandardHypoTestDemo.C:233
 StandardHypoTestDemo.C:234
 StandardHypoTestDemo.C:235
 StandardHypoTestDemo.C:236
 StandardHypoTestDemo.C:237
 StandardHypoTestDemo.C:238
 StandardHypoTestDemo.C:239
 StandardHypoTestDemo.C:240
 StandardHypoTestDemo.C:241
 StandardHypoTestDemo.C:242
 StandardHypoTestDemo.C:243
 StandardHypoTestDemo.C:244
 StandardHypoTestDemo.C:245
 StandardHypoTestDemo.C:246
 StandardHypoTestDemo.C:247
 StandardHypoTestDemo.C:248
 StandardHypoTestDemo.C:249
 StandardHypoTestDemo.C:250
 StandardHypoTestDemo.C:251
 StandardHypoTestDemo.C:252
 StandardHypoTestDemo.C:253
 StandardHypoTestDemo.C:254
 StandardHypoTestDemo.C:255
 StandardHypoTestDemo.C:256
 StandardHypoTestDemo.C:257
 StandardHypoTestDemo.C:258
 StandardHypoTestDemo.C:259
 StandardHypoTestDemo.C:260
 StandardHypoTestDemo.C:261
 StandardHypoTestDemo.C:262
 StandardHypoTestDemo.C:263
 StandardHypoTestDemo.C:264
 StandardHypoTestDemo.C:265
 StandardHypoTestDemo.C:266
 StandardHypoTestDemo.C:267
 StandardHypoTestDemo.C:268
 StandardHypoTestDemo.C:269
 StandardHypoTestDemo.C:270
 StandardHypoTestDemo.C:271
 StandardHypoTestDemo.C:272
 StandardHypoTestDemo.C:273
 StandardHypoTestDemo.C:274
 StandardHypoTestDemo.C:275
 StandardHypoTestDemo.C:276
 StandardHypoTestDemo.C:277
 StandardHypoTestDemo.C:278
 StandardHypoTestDemo.C:279
 StandardHypoTestDemo.C:280
 StandardHypoTestDemo.C:281
 StandardHypoTestDemo.C:282
 StandardHypoTestDemo.C:283
 StandardHypoTestDemo.C:284
 StandardHypoTestDemo.C:285
 StandardHypoTestDemo.C:286
 StandardHypoTestDemo.C:287
 StandardHypoTestDemo.C:288
 StandardHypoTestDemo.C:289
 StandardHypoTestDemo.C:290
 StandardHypoTestDemo.C:291
 StandardHypoTestDemo.C:292
 StandardHypoTestDemo.C:293
 StandardHypoTestDemo.C:294
 StandardHypoTestDemo.C:295
 StandardHypoTestDemo.C:296
 StandardHypoTestDemo.C:297
 StandardHypoTestDemo.C:298
 StandardHypoTestDemo.C:299
 StandardHypoTestDemo.C:300
 StandardHypoTestDemo.C:301
 StandardHypoTestDemo.C:302
 StandardHypoTestDemo.C:303
 StandardHypoTestDemo.C:304
 StandardHypoTestDemo.C:305
 StandardHypoTestDemo.C:306
 StandardHypoTestDemo.C:307
 StandardHypoTestDemo.C:308
 StandardHypoTestDemo.C:309
 StandardHypoTestDemo.C:310
 StandardHypoTestDemo.C:311
 StandardHypoTestDemo.C:312
 StandardHypoTestDemo.C:313
 StandardHypoTestDemo.C:314
 StandardHypoTestDemo.C:315
 StandardHypoTestDemo.C:316
 StandardHypoTestDemo.C:317
 StandardHypoTestDemo.C:318
 StandardHypoTestDemo.C:319
 StandardHypoTestDemo.C:320
 StandardHypoTestDemo.C:321
 StandardHypoTestDemo.C:322
 StandardHypoTestDemo.C:323
 StandardHypoTestDemo.C:324
 StandardHypoTestDemo.C:325
 StandardHypoTestDemo.C:326
 StandardHypoTestDemo.C:327
 StandardHypoTestDemo.C:328
 StandardHypoTestDemo.C:329
 StandardHypoTestDemo.C:330
 StandardHypoTestDemo.C:331
 StandardHypoTestDemo.C:332
 StandardHypoTestDemo.C:333
 StandardHypoTestDemo.C:334
 StandardHypoTestDemo.C:335
 StandardHypoTestDemo.C:336
 StandardHypoTestDemo.C:337
 StandardHypoTestDemo.C:338
 StandardHypoTestDemo.C:339
 StandardHypoTestDemo.C:340
 StandardHypoTestDemo.C:341
 StandardHypoTestDemo.C:342
 StandardHypoTestDemo.C:343
 StandardHypoTestDemo.C:344
 StandardHypoTestDemo.C:345
 StandardHypoTestDemo.C:346
 StandardHypoTestDemo.C:347
 StandardHypoTestDemo.C:348
 StandardHypoTestDemo.C:349
 StandardHypoTestDemo.C:350
 StandardHypoTestDemo.C:351
 StandardHypoTestDemo.C:352
 StandardHypoTestDemo.C:353
 StandardHypoTestDemo.C:354
 StandardHypoTestDemo.C:355
 StandardHypoTestDemo.C:356
 StandardHypoTestDemo.C:357
 StandardHypoTestDemo.C:358
 StandardHypoTestDemo.C:359
 StandardHypoTestDemo.C:360
 StandardHypoTestDemo.C:361
 StandardHypoTestDemo.C:362
 StandardHypoTestDemo.C:363
 StandardHypoTestDemo.C:364
 StandardHypoTestDemo.C:365
 StandardHypoTestDemo.C:366
 StandardHypoTestDemo.C:367
 StandardHypoTestDemo.C:368
 StandardHypoTestDemo.C:369
 StandardHypoTestDemo.C:370
 StandardHypoTestDemo.C:371
 StandardHypoTestDemo.C:372
 StandardHypoTestDemo.C:373
 StandardHypoTestDemo.C:374
 StandardHypoTestDemo.C:375
 StandardHypoTestDemo.C:376
 StandardHypoTestDemo.C:377
 StandardHypoTestDemo.C:378
 StandardHypoTestDemo.C:379
 StandardHypoTestDemo.C:380
 StandardHypoTestDemo.C:381
 StandardHypoTestDemo.C:382
 StandardHypoTestDemo.C:383
 StandardHypoTestDemo.C:384
 StandardHypoTestDemo.C:385
 StandardHypoTestDemo.C:386
 StandardHypoTestDemo.C:387
 StandardHypoTestDemo.C:388
 StandardHypoTestDemo.C:389
 StandardHypoTestDemo.C:390
 StandardHypoTestDemo.C:391
 StandardHypoTestDemo.C:392
 StandardHypoTestDemo.C:393
 StandardHypoTestDemo.C:394
 StandardHypoTestDemo.C:395
 StandardHypoTestDemo.C:396
 StandardHypoTestDemo.C:397
 StandardHypoTestDemo.C:398
 StandardHypoTestDemo.C:399
 StandardHypoTestDemo.C:400
 StandardHypoTestDemo.C:401
 StandardHypoTestDemo.C:402
 StandardHypoTestDemo.C:403
 StandardHypoTestDemo.C:404
 StandardHypoTestDemo.C:405
 StandardHypoTestDemo.C:406
 StandardHypoTestDemo.C:407
 StandardHypoTestDemo.C:408
 StandardHypoTestDemo.C:409
 StandardHypoTestDemo.C:410
 StandardHypoTestDemo.C:411
 StandardHypoTestDemo.C:412
 StandardHypoTestDemo.C:413
 StandardHypoTestDemo.C:414