// @(#)root/tmva $Id$
// Author: Andreas Hoecker, Peter Speckmayer, Matt Jachowski, Jan Therhaag, Jiahang Zhong

/**********************************************************************************
 * Project: TMVA - a Root-integrated toolkit for multivariate data analysis       *
 * Package: TMVA                                                                  *
 * Class  : MethodANNBase                                                         *
 * Web    : http://tmva.sourceforge.net                                           *
 *                                                                                *
 * Description:                                                                   *
 *      Artificial neural network base class for the discrimination of signal     *
 *      from background.                                                          *
 *                                                                                *
 * Authors (alphabetical):                                                        *
 *      Krzysztof Danielowski <danielow@cern.ch>       - IFJ & AGH, Poland        *
 *      Andreas Hoecker       <Andreas.Hocker@cern.ch> - CERN, Switzerland        *
 *      Matt Jachowski        <jachowski@stanford.edu> - Stanford University, USA *
 *      Kamil Kraszewski      <kalq@cern.ch>           - IFJ & UJ, Poland         *
 *      Maciej Kruk           <mkruk@cern.ch>          - IFJ & AGH, Poland        *
 *      Peter Speckmayer      <peter.speckmayer@cern.ch> - CERN, Switzerland      *
 *      Joerg Stelzer         <stelzer@cern.ch>        - DESY, Germany            *
 *      Jan Therhaag       <Jan.Therhaag@cern.ch>     - U of Bonn, Germany        *
 *      Jiahang Zhong         <Jiahang.Zhong@cern.ch>  - Academia Sinica, Taipei  *
 *                                                                                *
 * Copyright (c) 2005-2011:                                                       *
 *      CERN, Switzerland                                                         *
 *      U. of Bonn, Germany                                                       *
 *                                                                                *
 * Redistribution and use in source and binary forms, with or without             *
 * modification, are permitted according to the terms listed in LICENSE           *
 * (http://tmva.sourceforge.net/LICENSE)                                          *
 **********************************************************************************/

//_______________________________________________________________________
//
// Base class for all TMVA methods using artificial neural networks
//
//_______________________________________________________________________

#include <vector>
#include <cstdlib>
#include <stdexcept>
#if __cplusplus > 199711L
#include <atomic>
#endif

#include "TString.h"
#include "TTree.h"
#include "TDirectory.h"
#include "Riostream.h"
#include "TRandom3.h"
#include "TH2F.h"
#include "TH1.h"
#include "TMath.h"

#include "TMVA/MethodBase.h"
#include "TMVA/MethodANNBase.h"
#include "TMVA/TNeuron.h"
#include "TMVA/TSynapse.h"
#include "TMVA/TActivationChooser.h"
#include "TMVA/Types.h"
#include "TMVA/Tools.h"
#include "TMVA/TNeuronInputChooser.h"
#include "TMVA/Ranking.h"

using std::vector;

ClassImp(TMVA::MethodANNBase)

//______________________________________________________________________________
TMVA::MethodANNBase::MethodANNBase( const TString& jobName,
                                    Types::EMVA methodType,
                                    const TString& methodTitle,
                                    DataSetInfo& theData,
                                    const TString& theOption,
                                    TDirectory* theTargetDir )
   : TMVA::MethodBase( jobName, methodType, methodTitle, theData, theOption, theTargetDir )
   , fEstimator(kMSE)
   , fUseRegulator(kFALSE)
   , fRandomSeed(0)
{
   // standard constructor
   // Note: Right now it is an option to choose the neuron input function,
   // but only the input function "sum" leads to weight convergence --
   // otherwise the weights go to nan and lead to an ABORT.
   InitANNBase();

   DeclareOptions();
}

//______________________________________________________________________________
TMVA::MethodANNBase::MethodANNBase( Types::EMVA methodType,
                                    DataSetInfo& theData,
                                    const TString& theWeightFile,
                                    TDirectory* theTargetDir )
   : TMVA::MethodBase( methodType, theData, theWeightFile, theTargetDir )
   , fEstimator(kMSE)
   , fUseRegulator(kFALSE)
   , fRandomSeed(0)
{
   // construct the Method from the weight file
   InitANNBase();

   DeclareOptions();
}

//______________________________________________________________________________
void TMVA::MethodANNBase::DeclareOptions()
{
   // define the options (their key words) that can be set in the option string
   // here the options valid for ALL MVA methods are declared.
   // know options: NCycles=xx              :the number of training cycles
   //               Normalize=kTRUE,kFALSe  :if normalised in put variables should be used
   //               HiddenLayser="N-1,N-2"  :the specification of the hidden layers
   //               NeuronType=sigmoid,tanh,radial,linar  : the type of activation function
   //                                                       used at the neuronn
   //

   DeclareOptionRef( fNcycles    = 500,       "NCycles",         "Number of training cycles" );
   DeclareOptionRef( fLayerSpec  = "N,N-1",   "HiddenLayers",    "Specification of hidden layer architecture" );
   DeclareOptionRef( fNeuronType = "sigmoid", "NeuronType",      "Neuron activation function type" );
   DeclareOptionRef( fRandomSeed = 1, "RandomSeed", "Random seed for initial synapse weights (0 means unique seed for each run; default value '1')");

   DeclareOptionRef(fEstimatorS="MSE", "EstimatorType",
                    "MSE (Mean Square Estimator) for Gaussian Likelihood or CE(Cross-Entropy) for Bernoulli Likelihood" ); //zjh
   AddPreDefVal(TString("MSE"));  //zjh
   AddPreDefVal(TString("CE"));   //zjh


   TActivationChooser aChooser;
   std::vector<TString>* names = aChooser.GetAllActivationNames();
   Int_t nTypes = names->size();
   for (Int_t i = 0; i < nTypes; i++)
      AddPreDefVal(names->at(i));
   delete names;

   DeclareOptionRef(fNeuronInputType="sum", "NeuronInputType","Neuron input function type");
   TNeuronInputChooser iChooser;
   names = iChooser.GetAllNeuronInputNames();
   nTypes = names->size();
   for (Int_t i = 0; i < nTypes; i++) AddPreDefVal(names->at(i));
   delete names;
}


