Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members | Related Pages

RivetInterface.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Rivet.hh"
00003 #include "Rivet/AnalysisHandler.hh"
00004 #include "Rivet/Tools/Commandline.hh"
00005 #include "Rivet/Tools/Logging.hh"
00006 using namespace Rivet;
00007 using namespace TCLAP;
00008 
00009 #include "HepMC/IO_Ascii.h"
00010 #include "HepMC/GenEvent.h"
00011 using namespace HepMC;
00012 
00013 
00014 int main(int argc, char* argv[]) {
00015 
00016   // Configuration variables
00017   set<string> cfgAnalyses;
00018   string cfgHistoFileName;
00019   HistoFormat cfgHistoFormat;
00020   Log::LevelMap cfgLogLevels;
00021   unsigned int cfgNumEvents;
00022   string cfgEventFile;
00023 
00024   try {
00025     // Set up command line arguments
00026     CmdLine cmd("Analyse HepMC events from a file using the Rivet system", ' ', "0.1"); 
00027     //
00028     MultiArg<string>* analysesArg(0);
00029     SwitchArg* analysesAllArg(0);
00030     ValuesConstraint<string>* anaNameConstraint(0);
00031     Commandline::addAnalysisArgs(cmd, analysesArg, analysesAllArg, anaNameConstraint);
00032     //
00033     ValueArg<string>* histoNameArg(0);
00034     ValueArg<string>* histoTypeArg(0);
00035     ValuesConstraint<string>* histoTypeConstraint(0);
00036     Commandline::addHistoArgs(cmd, histoNameArg, histoTypeArg, histoTypeConstraint);
00037     //
00038     MultiArg<string>* logsArg(0);
00039     Commandline::addLoggingArgs(cmd, logsArg);
00040     //
00041     ValueArg<unsigned int> 
00042       numEventsArg("n", "numevents", "Max number of events to read (100000 by default)", 
00043                    false, 100000, "num", cmd);
00044     
00045     UnlabeledValueArg<string>
00046       eventFileArg("eventfile", "File containing ASCII format HepMC events", true, "-", "filename", cmd);      
00047     
00048     // Parse command line args
00049     cmd.parse(argc, argv);
00050     Commandline::useAnalysisArgs(cmd, analysesArg, analysesAllArg, cfgAnalyses);
00051     Commandline::useHistoArgs(cmd, histoNameArg, histoTypeArg, cfgHistoFileName, cfgHistoFormat);
00052     Commandline::useLoggingArgs(cmd, logsArg, cfgLogLevels);
00053     cfgNumEvents = numEventsArg.getValue();
00054     cfgEventFile = eventFileArg.getValue();
00055 
00056   } catch (ArgException& e) { 
00057     cerr << "Command line error: " << e.error() << " for arg " << e.argId() << endl; 
00058     return EXIT_FAILURE;
00059   } catch (runtime_error& e) { 
00060     cerr << "Error: " << e.what() << endl; 
00061     return EXIT_FAILURE;
00062   }
00063 
00064 
00065   /////////////////////////////////////////////////////////////////////////////
00066 
00067 
00068   // Add all the log levels from the command line into the logging framework
00069   Log::setDefaultLevels(cfgLogLevels);
00070   Log& log = Log::getLog("Rivet.Main");
00071 
00072 
00073   // Make a handler and add analyses
00074   AnalysisHandler handler;
00075   for (set<string>::const_iterator an = cfgAnalyses.begin(); an != cfgAnalyses.end(); ++an) {
00076     Analysis* a = AnalysisLoader::getAnalysis(*an);
00077     BeamPair beams = a->getBeams();
00078     log << Log::INFO << "Analysis name: " << a->getName() << endl;
00079     log << Log::DEBUG << a->getBeams() << a->isCompatible(PROTON, PROTON) << endl;
00080     log << Log::INFO << "Cuts:" << a->getCuts() << endl;
00081     handler.addAnalysis(a->getName());
00082   }
00083   handler.init();
00084 
00085 
00086   // Make a HepMC input
00087   IO_Ascii hepmcIn(cfgEventFile.c_str(), std::ios::in);
00088   if (hepmcIn.rdstate() != 0) {
00089     log << Log::ERROR << "Couldn't read HepMC events from file: " << cfgEventFile << endl;
00090   }
00091 
00092 
00093   // Event loop
00094   unsigned int num = 0;
00095   while (num < cfgNumEvents) {
00096     HepMC::GenEvent myevent;
00097     if (! hepmcIn.fill_next_event(&myevent)) break;
00098     log << Log::INFO << "Event number " << num + 1 << endl;
00099     handler.analyze(myevent);
00100     ++num;
00101   }
00102 
00103 
00104   // Finish up
00105   log << Log::INFO << "Finished!"  << endl;
00106   handler.finalize();
00107   return EXIT_SUCCESS;
00108 }