== Library Interface == Here is an explanation of how Delphes' library can be called from a program. == Installation == Commands to download and build the library: {{{ setup ROOT environment variables wget --no-check-certificate https://cp3.irmp.ucl.ac.be/projects/delphes/raw-attachment/wiki/WikiStart/Delphes_V_3.0.0.tar.gz tar -zxf Delphes_V_3.0.0.tar.gz cd Delphes_V_3.0.0 make -j 4 # library ls libDelphes.so # headers ls modules/Delphes.h ls classes/DelphesClasses.h ls classes/DelphesFactory.h # config file ls examples/delphes_card_CMS.tcl }}} == Complete example == A fully functional application using the Delphes library can be found in {{{ readers/DelphesSTDHEP.cpp }}} == Simplified example == {{{ReadEvent}}}, {{{ReadParticle}}}, {{{ConvertInput}}} and {{{ConvertOutput}}} functions should be implemented by the library's user. {{{ #include #include #include "TROOT.h" #include "TApplication.h" #include "TObjArray.h" #include "TDatabasePDG.h" #include "TParticlePDG.h" #include "TLorentzVector.h" #include "modules/Delphes.h" #include "classes/DelphesClasses.h" #include "classes/DelphesFactory.h" using namespace std; void ConvertInput(DelphesFactory *factory, TObjArray *particleOutputArray, TObjArray *candidateOutputArray, TObjArray *partonOutputArray); void ConvertOutput(Delphes *modularDelphes); bool ReadEvent() { return true; } bool ReadParticle() { return true; } int main() { // Declaration of variables ExRootConfReader *confReader; Delphes *modularDelphes; DelphesFactory *factory; TObjArray *particleOutputArray; TObjArray *candidateOutputArray; TObjArray *partonOutputArray; gROOT->SetBatch(); int appargc = 1; char *appName = "DelphesExample"; char *appargv[] = {appName}; TApplication app(appName, &appargc, appargv); try { // Initialization confReader = new ExRootConfReader; confReader->ReadFile("examples/delphes_card_CMS.tcl"); modularDelphes = new Delphes("Delphes"); modularDelphes->SetConfReader(confReader); factory = modularDelphes->GetFactory(); particleOutputArray = modularDelphes->ExportArray("particles"); candidateOutputArray = modularDelphes->ExportArray("candidates"); partonOutputArray = modularDelphes->ExportArray("partons"); modularDelphes->InitTask(); // Event loop while(ReadEvent()) { modularDelphes->Clear(); ConvertInput(factory, particleOutputArray, candidateOutputArray, partonOutputArray); modularDelphes->ProcessTask(); ConvertOutput(modularDelphes); } // Finalization modularDelphes->FinishTask(); delete modularDelphes; delete confReader; return 0; } catch(runtime_error &e) { cerr << "** ERROR: " << e.what() << endl; return 1; } } void ConvertInput(DelphesFactory *factory, TObjArray *particleOutputArray, TObjArray *candidateOutputArray, TObjArray *partonOutputArray) { Candidate *candidate; TDatabasePDG *pdg; TParticlePDG *pdgParticle; Int_t pid, status; Double_t px, py, pz, e; Double_t x, y, z, t; pdg = TDatabasePDG::Instance() while(ReadParticle()) { if(status == 1 || status == 2) { candidate = factory->NewCandidate(); candidate->PID = pid; candidate->Status = status; pdgParticle = pdg->GetParticle(fPID); candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999; candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9; candidate->Momentum.SetPxPyPzE(px, py, pz, e); candidate->Position.SetXYZT(x, y, z, t); particleOutputArray->Add(candidate); if(!pdgParticle) return; if(status == 1) { candidateOutputArray->Add(candidate); } else if(status == 2) { partonOutputArray->Add(candidate); } } } } void ConvertOutput(Delphes *modularDelphes) { const TObjArray *arrayJets = modularDelphes->ImportArray("FastJetFinder/jets"); TIter iteratorJets(arrayJets); Candidate *candidate; iteratorJets.Reset(); while((candidate = static_cast(iteratorJets.Next()))) { const TLorentzVector &momentum = candidate->GetMomentum(); } } }}}