//______________________________________________________________________________
void TMVA::MethodANNBase::ProcessOptions()
{
   // do nothing specific at this moment
   if      ( DoRegression() || DoMulticlass())  fEstimatorS = "MSE";    //zjh
   if      (fEstimatorS == "MSE" )  fEstimator = kMSE;   
   else if (fEstimatorS == "CE")    fEstimator = kCE;      //zjh
   std::vector<Int_t>* layout = ParseLayoutString(fLayerSpec);
   BuildNetwork(layout);
   delete layout;
}

//______________________________________________________________________________
std::vector<Int_t>* TMVA::MethodANNBase::ParseLayoutString(TString layerSpec)
{
   // parse layout specification string and return a vector, each entry
   // containing the number of neurons to go in each successive layer
   std::vector<Int_t>* layout = new std::vector<Int_t>();
   layout->push_back((Int_t)GetNvar());
   while(layerSpec.Length()>0) {
      TString sToAdd="";
      if (layerSpec.First(',')<0) {
         sToAdd = layerSpec;
         layerSpec = "";
      }
      else {
         sToAdd = layerSpec(0,layerSpec.First(','));
         layerSpec = layerSpec(layerSpec.First(',')+1,layerSpec.Length());
      }
      int nNodes = 0;
      if (sToAdd.BeginsWith("n") || sToAdd.BeginsWith("N")) { sToAdd.Remove(0,1); nNodes = GetNvar(); }
      nNodes += atoi(sToAdd);
      layout->push_back(nNodes);
   }
   if( DoRegression() )
      layout->push_back( DataInfo().GetNTargets() );  // one output node for each target
   else if( DoMulticlass() )
      layout->push_back( DataInfo().GetNClasses() );  // one output node for each class
   else
      layout->push_back(1);  // one output node (for signal/background classification)

   int n = 0;
   for( std::vector<Int_t>::iterator it = layout->begin(); it != layout->end(); it++ ){
      n++;
   }

   return layout;
}

//______________________________________________________________________________
void TMVA::MethodANNBase::InitANNBase()
{
   // initialize ANNBase object
   fNetwork         = NULL;
   frgen            = NULL;
   fActivation      = NULL;
   fOutput          = NULL; //zjh
   fIdentity        = NULL;
   fInputCalculator = NULL;
   fSynapses        = NULL;
   fEstimatorHistTrain = NULL;
   fEstimatorHistTest  = NULL;

   // reset monitorign histogram vectors
   fEpochMonHistS.clear();
   fEpochMonHistB.clear();
   fEpochMonHistW.clear();

   // these will be set in BuildNetwork()
   fInputLayer = NULL;
   fOutputNeurons.clear();

   frgen = new TRandom3(fRandomSeed);

   fSynapses = new TObjArray();
}

//______________________________________________________________________________
TMVA::MethodANNBase::~MethodANNBase()
{
   // destructor
   DeleteNetwork();
}

//______________________________________________________________________________
void TMVA::MethodANNBase::DeleteNetwork()
{
   // delete/clear network
   if (fNetwork != NULL) {
      TObjArray *layer;
      Int_t numLayers = fNetwork->GetEntriesFast();
      for (Int_t i = 0; i < numLayers; i++) {
         layer = (TObjArray*)fNetwork->At(i);
         DeleteNetworkLayer(layer);
      }
      delete fNetwork;
   }

   if (frgen != NULL)            delete frgen;
   if (fActivation != NULL)      delete fActivation;
   if (fOutput != NULL)          delete fOutput;  //zjh
   if (fIdentity != NULL)        delete fIdentity;
   if (fInputCalculator != NULL) delete fInputCalculator;
   if (fSynapses != NULL)        delete fSynapses;

   fNetwork         = NULL;
   frgen            = NULL;
   fActivation      = NULL;
   fOutput          = NULL; //zjh
   fIdentity        = NULL;
   fInputCalculator = NULL;
   fSynapses        = NULL;
}

//______________________________________________________________________________
void TMVA::MethodANNBase::DeleteNetworkLayer( TObjArray*& layer )
{
   // delete a network layer
   TNeuron* neuron;
   Int_t numNeurons = layer->GetEntriesFast();
   for (Int_t i = 0; i < numNeurons; i++) {
      neuron = (TNeuron*)layer->At(i);
      neuron->DeletePreLinks();
      delete neuron;
   }
   delete layer;
}

