Fork me on GitHub

source: git/modules/LLPFilter.cc@ d1ab205

Last change on this file since d1ab205 was d1ab205, checked in by christinaw97 <christina.wang@…>, 3 years ago

updated comments

  • Property mode set to 100644
File size: 6.5 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19/** \class LLPFilter
20 *
21 * Filter LLPs with particular PDG ID/status and calculate the EM and hadronic energy of LLP based on decay particles
22 * The classification of EM and hadronic energy of LLP is based on instructions from the HEPData entry for the CMS paper searching
23 * for neutral LLPs in the CMS endcap muon detectors: https://www.hepdata.net/record/104408
24 *
25 * \author Christina Wang
26 *
27 */
28
29#include "modules/LLPFilter.h"
30
31#include "classes/DelphesClasses.h"
32#include "classes/DelphesFactory.h"
33#include "classes/DelphesFormula.h"
34
35#include "ExRootAnalysis/ExRootClassifier.h"
36#include "ExRootAnalysis/ExRootFilter.h"
37#include "ExRootAnalysis/ExRootResult.h"
38
39#include "TDatabasePDG.h"
40#include "TFormula.h"
41#include "TLorentzVector.h"
42#include "TMath.h"
43#include "TObjArray.h"
44#include "TRandom3.h"
45#include "TString.h"
46
47#include <algorithm>
48#include <iostream>
49#include <sstream>
50#include <stdexcept>
51
52using namespace std;
53
54//------------------------------------------------------------------------------
55
56LLPFilter::LLPFilter() :
57 fItInputArray(0)
58{
59}
60
61//------------------------------------------------------------------------------
62
63LLPFilter::~LLPFilter()
64{
65}
66
67//------------------------------------------------------------------------------
68
69void LLPFilter::Init()
70{
71
72 ExRootConfParam param;
73 Size_t i, size;
74
75 // PT threshold
76 fPTMin = GetDouble("PTMin", 0.0);
77
78 fInvert = GetBool("Invert", false);
79 fDecayRegion = GetInt("DecayRegion", 0);
80 fDaughterNumber = GetInt("DaughterNumber", 0);
81
82
83
84 // no pileup
85 fRequireNotPileup = GetBool("RequireNotPileup", false);
86
87 fRequireStatus = GetBool("RequireStatus", false);
88 fStatus = GetInt("Status", 1);
89
90 fRequireCharge = GetBool("RequireCharge", false);
91 fCharge = GetInt("Charge", 1);
92
93 // import input array
94 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
95 fItInputArray = fInputArray->MakeIterator();
96
97 fParticleInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
98 fItParticleInputArray = fParticleInputArray->MakeIterator();
99
100
101 param = GetParam("PdgCode");
102 size = param.GetSize();
103
104 // read PdgCodes to be filtered out from the data card
105
106 fPdgCodes.clear();
107 for(i = 0; i < size; ++i)
108 {
109 fPdgCodes.push_back(param[i].GetInt());
110 }
111
112 // create output array
113 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
114}
115
116//------------------------------------------------------------------------------
117
118void LLPFilter::Finish()
119{
120 if(fItInputArray) delete fItInputArray;
121}
122
123//------------------------------------------------------------------------------
124
125void LLPFilter::Process()
126{
127
128 Candidate *candidate;
129 Int_t pdgCode;
130 Bool_t pass;
131 Double_t pt, eta;
132 Candidate *tempCandidate;
133
134 Candidate *daughter;
135 Int_t daughterPdg;
136
137
138 // loop over particles to find LLP
139 fItInputArray->Reset();
140 int index = -1;
141 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
142 {
143 index++;
144
145 //all distance units are in mm
146 pdgCode = candidate->PID;
147 const TLorentzVector &candidateMomentum = candidate->Momentum;
148 const TLorentzVector &candidateProdPosition = candidate->Position;
149 const TLorentzVector &candidateDecayPosition = candidate->DecayPosition;
150 // TLorentzVector candidateDecayPosition;
151 pt = candidateMomentum.Pt();
152 eta = candidateMomentum.Eta();
153 if(pt < fPTMin) continue;
154 if (fDaughterNumber > 0)
155 {
156 if (candidate->D2-candidate->D1 != fDaughterNumber) continue;//require at least fDaughterNumber daughters
157
158 }
159 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) == fPdgCodes.end()) continue; //require pdgID is one of the LLP id
160 if(fRequireStatus && (candidate->Status != fStatus)) continue;
161
162 // loop over particles to find LLP daughters and assign EM and hadronic energy
163 candidate->Eem = 0.0;
164 candidate->Ehad = 0.0;
165 fItParticleInputArray->Reset();
166
167 while((daughter = static_cast<Candidate *>(fItParticleInputArray->Next())))
168 {
169
170 daughterPdg = daughter->PID;
171 if (daughter->Status != 1)continue;
172 if (daughter->IsPU)continue;
173 if (abs(daughterPdg)==12 || abs(daughterPdg)==14 || abs(daughterPdg)==16 || abs(daughterPdg)==13)continue; // ignore neutrinos and muons
174 if (abs(daughterPdg) > 1000000) continue;//ignore BSM particles
175
176 const TLorentzVector &daughterMomentum = daughter->Momentum;
177
178 // look for mother until find LLP or reach the top of the tree
179 tempCandidate = daughter;
180 while(tempCandidate->M1 != -1 && tempCandidate->M1 != index)
181 {
182 tempCandidate = static_cast<Candidate *>(fParticleInputArray->At(tempCandidate->M1));
183 }
184 if (tempCandidate->M1 == -1) continue;
185
186 // candidateDecayPosition = daughter->Position;
187
188 if (abs(daughterPdg)==11 || abs(daughterPdg)==22 || abs(daughterPdg)==111)candidate->Eem += daughterMomentum.E();
189 else candidate->Ehad += daughterMomentum.E();
190 }
191
192 // used detector geometry in Figure 4.1.1, page141 from CERN-LHCC-97-032: https://cds.cern.ch/record/343814?ln=en
193 // decayRegion = 0: no cuts on decay region
194 // decayRegion = 1: select LLP that decays in CSC volume
195 // decayRegion = 2: select LLP that decays outside of calorimeters
196 if (fDecayRegion == 1)
197 {
198 if (abs(eta) < 2
199 && abs(candidateDecayPosition.Z())<11000 && abs(candidateDecayPosition.Z())>4000
200 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) < 6955)
201 {
202 fOutputArray->Add(candidate);
203 }
204
205 }
206 else if(fDecayRegion == 2)
207 {
208 if (abs(candidateDecayPosition.Z()) > 5680 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) > 3000)
209 {
210 fOutputArray->Add(candidate);
211 }
212 }
213 else{
214 fOutputArray->Add(candidate);
215 }
216 }//end of while loop
217}
Note: See TracBrowser for help on using the repository browser.