ATLAS Offline Software
Loading...
Searching...
No Matches
InvariantMassDeltaPhiInclusive2Charge.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4/*********************************
5 * InvariantMassDeltaPhiInclusive2Charge.cxx
6 * Based on V Sorin 2014 implementation of InvariantMassInclusive2. For questions contact atlas-trig-l1topo-algcom@cern.ch.
7 *
8 * @brief algorithm calculates the sqr of the INVMASS between two lists and applies invmass criteria. For pairs passing the INVMASS cut a further requirement based on DeltaPhi
9 * is applied, addressing ATR-19377.
10 * Events containing a pair of TGC muons with same charge are rejected
11 *
12**********************************/
13// TO DO size of the input list to be possbly refined
14
19
20#include <cmath>
21
22REGISTER_ALG_TCS(InvariantMassDeltaPhiInclusive2Charge)
23
24
26{
27 defineParameter("InputWidth1", 9);
28 defineParameter("InputWidth2", 9);
29 defineParameter("MaxTob1", 0);
30 defineParameter("MaxTob2", 0);
31 defineParameter("NumResultBits", 6);
32 defineParameter("MinMSqr", 0, 0);
33 defineParameter("MaxMSqr", 999, 0);
34 defineParameter("MinMSqr", 0, 1);
35 defineParameter("MaxMSqr", 999, 1);
36 defineParameter("MinMSqr", 0, 2);
37 defineParameter("MaxMSqr", 999, 2);
38 defineParameter("MinMSqr", 0, 3);
39 defineParameter("MaxMSqr", 999, 3);
40 defineParameter("MinMSqr", 0, 4);
41 defineParameter("MaxMSqr", 999, 4);
42 defineParameter("MinMSqr", 0, 5);
43 defineParameter("MaxMSqr", 999, 5);
44 defineParameter("MinET1",0,0);
45 defineParameter("MinET2",0,0);
46 defineParameter("MinET1",0,1);
47 defineParameter("MinET2",0,1);
48 defineParameter("MinET1",0,2);
49 defineParameter("MinET2",0,2);
50 defineParameter("MinET1",0,3);
51 defineParameter("MinET2",0,3);
52 defineParameter("MinET1",0,4);
53 defineParameter("MinET2",0,4);
54 defineParameter("MinET1",0,5);
55 defineParameter("MinET2",0,5);
56 defineParameter("ApplyEtaCut", 0);
57 defineParameter("MinEta1", 0);
58 defineParameter("MaxEta1", 31);
59 defineParameter("MinEta2", 0);
60 defineParameter("MaxEta2", 31);
61 defineParameter("MinDeltaPhi", 0, 0);
62 defineParameter("MaxDeltaPhi", 31, 0);
63 defineParameter("MinDeltaPhi", 0, 1);
64 defineParameter("MaxDeltaPhi", 31, 1);
65 defineParameter("MinDeltaPhi", 0, 2);
66 defineParameter("MaxDeltaPhi", 31, 2);
67 defineParameter("MinDeltaPhi", 0, 3);
68 defineParameter("MaxDeltaPhi", 31, 3);
69 defineParameter("MinDeltaPhi", 0, 4);
70 defineParameter("MaxDeltaPhi", 31, 4);
71 defineParameter("MinDeltaPhi", 0, 5);
72 defineParameter("MaxDeltaPhi", 31, 5);
73 //does this need to change?
75}
76
78
79
82 p_NumberLeading1 = parameter("InputWidth1").value();
83 p_NumberLeading2 = parameter("InputWidth2").value();
84 if(parameter("MaxTob1").value() > 0) p_NumberLeading1 = parameter("MaxTob1").value();
85 if(parameter("MaxTob2").value() > 0) p_NumberLeading2 = parameter("MaxTob2").value();
86
87 for(unsigned int i=0; i<numberOutputBits(); ++i) {
88 p_InvMassMin[i] = parameter("MinMSqr", i).value();
89 p_InvMassMax[i] = parameter("MaxMSqr", i).value();
90 p_DeltaPhiMin[i] = parameter("MinDeltaPhi", i).value();
91 p_DeltaPhiMax[i] = parameter("MaxDeltaPhi", i).value();
92 p_MinET1[i] = parameter("MinET1",i).value();
93 p_MinET2[i] = parameter("MinET2",i).value();
94 }
95
96 TRG_MSG_INFO("NumberLeading1 : " << p_NumberLeading1);
97 TRG_MSG_INFO("NumberLeading2 : " << p_NumberLeading2);
98 for(unsigned int i=0; i<numberOutputBits(); ++i) {
99 TRG_MSG_INFO("InvMassMin : " << p_InvMassMin[i]);
100 TRG_MSG_INFO("InvMassMax : " << p_InvMassMax[i]);
101 TRG_MSG_INFO("DeltaPhiMin : " << p_DeltaPhiMin[i]);
102 TRG_MSG_INFO("DeltaPhiMax : " << p_DeltaPhiMax[i]);
103 TRG_MSG_INFO("MinET1 : " << p_MinET1[i]);
104 TRG_MSG_INFO("MinET2 : " << p_MinET2[i]);
105 }
106
107 p_ApplyEtaCut = parameter("ApplyEtaCut").value();
108 p_MinEta1 = parameter("MinEta1" ).value();
109 p_MaxEta1 = parameter("MaxEta1" ).value();
110 p_MinEta2 = parameter("MinEta2" ).value();
111 p_MaxEta2 = parameter("MaxEta2" ).value();
112 TRG_MSG_INFO("ApplyEtaCut : "<<p_ApplyEtaCut );
113 TRG_MSG_INFO("MinEta1 : "<<p_MinEta1 );
114 TRG_MSG_INFO("MaxEta1 : "<<p_MaxEta1 );
115 TRG_MSG_INFO("MinEta2 : "<<p_MinEta2 );
116 TRG_MSG_INFO("MaxEta2 : "<<p_MaxEta2 );
117
118 TRG_MSG_INFO("number output : " << numberOutputBits());
119
120 // book histograms
121 for(unsigned int i=0; i<numberOutputBits(); ++i) {
122 std::string hname_accept = "hInvariantMassDeltaPhiInclusive2Charge_accept_bit"+std::to_string(static_cast<int>(i));
123 std::string hname_reject = "hInvariantMassDeltaPhiInclusive2Charge_reject_bit"+std::to_string(static_cast<int>(i));
124 // mass
125 bookHist(m_histAcceptM, hname_accept, "INVM vs DPHI", 100, std::sqrt(p_InvMassMin[i]), std::sqrt(p_InvMassMax[i]), 100, p_DeltaPhiMin[i], p_DeltaPhiMax[i]);
126 bookHist(m_histRejectM, hname_reject, "INVM vs DPHI", 100, std::sqrt(p_InvMassMin[i]), std::sqrt(p_InvMassMax[i]), 100, p_DeltaPhiMin[i], p_DeltaPhiMax[i]);
127 // eta2 vs. eta1
128 bookHist(m_histAcceptEta1Eta2, hname_accept, "ETA vs ETA", 100, p_MinEta1, p_MaxEta1, 100, p_MinEta2, p_MaxEta2);
129 bookHist(m_histRejectEta1Eta2, hname_reject, "ETA vs ETA", 100, p_MinEta1, p_MaxEta1, 100, p_MinEta2, p_MaxEta2);
130 }
131 return StatusCode::SUCCESS;
132}
133
134
135
137TCS::InvariantMassDeltaPhiInclusive2Charge::processBitCorrect( const std::vector<TCS::TOBArray const *> & input,
138 const std::vector<TCS::TOBArray *> & output,
139 Decision & decision )
140{
141
142 if( input.size() == 2) {
143 for( TOBArray::const_iterator tob1 = input[0]->begin();
144 tob1 != input[0]->end() && distance(input[0]->begin(), tob1) < p_NumberLeading1;
145 ++tob1)
146 {
147
148
149 for( TCS::TOBArray::const_iterator tob2 = input[1]->begin();
150 tob2 != input[1]->end() && distance(input[1]->begin(), tob2) < p_NumberLeading2;
151 ++tob2) {
152 // Inv Mass calculation
153 unsigned int invmass2 = calcInvMassBW( *tob1, *tob2 );
154 // test DeltaPhiMin, DeltaPhiMax
155 unsigned int deltaPhi = calcDeltaPhiBW( *tob1, *tob2 );
156 const int eta1 = (*tob1)->eta();
157 const int eta2 = (*tob2)->eta();
158 const unsigned int aeta1 = std::abs(eta1);
159 const unsigned int aeta2 = std::abs(eta2);
160 // Charge cut ( 1 = positive, -1 = negative, 0 = undefined (RPC) )
161 int charge1 = (*tob1)->charge();
162 int charge2 = (*tob2)->charge();
163 int totalCharge = charge1 + charge2;
164 bool acceptCharge = true;
165 if ( std::abs(totalCharge) == 2 ) { acceptCharge = false; }
166 for(unsigned int i=0; i<numberOutputBits(); ++i) {
167 bool accept = false;
168 if( static_cast<parType_t>((*tob1)->Et()) <= p_MinET1[i]) continue; // ET cut
169 if( static_cast<parType_t>((*tob2)->Et()) <= p_MinET2[i]) continue; // ET cut
170 if(p_ApplyEtaCut &&
171 ((aeta1 < p_MinEta1 || aeta1 > p_MaxEta1 ) ||
172 (aeta2 < p_MinEta2 || aeta2 > p_MaxEta2 ) )) continue;
173 accept = invmass2 >= p_InvMassMin[i] && invmass2 <= p_InvMassMax[i] && deltaPhi >= p_DeltaPhiMin[i] && deltaPhi <= p_DeltaPhiMax[i] && acceptCharge;
174 const bool fillAccept = fillHistos() and (fillHistosBasedOnHardware() ? getDecisionHardwareBit(i) : accept);
175 const bool fillReject = fillHistos() and not fillAccept;
176 const bool alreadyFilled = decision.bit(i);
177 if( accept ) {
178 decision.setBit(i, true);
179 output[i]->push_back( TCS::CompositeTOB(*tob1, *tob2) );
180 }
181 if(fillAccept and not alreadyFilled) {
182 fillHist2D(m_histAcceptM[i],std::sqrt(static_cast<float>(invmass2)),static_cast<float>(deltaPhi));
183 fillHist2D(m_histAcceptEta1Eta2[i],eta1, eta2);
184 } else if(fillReject) {
185 fillHist2D(m_histRejectM[i],std::sqrt(static_cast<float>(invmass2)),static_cast<float>(deltaPhi));
186 fillHist2D(m_histRejectEta1Eta2[i],eta1, eta2);
187 }
188 TRG_MSG_DEBUG("Decision " << i << ": " << (accept?"pass":"fail") << " invmass2 = " << invmass2);
189 }
190 }
191 }
192 for (unsigned int i=0; i < numberOutputBits(); ++i) {
193 bool hasAmbiguousInputs = TSU::isAmbiguousTruncation(input[0], p_NumberLeading1, p_MinET1[i])
195 output[i]->setAmbiguityFlag(hasAmbiguousInputs);
196 }
197 } else {
198
199 TCS_EXCEPTION("InvariantMassDeltaPhiInclusive2Charge alg must have 2 inputs, but got " << input.size());
200
201 }
203
204}
205
207TCS::InvariantMassDeltaPhiInclusive2Charge::process( const std::vector<TCS::TOBArray const *> & input,
208 const std::vector<TCS::TOBArray *> & output,
209 Decision & decision )
210{
211
212
213 if( input.size() == 2) {
214 for( TOBArray::const_iterator tob1 = input[0]->begin();
215 tob1 != input[0]->end() && distance(input[0]->begin(), tob1) < p_NumberLeading1;
216 ++tob1)
217 {
218 for( TCS::TOBArray::const_iterator tob2 = input[1]->begin();
219 tob2 != input[1]->end() && distance(input[1]->begin(), tob2) < p_NumberLeading2;
220 ++tob2) {
221 // Inv Mass calculation
222 unsigned int invmass2 = calcInvMass( *tob1, *tob2 );
223 // test DeltaPhiMin, DeltaPhiMax
224 unsigned int deltaPhi = calcDeltaPhi( *tob1, *tob2 );
225 const int eta1 = (*tob1)->eta();
226 const int eta2 = (*tob2)->eta();
227 const unsigned int aeta1 = std::abs(eta1);
228 const unsigned int aeta2 = std::abs(eta2);
229 // Charge cut ( 1 = positive, -1 = negative, 0 = undefined (RPC) )
230 int charge1 = (*tob1)->charge();
231 int charge2 = (*tob2)->charge();
232 int totalCharge = charge1 + charge2;
233 bool acceptCharge = true;
234 if ( std::abs(totalCharge) == 2 ) { acceptCharge = true; }
235 for(unsigned int i=0; i<numberOutputBits(); ++i) {
236 if( static_cast<parType_t>((*tob1)->Et()) <= p_MinET1[i]) continue; // ET cut
237 if( static_cast<parType_t>((*tob2)->Et()) <= p_MinET2[i]) continue; // ET cut
238 if(p_ApplyEtaCut &&
239 ((aeta1 < p_MinEta1 || aeta1 > p_MaxEta1 ) ||
240 (aeta2 < p_MinEta2 || aeta2 > p_MaxEta2 ) )) continue;
241 bool accept = invmass2 >= p_InvMassMin[i] && invmass2 <= p_InvMassMax[i] && deltaPhi >= p_DeltaPhiMin[i] && deltaPhi <= p_DeltaPhiMax[i] && acceptCharge;
242 const bool fillAccept = fillHistos() and (fillHistosBasedOnHardware() ? getDecisionHardwareBit(i) : accept);
243 const bool fillReject = fillHistos() and not fillAccept;
244 const bool alreadyFilled = decision.bit(i);
245 if( accept ) {
246 decision.setBit(i, true);
247 output[i]->push_back( TCS::CompositeTOB(*tob1, *tob2) );
248 }
249 if(fillAccept and not alreadyFilled) {
250 fillHist2D(m_histAcceptM[i],std::sqrt(static_cast<float>(invmass2)),static_cast<float>(deltaPhi));
251 fillHist2D(m_histAcceptEta1Eta2[i],eta1, eta2);
252 } else if(fillReject) {
253 fillHist2D(m_histRejectM[i],std::sqrt(static_cast<float>(invmass2)),static_cast<float>(deltaPhi));
254 fillHist2D(m_histRejectEta1Eta2[i],eta1, eta2);
255 }
256 TRG_MSG_DEBUG("Decision " << i << ": " << (accept ?"pass":"fail") << " invmass2 = " << invmass2);
257 }
258 }
259 }
260 } else {
261 TCS_EXCEPTION("InvariantMassDeltaPhiInclusive2Charge alg must have 2 inputs, but got " << input.size());
262 }
264}
#define REGISTER_ALG_TCS(CLASS)
Definition AlgFactory.h:62
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
const Parameter & parameter(const std::string &parameterName) const
const std::string & name() const
void bookHist(std::vector< std::string > &regName, const std::string &name, const std::string &title, const int binx, const int xmin, const int xmax)
unsigned int calcDeltaPhiBW(const TCS::GenericTOB *tob1, const TCS::GenericTOB *tob2)
unsigned int calcInvMassBW(const TCS::GenericTOB *tob1, const TCS::GenericTOB *tob2)
void defineParameter(const std::string &name, TCS::parType_t value)
unsigned int calcDeltaPhi(const TCS::GenericTOB *tob1, const TCS::GenericTOB *tob2)
unsigned int calcInvMass(const TCS::GenericTOB *tob1, const TCS::GenericTOB *tob2)
void fillHist2D(const std::string &histName, double x, double y)
data_t::const_iterator const_iterator
std::vector< std::string > m_histAcceptEta1Eta2
Definition DecisionAlg.h:79
std::vector< std::string > m_histAcceptM
Definition DecisionAlg.h:75
void setNumberOutputBits(unsigned int numberOutputBits)
Definition DecisionAlg.h:40
DecisionAlg(const std::string &name)
Definition DecisionAlg.h:25
std::vector< std::string > m_histRejectEta1Eta2
Definition DecisionAlg.h:80
bool fillHistosBasedOnHardware() const
! getter
bool fillHistos() const
whether the monitoring histograms should be filled
unsigned int numberOutputBits() const
Definition DecisionAlg.h:39
bool getDecisionHardwareBit(const unsigned int &bitNumber) const
! get one hardware decision bit from this algo
std::vector< std::string > m_histRejectM
Definition DecisionAlg.h:76
bool bit(unsigned int index) const
Definition Decision.h:40
void setBit(unsigned int index, bool value)
Definition Decision.cxx:12
virtual StatusCode processBitCorrect(const std::vector< TCS::TOBArray const * > &input, const std::vector< TCS::TOBArray * > &output, Decision &decison) override final
virtual StatusCode process(const std::vector< TCS::TOBArray const * > &input, const std::vector< TCS::TOBArray * > &output, Decision &decison) override final
uint32_t parType_t
Definition Parameter.h:22
bool isAmbiguousTruncation(TCS::TOBArray const *tobs, size_t pos, unsigned minEt=0)
STL namespace.