//______________________________________________________________________________
void TMVA::MethodANNBase::BuildNetwork( std::vector<Int_t>* layout, std::vector<Double_t>* weights, Bool_t fromFile )
{
   // build network given a layout (number of neurons in each layer)
   // and optional weights array

   if (fEstimatorS == "MSE")  fEstimator = kMSE;    //zjh
   else if (fEstimatorS == "CE")    fEstimator = kCE;      //zjh
   else Log()<<kWARNING<<"fEstimator="<<fEstimator<<"\tfEstimatorS="<<fEstimatorS<<Endl;
   if (fEstimator!=kMSE && fEstimator!=kCE) Log()<<kWARNING<<"Estimator type unspecified \t"<<Endl; //zjh

   Log() << kINFO << "Building Network" << Endl;

   DeleteNetwork();
   InitANNBase();

   // set activation and input functions
   TActivationChooser aChooser;
   fActivation = aChooser.CreateActivation(fNeuronType);
   fIdentity   = aChooser.CreateActivation("linear");
   if (fEstimator==kMSE)  fOutput = aChooser.CreateActivation("linear");  //zjh
   else if (fEstimator==kCE)   fOutput = aChooser.CreateActivation("sigmoid"); //zjh
   TNeuronInputChooser iChooser;
   fInputCalculator = iChooser.CreateNeuronInput(fNeuronInputType);

   fNetwork = new TObjArray();
   fRegulatorIdx.clear();		//zjh
   fRegulators.clear();			//zjh
   BuildLayers( layout, fromFile );

   // cache input layer and output neuron for fast access
   fInputLayer   = (TObjArray*)fNetwork->At(0);
   TObjArray* outputLayer = (TObjArray*)fNetwork->At(fNetwork->GetEntriesFast()-1);
   fOutputNeurons.clear();
   for (Int_t i = 0; i < outputLayer->GetEntries(); i++) {
      fOutputNeurons.push_back( (TNeuron*)outputLayer->At(i) );
   }

   if (weights == NULL) InitWeights();
   else                 ForceWeights(weights);
}




