Fork me on GitHub

source: git/modules/LLPFilter.cc

Last change on this file was cc8716b, checked in by GitHub <noreply@…>, 3 years ago

Update to handle CMS endcap muon detector showers for long-lived particles (#103)

Co-authored-by: christinaw97 <christina.wang@…>

  • Property mode set to 100644
File size: 6.7 KB
RevLine 
[cc8716b]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 * Muons and neutrinos are ignored. Photons, electrons, and pi0 are EM energy and everything else is hadronic energy.
25 *
26 * \author Christina Wang
27 *
28 */
29
30#include "modules/LLPFilter.h"
31
32#include "classes/DelphesClasses.h"
33#include "classes/DelphesFactory.h"
34#include "classes/DelphesFormula.h"
35
36#include "ExRootAnalysis/ExRootClassifier.h"
37#include "ExRootAnalysis/ExRootFilter.h"
38#include "ExRootAnalysis/ExRootResult.h"
39
40#include "TDatabasePDG.h"
41#include "TFormula.h"
42#include "TLorentzVector.h"
43#include "TMath.h"
44#include "TObjArray.h"
45#include "TRandom3.h"
46#include "TString.h"
47
48#include <algorithm>
49#include <iostream>
50#include <sstream>
51#include <stdexcept>
52
53using namespace std;
54
55//------------------------------------------------------------------------------
56
57LLPFilter::LLPFilter() :
58 fItInputArray(0)
59{
60}
61
62//------------------------------------------------------------------------------
63
64LLPFilter::~LLPFilter()
65{
66}
67
68//------------------------------------------------------------------------------
69
70void LLPFilter::Init()
71{
72
73 ExRootConfParam param;
74 Size_t i, size;
75
76 // PT threshold
77 fPTMin = GetDouble("PTMin", 0.0);
78
79 fInvert = GetBool("Invert", false);
80 fDaughterNumber = GetInt("DaughterNumber", 0);
81 fRequireDecayRegion = GetBool("RequireDecayRegion", 0);
82
83 fDecayRegionRMax = GetDouble("DecayRegionRMax", 0.0); //mm
84 fDecayRegionRMin = GetDouble("DecayRegionRMin", 0.0); //mm
85 fDecayRegionZMax = GetDouble("DecayRegionZMax", 0.0); //mm
86 fDecayRegionZMin = GetDouble("DecayRegionZMin", 0.0); //mm
87 fDecayRegionEtaMax = GetDouble("DecayRegionEtaMax", 0.0); // requirement on abs(eta)
88 fDecayRegionEtaMin = GetDouble("DecayRegionEtaMin", 0.0); //requirement on abs(eta)
89
90
91 // no pileup
92 fRequireNotPileup = GetBool("RequireNotPileup", false);
93
94 fRequireStatus = GetBool("RequireStatus", false);
95 fStatus = GetInt("Status", 1);
96
97 fRequireCharge = GetBool("RequireCharge", false);
98 fCharge = GetInt("Charge", 1);
99
100 // import input array
101 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
102 fItInputArray = fInputArray->MakeIterator();
103
104 fParticleInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
105 fItParticleInputArray = fParticleInputArray->MakeIterator();
106
107
108 param = GetParam("PdgCode");
109 size = param.GetSize();
110
111 // read PdgCodes to be filtered out from the data card
112
113 fPdgCodes.clear();
114 for(i = 0; i < size; ++i)
115 {
116 fPdgCodes.push_back(param[i].GetInt());
117 }
118
119 // create output array
120 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
121}
122
123//------------------------------------------------------------------------------
124
125void LLPFilter::Finish()
126{
127 if(fItInputArray) delete fItInputArray;
128}
129
130//------------------------------------------------------------------------------
131
132void LLPFilter::Process()
133{
134
135 Candidate *candidate;
136 Int_t pdgCode;
137 Bool_t pass;
138 Double_t pt, eta;
139 Candidate *tempCandidate;
140
141 Candidate *daughter;
142 Int_t daughterPdg;
143
144
145 // loop over particles to find LLP
146 fItInputArray->Reset();
147 int index = -1;
148 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
149 {
150 index++;
151
152 //all distance units are in mm
153 pdgCode = candidate->PID;
154 const TLorentzVector &candidateMomentum = candidate->Momentum;
155 const TLorentzVector &candidateProdPosition = candidate->Position;
156 const TLorentzVector &candidateDecayPosition = candidate->DecayPosition;
157 pt = candidateMomentum.Pt();
158 eta = candidateMomentum.Eta();
159 if(pt < fPTMin) continue;
160 if (fDaughterNumber > 0)
161 {
162 if (candidate->D2-candidate->D1 != fDaughterNumber) continue;//require at least fDaughterNumber daughters
163
164 }
165 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) == fPdgCodes.end()) continue; //require pdgID is one of the LLP id
166 if(fRequireStatus && (candidate->Status != fStatus)) continue;
167
168 // loop over particles to find LLP daughters and assign EM and hadronic energy
169 candidate->Eem = 0.0;
170 candidate->Ehad = 0.0;
171 fItParticleInputArray->Reset();
172
173 while((daughter = static_cast<Candidate *>(fItParticleInputArray->Next())))
174 {
175
176 daughterPdg = daughter->PID;
177 if (daughter->Status != 1)continue;
178 if (daughter->IsPU)continue;
179 if (abs(daughterPdg)==12 || abs(daughterPdg)==14 || abs(daughterPdg)==16 || abs(daughterPdg)==13)continue; // ignore neutrinos and muons
180 if (abs(daughterPdg) > 1000000) continue;//ignore BSM particles
181
182 const TLorentzVector &daughterMomentum = daughter->Momentum;
183
184 // look for mother until find LLP or reach the top of the tree
185 tempCandidate = daughter;
186 while(tempCandidate->M1 != -1 && tempCandidate->M1 != index)
187 {
188 tempCandidate = static_cast<Candidate *>(fParticleInputArray->At(tempCandidate->M1));
189 }
190 if (tempCandidate->M1 == -1) continue;
191
192 // assign LLP EM or hadronic energy, depending on the daughter ID
193 if (abs(daughterPdg)==11 || abs(daughterPdg)==22 || abs(daughterPdg)==111)candidate->Eem += daughterMomentum.E();
194 else candidate->Ehad += daughterMomentum.E();
195 }
196
197 if (fRequireDecayRegion)
198 {
199 if (abs(eta) < fDecayRegionEtaMax && abs(eta) > fDecayRegionEtaMin
200 && abs(candidateDecayPosition.Z()) < fDecayRegionZMax && abs(candidateDecayPosition.Z()) > fDecayRegionZMin
201 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) < fDecayRegionRMax
202 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) > fDecayRegionRMin)
203 {
204 fOutputArray->Add(candidate);
205 }
206
207 }
208 else{
209 fOutputArray->Add(candidate);
210 }
211 }//end of while loop
212}
Note: See TracBrowser for help on using the repository browser.