00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #include <stdlib.h>
00018 #include <time.h>
00019
00020 #include <stdio.h>
00021 #include <iostream>
00022 #include <fstream>
00023 #include <sstream>
00024 #include "HepMC/PythiaWrapper.h"
00025 #include "HepMC/IO_HEPEVT.h"
00026 #include "HepMC/IO_Ascii.h"
00027 #include "HepMC/GenEvent.h"
00028 #include "HepMC/SimpleVector.h"
00029 #include "MyPythia6Wrapper.h"
00030
00031 #include "TH1.h"
00032 #include "TFile.h"
00033 #include "TMath.h"
00034 #include <TCanvas.h>
00035
00036 #include "TF1.h"
00037 #include <TH2.h>
00038 #include <TStyle.h>
00039 #include "TLegend.h"
00040 #include <TTree.h>
00041 #include <TString.h>
00042
00043 using namespace std;
00044
00045
00046 #include "../../include/Configuration.h"
00047 #include "../../include/JetFinder.h"
00048 #include "../../include/TopAnalysis.h"
00049 #include "../../include/TauAnalysis.h"
00050 #include "../../include/DiJetAnalysis.h"
00051 #include "../../include/WplusJetAnalysis.h"
00052 #include "../../include/ZAnalysis.h"
00053 #include "../../include/UEAnalysis.h"
00054 #include "../../include/EtMissAnalysis.h"
00055 #include "../../include/UserAnalysis.h"
00056 #include "readpythiaConfigFile.h"
00057
00058
00059
00060 int main(int argc, const char* argv[]) {
00061
00062 if(argc <= 1 ) {
00063 std::cout << "Usage: " << argv[0] << " <Configuration File>" << std::endl;
00064 return 0;
00065 }
00066
00067
00068 Configuration* currentconfig = new Configuration();
00069 if(currentconfig->read(argv[1]))
00070 {
00071 std::cout << " ... skip !" << std::endl;
00072 return 0;
00073 }
00074
00075 JetFinder FindJet;
00076 FindJet.Init(currentconfig->max_eta,currentconfig->min_pt);
00077 FindJet.InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00078
00079
00080 std::vector<baseAnalysis*> analysis;
00081
00082
00083 if (currentconfig->dijet) {
00084 std::cout << "Adding module DiJetAnalysis" << std::endl;
00085 analysis.push_back(new DiJetAnalysis);
00086 }
00087 if (currentconfig->wplusjet) {
00088 std::cout << "Adding module WplusJetAnalysis" << std::endl;
00089 analysis.push_back(new WplusJetAnalysis);
00090 }
00091 if(currentconfig->z){
00092 std::cout << "Adding module ZAnalysis" << std::endl;
00093 analysis.push_back(new ZAnalysis);
00094 }
00095 if(currentconfig->tau){
00096 std::cout << "Adding module TauAnalysis" << std::endl;
00097 analysis.push_back(new TauAnalysis);
00098 }
00099 if(currentconfig->top){
00100 std::cout << "Adding module TopAnalysis" << std::endl;
00101 analysis.push_back(new TopAnalysis);
00102 }
00103 if(currentconfig->ue){
00104 std::cout << "Adding module UEAnalysis" << std::endl;
00105 analysis.push_back(new UEAnalysis);
00106 }
00107 if(currentconfig->etmiss){
00108 std::cout << "Adding module EtMissAnalysis" << std::endl;
00109 analysis.push_back(new EtMissAnalysis);
00110 }
00111 if(currentconfig->user){
00112 std::cout << "Adding module UserAnalysis" << std::endl;
00113 analysis.push_back(new UserAnalysis);
00114 }
00115
00116
00117
00118 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00119 (*i)->Init(currentconfig->max_eta,currentconfig->min_pt);
00120 (*i)->InitJetFinder(currentconfig->coneRadius,currentconfig->overlapThreshold,currentconfig->jet_ptmin,currentconfig->lepton_ptmin,currentconfig->DeltaR_lepton_track);
00121 }
00122
00123
00124
00125
00126
00127 pyjets.n=0;
00128 call_pyhepc(1);
00129
00130 HepMC::HEPEVT_Wrapper::set_max_number_entries(10000);
00131 HepMC::HEPEVT_Wrapper::set_sizeof_real(8);
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142 pydatr.mrpy[1-1] = currentconfig->rseed ;
00143
00144
00145
00146
00147 std::cout << "Configuration Parameter List:" << std::endl;
00148
00149 for (unsigned int files=0; files<(currentconfig->ConfigFileNames).size(); files++)
00150 readConfigFile((currentconfig->ConfigFileNames)[files]);
00151
00152
00153
00154
00155
00156
00157
00158
00159 call_pyinit("CMS", "p", "p", 14000.);
00160
00161
00162
00163
00164 HepMC::IO_HEPEVT hepevtio;
00165
00166
00167
00168
00169
00170
00171
00172 int nShow = 50;
00173 int nShowPace = max(1,int(currentconfig->nevents)/nShow);
00174
00175 for (unsigned int iEvent = 1 ; iEvent <= currentconfig->nevents; iEvent++) {
00176
00177
00178 if (iEvent%nShowPace == 0) cout << "pythia6: Now begin event "
00179 << iEvent << endl;
00180
00181 call_pyevnt();
00182
00183 call_pyhepc(1);
00184 HepMC::GenEvent* evt = hepevtio.read_next_event();
00185
00186 evt->set_event_number(iEvent);
00187 evt->set_signal_process_id(20);
00188
00189
00190
00191
00192 vector<fastjet::PseudoJet> inclusive_jets = FindJet.GetJet(evt);
00193
00194
00195 int ret_all = true;
00196 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00197 (*i)->SetJet(&inclusive_jets);
00198 int ret = (*i)->Process(evt);
00199 if (ret == false) ret_all = false;
00200 }
00201
00202 if (ret_all==false) {
00203 std::cout << "Processing event number " << iEvent << std::endl;
00204 std::cout << "throw FAILURE"<<std::endl;
00205
00206
00207
00208
00209
00210 std::cout<<std::endl;
00211 }
00212
00213
00214 FindJet.ClearEvent(evt);
00215
00216
00217 delete evt;
00218 }
00219
00220
00221
00222 call_pystat(1);
00223
00224 TFile* f = new TFile(currentconfig->OutputFileName.c_str(),"RECREATE");
00225 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00226 (*i)->Finalize(f);
00227 }
00228
00229
00230 for (std::vector<baseAnalysis*>::const_iterator i(analysis.begin()); i != analysis.end(); ++i) {
00231 delete (*i);
00232 }
00233 analysis.resize(0);
00234
00235
00236 delete currentconfig;
00237
00238 std::cout<<"pythia6: program ends"<<std::endl;
00239 f->Close();
00240
00241 std::cout<<"Run successfully!"<<std::endl;
00242 return 0;
00243
00244 }
00245