//______________________________________________________________________________
void TMVA::MethodANNBase::BuildLayers( std::vector<Int_t>* layout, Bool_t fromFile )
{
   // build the network layers

   TObjArray* curLayer;
   TObjArray* prevLayer = NULL;

   Int_t numLayers = layout->size();

   for (Int_t i = 0; i < numLayers; i++) {
      curLayer = new TObjArray();
      BuildLayer(layout->at(i), curLayer, prevLayer, i, numLayers, fromFile);
      prevLayer = curLayer;
      fNetwork->Add(curLayer);
   }

   // cache pointers to synapses for fast access, the order matters
   for (Int_t i = 0; i < numLayers; i++) {
      TObjArray* layer = (TObjArray*)fNetwork->At(i);
      Int_t numNeurons = layer->GetEntriesFast();
      if (i!=0 && i!=numLayers-1) fRegulators.push_back(0.);  //zjh
      for (Int_t j = 0; j < numNeurons; j++) {
         if (i==0) fRegulators.push_back(0.);			//zjh
         TNeuron* neuron = (TNeuron*)layer->At(j);
         Int_t numSynapses = neuron->NumPostLinks();
         for (Int_t k = 0; k < numSynapses; k++) {
            TSynapse* synapse = neuron->PostLinkAt(k);
            fSynapses->Add(synapse);
            fRegulatorIdx.push_back(fRegulators.size()-1);	//zjh
         }
      }
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::BuildLayer( Int_t numNeurons, TObjArray* curLayer, 
                                      TObjArray* prevLayer, Int_t layerIndex, 
                                      Int_t numLayers, Bool_t fromFile )
{
   // build a single layer with neurons and synapses connecting this
   // layer to the previous layer

   TNeuron* neuron;
   for (Int_t j = 0; j < numNeurons; j++) {
      if (fromFile && (layerIndex != numLayers-1) && (j==numNeurons-1)){
         neuron = new TNeuron();
         neuron->SetActivationEqn(fIdentity);
         neuron->SetBiasNeuron();
         neuron->ForceValue(1.0);
         curLayer->Add(neuron);  
      }
      else {
         neuron = new TNeuron();
         neuron->SetInputCalculator(fInputCalculator);
         
         // input layer
         if (layerIndex == 0) {
            neuron->SetActivationEqn(fIdentity);
            neuron->SetInputNeuron();
         }
         else {
            // output layer
            if (layerIndex == numLayers-1) {
               neuron->SetOutputNeuron();
               neuron->SetActivationEqn(fOutput);     //zjh
         }
            // hidden layers
            else neuron->SetActivationEqn(fActivation);
            AddPreLinks(neuron, prevLayer);
         }
         
         curLayer->Add(neuron);
      }
   }

   // add bias neutron (except to output layer)
   if(!fromFile){
      if (layerIndex != numLayers-1) {
         neuron = new TNeuron();
         neuron->SetActivationEqn(fIdentity);
         neuron->SetBiasNeuron();
         neuron->ForceValue(1.0);
         curLayer->Add(neuron);
      }
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::AddPreLinks(TNeuron* neuron, TObjArray* prevLayer)
{
   // add synapses connecting a neuron to its preceding layer

   TSynapse* synapse;
   int numNeurons = prevLayer->GetEntriesFast();
   TNeuron* preNeuron;

   for (Int_t i = 0; i < numNeurons; i++) {
      preNeuron = (TNeuron*)prevLayer->At(i);
      synapse = new TSynapse();
      synapse->SetPreNeuron(preNeuron);
      synapse->SetPostNeuron(neuron);
      preNeuron->AddPostLink(synapse);
      neuron->AddPreLink(synapse);
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::InitWeights()
{
   // initialize the synapse weights randomly
   PrintMessage("Initializing weights");
   
   // init synapse weights
   Int_t numSynapses = fSynapses->GetEntriesFast();
   TSynapse* synapse;
   for (Int_t i = 0; i < numSynapses; i++) {
      synapse = (TSynapse*)fSynapses->At(i);
      synapse->SetWeight(4.0*frgen->Rndm() - 2.0);
   }
}

//_______________________________________________________________________
void TMVA::MethodANNBase::ForceWeights(std::vector<Double_t>* weights)
{
   // force the synapse weights
   PrintMessage("Forcing weights");

   Int_t numSynapses = fSynapses->GetEntriesFast();
   TSynapse* synapse;
   for (Int_t i = 0; i < numSynapses; i++) {
      synapse = (TSynapse*)fSynapses->At(i);
      synapse->SetWeight(weights->at(i));
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::ForceNetworkInputs( const Event* ev, Int_t ignoreIndex)
{
   // force the input values of the input neurons
   // force the value for each input neuron

   Double_t x;
   TNeuron* neuron;

   //   const Event* ev = GetEvent();
   for (UInt_t j = 0; j < GetNvar(); j++) {

      x = (j != (UInt_t)ignoreIndex)?ev->GetValue(j):0;

      neuron = GetInputNeuron(j);
      neuron->ForceValue(x);
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::ForceNetworkCalculations()
{
   // calculate input values to each neuron

   TObjArray* curLayer;
   TNeuron* neuron;
   Int_t numLayers = fNetwork->GetEntriesFast();
   Int_t numNeurons;

   for (Int_t i = 0; i < numLayers; i++) {
      curLayer = (TObjArray*)fNetwork->At(i);
      numNeurons = curLayer->GetEntriesFast();

      for (Int_t j = 0; j < numNeurons; j++) {
         neuron = (TNeuron*) curLayer->At(j);
         neuron->CalculateValue();
         neuron->CalculateActivationValue();

      }
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::PrintMessage(TString message, Bool_t force) const
{
   // print messages, turn off printing by setting verbose and debug flag appropriately
   if (Verbose() || Debug() || force) Log() << kINFO << message << Endl;
}

//______________________________________________________________________________
void TMVA::MethodANNBase::WaitForKeyboard()
{
   // wait for keyboard input, for debugging
   std::string dummy;
   Log() << kINFO << "***Type anything to continue (q to quit): ";
   std::getline(std::cin, dummy);
   if (dummy == "q" || dummy == "Q") {
      PrintMessage( "quit" );
      delete this;
      exit(0);
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::PrintNetwork() const
{
   // print network representation, for debugging
   if (!Debug()) return;

   Log() << kINFO << Endl;
   PrintMessage( "Printing network " );
   Log() << kINFO << "-------------------------------------------------------------------" << Endl;

   TObjArray* curLayer;
   Int_t numLayers = fNetwork->GetEntriesFast();

   for (Int_t i = 0; i < numLayers; i++) {

      curLayer = (TObjArray*)fNetwork->At(i);
      Int_t numNeurons = curLayer->GetEntriesFast();

      Log() << kINFO << "Layer #" << i << " (" << numNeurons << " neurons):" << Endl;
      PrintLayer( curLayer );
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::PrintLayer(TObjArray* layer) const
{
   // print a single layer, for debugging

   Int_t numNeurons = layer->GetEntriesFast();
   TNeuron* neuron;
  
   for (Int_t j = 0; j < numNeurons; j++) {
      neuron = (TNeuron*) layer->At(j);
      Log() << kINFO << "\tNeuron #" << j << " (LinksIn: " << neuron->NumPreLinks() 
              << " , LinksOut: " << neuron->NumPostLinks() << ")" << Endl;
      PrintNeuron( neuron );
   }
}

//______________________________________________________________________________
void TMVA::MethodANNBase::PrintNeuron(TNeuron* neuron) const
{
   // print a neuron, for debugging
   Log() << kINFO 
           << "\t\tValue:\t"     << neuron->GetValue()
           << "\t\tActivation: " << neuron->GetActivationValue()
           << "\t\tDelta: "      << neuron->GetDelta() << Endl;
   Log() << kINFO << "\t\tActivationEquation:\t";
   neuron->PrintActivationEqn();
   Log() << kINFO << "\t\tLinksIn:" << Endl;
   neuron->PrintPreLinks();
   Log() << kINFO << "\t\tLinksOut:" << Endl;
   neuron->PrintPostLinks();
}

//_______________________________________________________________________
Double_t TMVA::MethodANNBase::GetMvaValue( Double_t* err, Double_t* errUpper )
{
   // get the mva value generated by the NN
   TNeuron* neuron;

   TObjArray* inputLayer = (TObjArray*)fNetwork->At(0);

   const Event * ev = GetEvent();

   for (UInt_t i = 0; i < GetNvar(); i++) {
      neuron = (TNeuron*)inputLayer->At(i);
      neuron->ForceValue( ev->GetValue(i) );
   }
   ForceNetworkCalculations();

   // check the output of the network
   TObjArray* outputLayer = (TObjArray*)fNetwork->At( fNetwork->GetEntriesFast()-1 );
   neuron = (TNeuron*)outputLayer->At(0);

   // cannot determine error
   NoErrorCalc(err, errUpper);

   return neuron->GetActivationValue();
}

//_______________________________________________________________________
const std::vector<Float_t> &TMVA::MethodANNBase::GetRegressionValues() 
{
   // get the regression value generated by the NN
   TNeuron* neuron;

   TObjArray* inputLayer = (TObjArray*)fNetwork->At(0);

   const Event * ev = GetEvent();

   for (UInt_t i = 0; i < GetNvar(); i++) {
      neuron = (TNeuron*)inputLayer->At(i);
      neuron->ForceValue( ev->GetValue(i) );
   }
   ForceNetworkCalculations();

   // check the output of the network
   TObjArray* outputLayer = (TObjArray*)fNetwork->At( fNetwork->GetEntriesFast()-1 );

   if (fRegressionReturnVal == NULL) fRegressionReturnVal = new std::vector<Float_t>();
   fRegressionReturnVal->clear();

   Event * evT = new Event(*ev);
   UInt_t ntgts = outputLayer->GetEntriesFast();
   for (UInt_t itgt = 0; itgt < ntgts; itgt++) {
      evT->SetTarget(itgt,((TNeuron*)outputLayer->At(itgt))->GetActivationValue());
   }

   const Event* evT2 = GetTransformationHandler().InverseTransform( evT );
   for (UInt_t itgt = 0; itgt < ntgts; itgt++) {
      fRegressionReturnVal->push_back( evT2->GetTarget(itgt) );
   }

   delete evT;

   return *fRegressionReturnVal;
}









//_______________________________________________________________________
const std::vector<Float_t> &TMVA::MethodANNBase::GetMulticlassValues()
{
   // get the multiclass classification values generated by the NN
   TNeuron* neuron;

   TObjArray* inputLayer = (TObjArray*)fNetwork->At(0);

   const Event * ev = GetEvent();
   
   for (UInt_t i = 0; i < GetNvar(); i++) {
      neuron = (TNeuron*)inputLayer->At(i);
      neuron->ForceValue( ev->GetValue(i) );
   }
   ForceNetworkCalculations();

   // check the output of the network
 
   if (fMulticlassReturnVal == NULL) fMulticlassReturnVal = new std::vector<Float_t>();
   fMulticlassReturnVal->clear();
   std::vector<Float_t> temp;

   UInt_t nClasses = DataInfo().GetNClasses();
   for (UInt_t icls = 0; icls < nClasses; icls++) {
      temp.push_back(GetOutputNeuron( icls )->GetActivationValue() );
   }
   
   for(UInt_t iClass=0; iClass<nClasses; iClass++){
      Double_t norm = 0.0;
      for(UInt_t j=0;j<nClasses;j++){
         if(iClass!=j)
            norm+=exp(temp[j]-temp[iClass]);
         }
      (*fMulticlassReturnVal).push_back(1.0/(1.0+norm));
   }


   
   return *fMulticlassReturnVal;
}


//_______________________________________________________________________
void TMVA::MethodANNBase::AddWeightsXMLTo( void* parent ) const 
{
   // create XML description of ANN classifier
   Int_t numLayers = fNetwork->GetEntriesFast();
   void* wght = gTools().xmlengine().NewChild(parent, 0, "Weights");
   void* xmlLayout = gTools().xmlengine().NewChild(wght, 0, "Layout");
   gTools().xmlengine().NewAttr(xmlLayout, 0, "NLayers", gTools().StringFromInt(fNetwork->GetEntriesFast()) );
   TString weights = "";
   for (Int_t i = 0; i < numLayers; i++) {
      TObjArray* layer = (TObjArray*)fNetwork->At(i);
      Int_t numNeurons = layer->GetEntriesFast();
      void* layerxml = gTools().xmlengine().NewChild(xmlLayout, 0, "Layer");
      gTools().xmlengine().NewAttr(layerxml, 0, "Index",    gTools().StringFromInt(i) );
      gTools().xmlengine().NewAttr(layerxml, 0, "NNeurons", gTools().StringFromInt(numNeurons) );
      for (Int_t j = 0; j < numNeurons; j++) {
         TNeuron* neuron = (TNeuron*)layer->At(j);
         Int_t numSynapses = neuron->NumPostLinks();
         void* neuronxml = gTools().AddChild(layerxml, "Neuron");
         gTools().AddAttr(neuronxml, "NSynapses", gTools().StringFromInt(numSynapses) );
         if(numSynapses==0) continue;
	     std::stringstream s("");
         s.precision( 16 );
         for (Int_t k = 0; k < numSynapses; k++) {
            TSynapse* synapse = neuron->PostLinkAt(k);
            s << std::scientific << synapse->GetWeight() << " ";
         }
         gTools().AddRawLine( neuronxml, s.str().c_str() );
      }
   }

   // if inverse hessian exists, write inverse hessian to weight file
   if( fInvHessian.GetNcols()>0 ){
      void* xmlInvHessian = gTools().xmlengine().NewChild(wght, 0, "InverseHessian");

      // get the matrix dimensions
      Int_t nElements = fInvHessian.GetNoElements();
      Int_t nRows     = fInvHessian.GetNrows();
      Int_t nCols     = fInvHessian.GetNcols();
      gTools().xmlengine().NewAttr(xmlInvHessian, 0, "NElements", gTools().StringFromInt(nElements) );
      gTools().xmlengine().NewAttr(xmlInvHessian, 0, "NRows", gTools().StringFromInt(nRows) );
      gTools().xmlengine().NewAttr(xmlInvHessian, 0, "NCols", gTools().StringFromInt(nCols) );

      // read in the matrix elements
      Double_t* elements = new Double_t[nElements+10];
      fInvHessian.GetMatrix2Array( elements );

      // store the matrix elements row-wise
      Int_t index = 0;
      for( Int_t row = 0; row < nRows; ++row ){
         void* xmlRow = gTools().xmlengine().NewChild(xmlInvHessian, 0, "Row");
         gTools().xmlengine().NewAttr(xmlRow, 0, "Index", gTools().StringFromInt(row) );

         // create the rows
	     std::stringstream s("");
         s.precision( 16 );
         for( Int_t col = 0; col < nCols; ++col ){
            s << std::scientific << (*(elements+index)) << " ";
            ++index;
         }
         gTools().xmlengine().AddRawLine( xmlRow, s.str().c_str() );
      }
      delete[] elements;
   }
}


//_______________________________________________________________________
void TMVA::MethodANNBase::ReadWeightsFromXML( void* wghtnode )
{
   // read MLP from xml weight file

   // build the layout first
   Bool_t fromFile = kTRUE;
   std::vector<Int_t>* layout = new std::vector<Int_t>();

   void* xmlLayout = NULL;
   xmlLayout = gTools().GetChild(wghtnode, "Layout");
   if( !xmlLayout )
      xmlLayout = wghtnode;

   UInt_t nLayers;
   gTools().ReadAttr( xmlLayout, "NLayers", nLayers );
   layout->resize( nLayers );

   void* ch = gTools().xmlengine().GetChild(xmlLayout);
   UInt_t index;
   UInt_t nNeurons;
   while (ch) {
      gTools().ReadAttr( ch, "Index",   index   );
      gTools().ReadAttr( ch, "NNeurons", nNeurons );
      layout->at(index) = nNeurons;
      ch = gTools().GetNextChild(ch);
   }

   BuildNetwork( layout, NULL, fromFile );
   // fill the weights of the synapses
   UInt_t nSyn;
   Float_t weight;
   ch = gTools().xmlengine().GetChild(xmlLayout);
   UInt_t iLayer = 0;
   while (ch) {  // layers
      TObjArray* layer = (TObjArray*)fNetwork->At(iLayer);
      gTools().ReadAttr( ch, "Index",   index   );
      gTools().ReadAttr( ch, "NNeurons", nNeurons );

      void* nodeN = gTools().GetChild(ch);
      UInt_t iNeuron = 0;
      while( nodeN ){ // neurons
         TNeuron *neuron = (TNeuron*)layer->At(iNeuron);
         gTools().ReadAttr( nodeN, "NSynapses", nSyn );
         if( nSyn > 0 ){
            const char* content = gTools().GetContent(nodeN);
            std::stringstream s(content);
            for (UInt_t iSyn = 0; iSyn<nSyn; iSyn++) { // synapses

               TSynapse* synapse = neuron->PostLinkAt(iSyn);
               s >> weight;
               //Log() << kWARNING << neuron << " " << weight <<  Endl;
               synapse->SetWeight(weight);
            }
         }
         nodeN = gTools().GetNextChild(nodeN);
         iNeuron++;
      }
      ch = gTools().GetNextChild(ch);
      iLayer++;
   }

   delete layout;

   void* xmlInvHessian = NULL;
   xmlInvHessian = gTools().GetChild(wghtnode, "InverseHessian");
   if( !xmlInvHessian )
      // no inverse hessian available
      return;

   fUseRegulator = kTRUE;

   Int_t nElements = 0;
   Int_t nRows     = 0;
   Int_t nCols     = 0;
   gTools().ReadAttr( xmlInvHessian, "NElements", nElements );
   gTools().ReadAttr( xmlInvHessian, "NRows", nRows );
   gTools().ReadAttr( xmlInvHessian, "NCols", nCols );

   // adjust the matrix dimensions
   fInvHessian.ResizeTo( nRows, nCols );

   // prepare an array to read in the values
   Double_t* elements;
   if (nElements > std::numeric_limits<int>::max()-100){
      Log() << kFATAL << "you tried to read a hessian matrix with " << nElements << " elements, --> too large, guess s.th. went wrong reading from the weight file" << Endl;
      return;
   } else {
      elements = new Double_t[nElements+10];
   }     



   void* xmlRow = gTools().xmlengine().GetChild(xmlInvHessian);
   Int_t row = 0;
   index = 0;
   while (xmlRow) {  // rows
      gTools().ReadAttr( xmlRow, "Index",   row   );

      const char* content = gTools().xmlengine().GetNodeContent(xmlRow);

      std::stringstream s(content);
      for (Int_t iCol = 0; iCol<nCols; iCol++) { // columns
	 s >> (*(elements+index));
	 ++index;
      }
      xmlRow = gTools().xmlengine().GetNext(xmlRow);
      ++row;
   }

   fInvHessian.SetMatrixArray( elements );

   delete[] elements;
}


//_______________________________________________________________________
void TMVA::MethodANNBase::ReadWeightsFromStream( std::istream & istr)
{
   // destroy/clear the network then read it back in from the weights file

   // delete network so we can reconstruct network from scratch

   TString dummy;

   // synapse weights
   Double_t weight;
   std::vector<Double_t>* weights = new std::vector<Double_t>();
   istr>> dummy;
   while (istr>> dummy >> weight) weights->push_back(weight); // use w/ slower write-out

   ForceWeights(weights);
   

   delete weights;
}

//_______________________________________________________________________
const TMVA::Ranking* TMVA::MethodANNBase::CreateRanking()
{
   // compute ranking of input variables by summing function of weights

   // create the ranking object
   fRanking = new Ranking( GetName(), "Importance" );

   TNeuron*  neuron;
   TSynapse* synapse;
   Double_t  importance, avgVal;
   TString varName;

   for (UInt_t ivar = 0; ivar < GetNvar(); ivar++) {

      neuron = GetInputNeuron(ivar);
      Int_t numSynapses = neuron->NumPostLinks();
      importance = 0;
      varName = GetInputVar(ivar); // fix this line

      // figure out average value of variable i
      Double_t meanS, meanB, rmsS, rmsB, xmin, xmax;
      Statistics( TMVA::Types::kTraining, varName, 
                  meanS, meanB, rmsS, rmsB, xmin, xmax );

      avgVal = (TMath::Abs(meanS) + TMath::Abs(meanB))/2.0;
      double meanrms = (TMath::Abs(rmsS) + TMath::Abs(rmsB))/2.;
      if (avgVal<meanrms) avgVal = meanrms;      
      if (IsNormalised()) avgVal = 0.5*(1 + gTools().NormVariable( avgVal, GetXmin( ivar ), GetXmax( ivar ))); 

      for (Int_t j = 0; j < numSynapses; j++) {
         synapse = neuron->PostLinkAt(j);
         importance += synapse->GetWeight() * synapse->GetWeight();
      }
      
      importance *= avgVal * avgVal;

      fRanking->AddRank( Rank( varName, importance ) );
   }

   return fRanking;
}

//_______________________________________________________________________
void TMVA::MethodANNBase::CreateWeightMonitoringHists( const TString& bulkname, 
                                                       std::vector<TH1*>* hv ) const
{
   TH2F* hist;
   Int_t numLayers = fNetwork->GetEntriesFast();

   for (Int_t i = 0; i < numLayers-1; i++) {
      
      TObjArray* layer1 = (TObjArray*)fNetwork->At(i);
      TObjArray* layer2 = (TObjArray*)fNetwork->At(i+1);
      Int_t numNeurons1 = layer1->GetEntriesFast();
      Int_t numNeurons2 = layer2->GetEntriesFast();
      
      TString name = Form("%s%i%i", bulkname.Data(), i, i+1);
      hist = new TH2F(name + "", name + "", 
                      numNeurons1, 0, numNeurons1, numNeurons2, 0, numNeurons2);
      
      for (Int_t j = 0; j < numNeurons1; j++) {
         
         TNeuron* neuron = (TNeuron*)layer1->At(j);
         Int_t numSynapses = neuron->NumPostLinks();

         for (Int_t k = 0; k < numSynapses; k++) {
            
            TSynapse* synapse = neuron->PostLinkAt(k);
            hist->SetBinContent(j+1, k+1, synapse->GetWeight());
            
         }
      }
      
      if (hv) hv->push_back( hist );
      else {
         hist->Write();
         delete hist;
      }
   }
}

//_______________________________________________________________________
void TMVA::MethodANNBase::WriteMonitoringHistosToFile() const
{
   // write histograms to file
   PrintMessage(Form("Write special histos to file: %s", BaseDir()->GetPath()), kTRUE);

   if (fEstimatorHistTrain) fEstimatorHistTrain->Write();
   if (fEstimatorHistTest ) fEstimatorHistTest ->Write();

   // histograms containing weights for architecture plotting (used in macro "network.C")
   CreateWeightMonitoringHists( "weights_hist" );

   // now save all the epoch-wise monitoring information
#if __cplusplus > 199711L
   static std::atomic<int> epochMonitoringDirectoryNumber{0};
#else
   static int epochMonitoringDirectoryNumber = 0;
#endif
   int epochVal = epochMonitoringDirectoryNumber++;
   TDirectory* epochdir = NULL;
   if( epochVal == 0 )
      epochdir = BaseDir()->mkdir( "EpochMonitoring" );
   else
      epochdir = BaseDir()->mkdir( Form("EpochMonitoring_%4d",epochVal) );

   epochdir->cd();
   for (std::vector<TH1*>::const_iterator it = fEpochMonHistS.begin(); it != fEpochMonHistS.end(); it++) {
      (*it)->Write();
      delete (*it);
   }
   for (std::vector<TH1*>::const_iterator it = fEpochMonHistB.begin(); it != fEpochMonHistB.end(); it++) {
      (*it)->Write();
      delete (*it);
   }
   for (std::vector<TH1*>::const_iterator it = fEpochMonHistW.begin(); it != fEpochMonHistW.end(); it++) {
      (*it)->Write();
      delete (*it);
   }
   BaseDir()->cd();
}

//_______________________________________________________________________
void TMVA::MethodANNBase::MakeClassSpecific( std::ostream& fout, const TString& className ) const
{
   // write specific classifier response
   Int_t numLayers = fNetwork->GetEntries();

   fout << std::endl;
   fout << "   double ActivationFnc(double x) const;" << std::endl;
   fout << "   double OutputActivationFnc(double x) const;" << std::endl;     //zjh
   fout << std::endl;
   fout << "   int fLayers;" << std::endl;
   fout << "   int fLayerSize["<<numLayers<<"];" << std::endl;
   int numNodesFrom = -1;
   for (Int_t lIdx = 0; lIdx < numLayers; lIdx++) {
      int numNodesTo = ((TObjArray*)fNetwork->At(lIdx))->GetEntries();
      if (numNodesFrom<0) { numNodesFrom=numNodesTo; continue; }
      fout << "   double fWeightMatrix" << lIdx-1  << "to" << lIdx << "[" << numNodesTo << "][" << numNodesFrom << "];";
      fout << "   // weight matrix from layer " << lIdx-1  << " to " << lIdx << std::endl;
      numNodesFrom = numNodesTo;
   }
   fout << std::endl;
   fout << "   double * fWeights["<<numLayers<<"];" << std::endl;
   fout << "};" << std::endl;

   fout << std::endl;

   fout << "inline void " << className << "::Initialize()" << std::endl;
   fout << "{" << std::endl;
   fout << "   // build network structure" << std::endl;
   fout << "   fLayers = " << numLayers << ";" << std::endl;
   for (Int_t lIdx = 0; lIdx < numLayers; lIdx++) {
      TObjArray* layer = (TObjArray*)fNetwork->At(lIdx);
      int numNodes = layer->GetEntries();
      fout << "   fLayerSize[" << lIdx << "] = " << numNodes << "; fWeights["<<lIdx<<"] = new double["<<numNodes<<"]; " << std::endl;
   }

   for (Int_t i = 0; i < numLayers-1; i++) {
      fout << "   // weight matrix from layer " << i  << " to " << i+1 << std::endl;
      TObjArray* layer = (TObjArray*)fNetwork->At(i);
      Int_t numNeurons = layer->GetEntriesFast();
      for (Int_t j = 0; j < numNeurons; j++) {
         TNeuron* neuron = (TNeuron*)layer->At(j);
         Int_t numSynapses = neuron->NumPostLinks();
         for (Int_t k = 0; k < numSynapses; k++) {
            TSynapse* synapse = neuron->PostLinkAt(k);
            fout << "   fWeightMatrix" << i  << "to" << i+1 << "[" << k << "][" << j << "] = " << synapse->GetWeight() << ";" << std::endl;
         }
      }
   }

   fout << "}" << std::endl;
   fout << std::endl;

   // writing of the GetMvaValue__ method
   fout << "inline double " << className << "::GetMvaValue__( const std::vector<double>& inputValues ) const" << std::endl;
   fout << "{" << std::endl;
   fout << "   if (inputValues.size() != (unsigned int)fLayerSize[0]-1) {" << std::endl;
   fout << "      std::cout << \"Input vector needs to be of size \" << fLayerSize[0]-1 << std::endl;" << std::endl;
   fout << "      return 0;" << std::endl;
   fout << "   }" << std::endl;
   fout << std::endl;
   fout << "   for (int l=0; l<fLayers; l++)" << std::endl;
   fout << "      for (int i=0; i<fLayerSize[l]; i++) fWeights[l][i]=0;" << std::endl;
   fout << std::endl;
   fout << "   for (int l=0; l<fLayers-1; l++)" << std::endl;
   fout << "      fWeights[l][fLayerSize[l]-1]=1;" << std::endl;
   fout << std::endl;
   fout << "   for (int i=0; i<fLayerSize[0]-1; i++)" << std::endl;
   fout << "      fWeights[0][i]=inputValues[i];" << std::endl;
   fout << std::endl;
   for (Int_t i = 0; i < numLayers-1; i++) {
      fout << "   // layer " << i << " to " << i+1 << std::endl;
      if (i+1 == numLayers-1) {
         fout << "   for (int o=0; o<fLayerSize[" << i+1 << "]; o++) {" << std::endl;
      } 
      else {
         fout << "   for (int o=0; o<fLayerSize[" << i+1 << "]-1; o++) {" << std::endl;
      }
      fout << "      for (int i=0; i<fLayerSize[" << i << "]; i++) {" << std::endl;
      fout << "         double inputVal = fWeightMatrix" << i << "to" << i+1 << "[o][i] * fWeights[" << i << "][i];" << std::endl;

      if ( fNeuronInputType == "sum") {
         fout << "         fWeights[" << i+1 << "][o] += inputVal;" << std::endl;
      } 
      else if ( fNeuronInputType == "sqsum") {
         fout << "         fWeights[" << i+1 << "][o] += inputVal*inputVal;" << std::endl;
      } 
      else { // fNeuronInputType == TNeuronInputChooser::kAbsSum
         fout << "         fWeights[" << i+1 << "][o] += fabs(inputVal);" << std::endl;
      }
      fout << "      }" << std::endl;
      if (i+1 != numLayers-1) // in the last layer no activation function is applied
         fout << "      fWeights[" << i+1 << "][o] = ActivationFnc(fWeights[" << i+1 << "][o]);" << std::endl;
      else	fout << "      fWeights[" << i+1 << "][o] = OutputActivationFnc(fWeights[" << i+1 << "][o]);" << std::endl; //zjh
      fout << "   }" << std::endl;
   }
   fout << std::endl;
   fout << "   return fWeights[" << numLayers-1 << "][0];" << std::endl;   
   fout << "}" << std::endl;

   fout << std::endl;
   TString fncName = className+"::ActivationFnc";
   fActivation->MakeFunction(fout, fncName);
   fncName = className+"::OutputActivationFnc";  	//zjh
   fOutput->MakeFunction(fout, fncName); 			//zjh

   fout << "   " << std::endl;
   fout << "// Clean up" << std::endl;
   fout << "inline void " << className << "::Clear() " << std::endl;
   fout << "{" << std::endl;
   fout << "   // clean up the arrays" << std::endl;
   fout << "   for (int lIdx = 0; lIdx < "<<numLayers<<"; lIdx++) {" << std::endl;
   fout << "      delete[] fWeights[lIdx];" << std::endl;
   fout << "   }" << std::endl;
   fout << "}" << std::endl;
}

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