Quick start with Delphes ======================== Commands to get the code: wget --no-check-certificate https://cp3.irmp.ucl.ac.be/projects/delphes/raw-attachment/wiki/WikiStart/Delphes-3.0.0.tar.gz tar -zxf Delphes-3.0.0.tar.gz Commands to compile the code: cd Delphes-3.0.0 make -j 4 Finally, we can run Delphes: ./DelphesHepMC and simulate some Z->ee events: ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output.root examples/z_ee.hep or curl -s http://cp3.irmp.ucl.ac.be/~demin/z_ee.hep.gz | gunzip | ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output.root Simple analysis using TTree::Draw ================================= Now we can start ROOT and look at the data stored on the tree Start ROOT and load shared library: root -l gSystem->Load("libDelphes"); Open ROOT tree file and do some basic analysis using Draw or TBrowser: TFile::Open("delphes_output.root"); Delphes->Draw("Electron.PT"); TBrowser browser; Note 1: Delphes - tree name, it can be learnt e.g. from TBrowser Note 2: Electron - branch name; PT - variable (leaf) of this branch Complete description of all branches can be found in Delphes/doc/RootTreeDescription.html This file is also available via web: http://cp3.irmp.ucl.ac.be/~demin/RootTreeDescription.html Macro-based analysis ==================== Analysis macro consists of histogram booking, event loop (histogram filling), histogram display Let's simulate some bbbar events: ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output.root examples/bbar.hep or curl -s http://cp3.irmp.ucl.ac.be/~demin/bbar.hep.gz | gunzip | ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output.root Start ROOT and load shared library: root -l gSystem->Load("libDelphes"); Basic analysis macro: { // Create chain of root trees TChain chain("Delphes"); chain.Add("delphes_output.root"); // Create object of class ExRootTreeReader ExRootTreeReader *treeReader = new ExRootTreeReader(&chain); Long64_t numberOfEntries = treeReader->GetEntries(); // Get pointers to branches used in this analysis TClonesArray *branchJet = treeReader->UseBranch("Jet"); // Book histograms TH1 *histJetPT = new TH1F("jet_pt", "jet P_{T}", 100, 0.0, 100.0); // Loop over all events for(Int_t entry = 0; entry < numberOfEntries; ++entry) { // Load selected branches with data from specified event treeReader->ReadEntry(entry); // If event contains at least 1 jet if(branchJet->GetEntries() > 0) { // Take first jet Jet *jet = (Jet*) branchJet->At(0); // Plot jet transverse momentum histJetPT->Fill(jet->PT); // Print jet transverse momentum cout << jet->PT << endl; } } // Show resulting histograms histJetPT->Draw(); } More advanced macro-based analysis ================================== Delphes/examples contains macros Example1.C and Example2.C Here are commands to run these macros: root -l gSystem->Load("libDelphes"); .X examples/Example2.C("delphes_output.root"); Modifying configuration file ============================ Open default configuration file in an editor: vi examples/delphes_card_CMS.tcl Replace the default b-tagging efficiency formulas with the the new ones: # default efficiency formula (misidentification rate) add EfficiencyFormula {0} {0.001} # efficiency formula for c-jets (misidentification rate) add EfficiencyFormula {4} { (pt <= 15.0) * (0.000) + \ (abs(eta) <= 1.2) * (pt > 15.0) * (0.2*tanh(pt*0.03 - 0.4)) + \ (abs(eta) > 1.2 && abs(eta) <= 2.5) * (pt > 15.0) * (0.1*tanh(pt*0.03 - 0.4)) + \ (abs(eta) > 2.5) * (0.000)} # efficiency formula for b-jets add EfficiencyFormula {5} { (pt <= 15.0) * (0.000) + \ (abs(eta) <= 1.2) * (pt > 15.0) * (0.5*tanh(pt*0.03 - 0.4)) + \ (abs(eta) > 1.2 && abs(eta) <= 2.5) * (pt > 15.0) * (0.4*tanh(pt*0.03 - 0.4)) + \ (abs(eta) > 2.5) * (0.000)} Now we can run Delphes with the modified configuration file: ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output_new.root examples/bbar.hep or curl -s http://cp3.irmp.ucl.ac.be/~demin/bbar.hep.gz | gunzip | ./DelphesSTDHEP examples/delphes_card_CMS.tcl delphes_output_new.root Open ROOT tree file and check the PT spectrum of the b-jets: root -l TFile *file = TFile::Open("delphes_output.root"); Delphes->Draw("Jet.PT", "Jet.BTag == 1"); file->Close(); TFile *file = TFile::Open("delphes_output_new.root"); Delphes->Draw("Jet.PT", "Jet.BTag == 1", "SAME"); htemp->SetLineColor(2); htemp->Draw("SAME"); file->Close();