Fork me on GitHub

source: git/modules/LLPFilter.cc@ d612dec

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

added CscCluster modules and classes

  • Property mode set to 100644
File size: 6.6 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 * Removes particles with specific PDG codes
22 *
23 * \author M. Selvaggi
24 *
25 */
26
27#include "modules/LLPFilter.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
36
37#include "TDatabasePDG.h"
38#include "TFormula.h"
39#include "TLorentzVector.h"
40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
44
45#include <algorithm>
46#include <iostream>
47#include <sstream>
48#include <stdexcept>
49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54LLPFilter::LLPFilter() :
55 fItInputArray(0)
56{
57}
58
59//------------------------------------------------------------------------------
60
61LLPFilter::~LLPFilter()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void LLPFilter::Init()
68{
69
70 ExRootConfParam param;
71 Size_t i, size;
72
73 // PT threshold
74 fPTMin = GetDouble("PTMin", 0.0);
75
76 fInvert = GetBool("Invert", false);
77 fDecayRegion = GetInt("DecayRegion", 0);
78 fDaughterNumber = GetInt("DaughterNumber", 0);
79
80
81
82 // no pileup
83 fRequireNotPileup = GetBool("RequireNotPileup", false);
84
85 fRequireStatus = GetBool("RequireStatus", false);
86 fStatus = GetInt("Status", 1);
87
88 fRequireCharge = GetBool("RequireCharge", false);
89 fCharge = GetInt("Charge", 1);
90
91 // import input array
92 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
93 fItInputArray = fInputArray->MakeIterator();
94
95 fParticleInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
96 fItParticleInputArray = fParticleInputArray->MakeIterator();
97
98
99 param = GetParam("PdgCode");
100 size = param.GetSize();
101
102 // read PdgCodes to be filtered out from the data card
103
104 fPdgCodes.clear();
105 for(i = 0; i < size; ++i)
106 {
107 fPdgCodes.push_back(param[i].GetInt());
108 }
109
110 // create output array
111 fOutputArray = ExportArray(GetString("OutputArray", "filteredParticles"));
112}
113
114//------------------------------------------------------------------------------
115
116void LLPFilter::Finish()
117{
118 if(fItInputArray) delete fItInputArray;
119}
120
121//------------------------------------------------------------------------------
122
123void LLPFilter::Process()
124{
125
126 Candidate *candidate;
127 Int_t pdgCode;
128 Bool_t pass;
129 Double_t pt, eta;
130 Candidate *tempCandidate;
131
132 Candidate *daughter;
133 Int_t daughterPdg;
134
135
136 // loop over particles to find LLP
137 fItInputArray->Reset();
138 int index = -1;
139 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
140 {
141 index++;
142
143 //all distance units are in mm
144 pdgCode = candidate->PID;
145 const TLorentzVector &candidateMomentum = candidate->Momentum;
146 const TLorentzVector &candidateProdPosition = candidate->Position;
147 const TLorentzVector &candidateDecayPosition = candidate->DecayPosition;
148 pt = candidateMomentum.Pt();
149 eta = candidateMomentum.Eta();
150 if(pt < fPTMin) continue;
151 if (fDaughterNumber > 0)
152 {
153 if (candidate->D2-candidate->D1 != fDaughterNumber) continue;//require 3 daughters
154
155 }
156 if(find(fPdgCodes.begin(), fPdgCodes.end(), pdgCode) == fPdgCodes.end()) continue; //require pdgID is one of the LLP id
157 if(fRequireStatus && (candidate->Status != fStatus)) continue;
158
159 // loop over particles to find LLP daughters and assign EM and hadronic energy
160 candidate->Eem = 0.0;
161 candidate->Ehad = 0.0;
162
163 fItParticleInputArray->Reset();
164
165 while((daughter = static_cast<Candidate *>(fItParticleInputArray->Next())))
166 {
167
168 // check mother is LLP
169 //check ID is 123456 or neutrinos or charged leptons
170 daughterPdg = daughter->PID;
171 if (daughter->Status != 1)continue;
172
173 if (daughter->IsPU)continue;
174
175 if (abs(daughterPdg)==12 || abs(daughterPdg)==14 || abs(daughterPdg)==16 || abs(daughterPdg)==13)continue; // ignore neutrinos and muons
176
177
178 if (abs(daughterPdg) > 1000000) continue;//ignore BSM particles
179 const TLorentzVector &daughterProdPosition = daughter->Position;
180 const TLorentzVector &daughterMomentum = daughter->Momentum;
181
182 const TLorentzVector distance = daughterProdPosition - candidateDecayPosition;
183
184 if (sqrt(pow(distance.X(), 2) + pow(distance.Y(), 2) + pow(distance.Z(), 2))>100) continue; //if vertices are close, then matched, this 1 m only works for those decay in muon system, would be too large for tracker related
185 tempCandidate = daughter;
186 while(tempCandidate->M1 != -1 && tempCandidate->M1 != index)
187 {
188 tempCandidate = static_cast<Candidate *>(fParticleInputArray->At(tempCandidate->M1));
189 // if (tempCandidate->PID == pdgCode)break;
190
191 }
192 if (tempCandidate->M1 == -1) continue;
193
194 if (abs(daughterPdg)==11 || abs(daughterPdg)==22 || abs(daughterPdg)==111)candidate->Eem += daughterMomentum.E();
195 else candidate->Ehad += daughterMomentum.E();
196
197
198 }
199
200
201
202 // used detector geometry in Figure 4.1.1, page141 from CERN-LHCC-97-032: https://cds.cern.ch/record/343814?ln=en
203 // decayRegion = 0: no cuts on decay region
204 // decayRegion = 1: select LLP that decays in CSC volume
205 // decayRegion = 2: select LLP that decays outside of calorimeters
206 if (fDecayRegion == 1)
207 {
208 if (abs(eta) < 2
209 && abs(candidateDecayPosition.Z())<11000 && abs(candidateDecayPosition.Z())>4000
210 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) < 6955)
211 {
212 fOutputArray->Add(candidate);
213 }
214
215 }
216 else if(fDecayRegion == 2)
217 {
218 if (abs(candidateDecayPosition.Z()) > 5680 && sqrt(pow(candidateDecayPosition.X(),2)+pow(candidateDecayPosition.Y(),2)) > 3000)
219 {
220 fOutputArray->Add(candidate);
221 }
222 }
223 else{
224 fOutputArray->Add(candidate);
225 }
226 }//end of while loop
227}
Note: See TracBrowser for help on using the repository browser.