ATLAS Offline Software
Loading...
Searching...
No Matches
MissingMassCalculator.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5/*
6------ MissingMassCalculator
7------ Author: Aliaksandr Pranko (appranko@lbl.gov)
8------ Code developers: David Rousseau (rousseau@lal.in2p3.fr), Dimitris Varouchas (Dimitris.Varouchas@cern.ch)
9MissingMassCalculator is designed to reconstruct mass in
10events where two particles decay into states with missing ET.
11*/
12#ifndef MissingMassCalculator_h
13#define MissingMassCalculator_h
14
15
16
17#if !defined (__CINT__) || defined (__MAKECINT__)
18
19#include <TRandom2.h>
20#include <TH1.h>
21#include <TGraph.h>
22#include <TF1.h>
23#include <TMath.h>
24#include <Math/Vector4D.h>
25#include <Math/Vector2D.h>
26#include <vector>
27#include <TObject.h>
28#include <TDirectory.h>
29
30#include <memory>
31#include <string>
32
37
39
40#endif
41
42
43namespace DiTauMassTools{
44 using ROOT::Math::PtEtaPhiMVector;
45
47
48 public:
49
50 private:
51
52 //---------------- structures
53 struct DitauStuff {
54 double Mditau_best{}; // best fitted M(ditau)
55 double Sign_best{}; // best significance of M(ditau) fit
56 PtEtaPhiMVector nutau1; // fitted 4-vec for neutrino from tau-1
57 PtEtaPhiMVector nutau2; // fitted 4-vec for neutrino from tau-2
58 PtEtaPhiMVector vistau1; // fitted 4-vec for visible tau-1
59 PtEtaPhiMVector vistau2; // fitted 4-vec for visible tau-2
60 double RMSoverMPV{};
61 };
62
63 TRandom2 m_randomGen;
64
66
67 bool m_fUseEfficiencyRecovery{}; // switch to turn ON/OFF re-fit in order to recover efficiency
68 bool m_fUseFloatStopping{}; // switch to turn ON/OFF floating stopping criterion
72
76 int m_rmsStop{};
77 double m_meanbinStop{};
78
79 // these are temporary vectors. Declared globally to avoid construction/destruction
80 std::vector<PtEtaPhiMVector> m_nuvecsol1;
81 std::vector<PtEtaPhiMVector> m_nuvecsol2;
82
83 std::vector<PtEtaPhiMVector> m_tauvecsol1;
84 std::vector<PtEtaPhiMVector> m_tauvecsol2;
85 std::vector<double> m_tauvecprob1;
86 std::vector<double> m_tauvecprob2;
87
88 std::vector<PtEtaPhiMVector> m_nuvec1_tmp;
89 std::vector<PtEtaPhiMVector> m_nuvec2_tmp;
90
91 PtEtaPhiMVector m_tautau_tmp;
92
94
97 double m_nsigma_METscan_lfv_lh{}, m_beamEnergy{}; // number of sigmas for MET-scan
98
100
101 double m_prob_tmp{};
102
104 double m_mtautauSum{};
105
107 int m_seed{};
108 // data member for the spaceWalker approach
109
110 int m_iter0{};
116 int m_nosol1{};
117 int m_nosol2{};
119 bool m_switch1{};
120 bool m_switch2{};
121
123
129
133
137
138 double m_mTau{},m_mTau2{};
140 double m_eTau1{}, m_eTau2{};
141 double m_eTau10{}, m_eTau20{};
150
156 std::vector<double> m_probFinalSolOldVec;
157 std::vector<double> m_mtautauFinalSolOldVec;
158 std::vector<PtEtaPhiMVector> m_nu1FinalSolOldVec;
159 std::vector<PtEtaPhiMVector> m_nu2FinalSolOldVec;
160
161 int m_nsol{};
162 std::vector<double> m_probFinalSolVec;
163 std::vector<double> m_mtautauFinalSolVec;
164 std::vector<PtEtaPhiMVector> m_nu1FinalSolVec;
165 std::vector<PtEtaPhiMVector> m_nu2FinalSolVec;
166
167
170 double m_walkWeight{};
172
174
175 PtEtaPhiMVector m_tauVec1,m_tauVec2;
181 double m_tauVec1E{};
182 double m_tauVec2E{};
183 double m_m2Nu1{};
184 double m_m2Nu2{};
185 double m_ET2v1{};
186 double m_ET2v2{};
187 double m_E2v1{};
188 double m_E2v2{};
189 double m_Ev2{};
190 double m_Ev1{};
191 double m_Mvis{},m_Meff{};
192
193 //--- define histograms for histogram method
194 //--- upper limits need to be revisied in the future!!! It may be not enough for some analyses
195 std::shared_ptr<TH1F> m_fMfit_all;
196 std::shared_ptr<TH1F> m_fMEtP_all;
197 std::shared_ptr<TH1F> m_fMEtL_all;
198 std::shared_ptr<TH1F> m_fMnu1_all;
199 std::shared_ptr<TH1F> m_fMnu2_all;
200 std::shared_ptr<TH1F> m_fPhi1_all;
201 std::shared_ptr<TH1F> m_fPhi2_all;
202 std::shared_ptr<TGraph> m_fMfit_allGraph;
203 std::shared_ptr<TH1F> m_fMfit_allNoWeight;
204
205 std::shared_ptr<TH1F> m_fPXfit1;
206 std::shared_ptr<TH1F> m_fPYfit1;
207 std::shared_ptr<TH1F> m_fPZfit1;
208 std::shared_ptr<TH1F> m_fPXfit2;
209 std::shared_ptr<TH1F> m_fPYfit2;
210 std::shared_ptr<TH1F> m_fPZfit2;
211
212 // these histograms are used for the floating stopping criterion
213 std::shared_ptr<TH1F> m_fMmass_split1;
214 std::shared_ptr<TH1F> m_fMEtP_split1;
215 std::shared_ptr<TH1F> m_fMEtL_split1;
216 std::shared_ptr<TH1F> m_fMnu1_split1;
217 std::shared_ptr<TH1F> m_fMnu2_split1;
218 std::shared_ptr<TH1F> m_fPhi1_split1;
219 std::shared_ptr<TH1F> m_fPhi2_split1;
220 std::shared_ptr<TH1F> m_fMmass_split2;
221 std::shared_ptr<TH1F> m_fMEtP_split2;
222 std::shared_ptr<TH1F> m_fMEtL_split2;
223 std::shared_ptr<TH1F> m_fMnu1_split2;
224 std::shared_ptr<TH1F> m_fMnu2_split2;
225 std::shared_ptr<TH1F> m_fPhi1_split2;
226 std::shared_ptr<TH1F> m_fPhi2_split2;
227
229
230 TH1F* m_fPhi1{};
231 TH1F* m_fPhi2{};
232 TH1F* m_fMnu1{};
233 TH1F* m_fMnu2{};
234 TH1F* m_fMetx{};
235 TH1F* m_fMety{};
236 TH1F* m_fTheta3D{};
237 TH1F* m_fTauProb{};
238
239 // for intermediate calc
240 PtEtaPhiMVector m_TLVdummy;
241
242 //---------------- protected variables
243 DitauStuff m_fDitauStuffFit; // results based on fit method
244 DitauStuff m_fDitauStuffHisto; // results based on histo method
245
246 int m_niter_fit1{}; // number of iterations for dR-dPhi scan
247 int m_niter_fit2{}; // number of iterations for MET-scan
248 int m_niter_fit3{}; // number of iterations for Mnu-scan
249 int m_NiterRandom{}; // number of random iterations (for lh, multiply or divide by 10 for ll and hh)
252 int m_RndmSeedAltering{}; // reset seed (not necessary by default)
253
254 double m_dRmax_tau{}; // maximum dR(nu-visTau)
255
256 double m_MnuScanRange{}; // range of M(nunu) scan; M(nunu) range can be affected by selection cuts
257
258
259 //---------------- protected functions
260 void ClearDitauStuff(DitauStuff &fStuff);
261 void DoOutputInfo();
262 void PrintOtherInput();
263 void PrintResults();
264
265 inline int NuPsolutionV3(const double & mNu1, const double & mNu2, const double & phi1, const double & phi2,
266 int & nsol1, int & nsol2);
267
268 inline int NuPsolutionLFV(const XYVector & met_vec, const PtEtaPhiMVector & tau,
269 const double & m_nu, std::vector<PtEtaPhiMVector> &nu_vec);
270
271
272 protected:
273 inline int CheckSolutions(PtEtaPhiMVector nu_vec, PtEtaPhiMVector vis_vec, int decayType);
274 inline int TailCleanUp(const PtEtaPhiMVector & vis1, const PtEtaPhiMVector & nu1,
275 const PtEtaPhiMVector & vis2, const PtEtaPhiMVector & nu2,
276 const double & mmc_mass, const double & vis_mass, const double & eff_mass, const double & dphiTT);
277
278
279 inline int refineSolutions ( const double & M_nu1, const double & M_nu2,
280 const int nsol1, const int nsol2,
281 const double & Mvis, const double & Meff);
282
283
284
285 inline void handleSolutions();
286
287 inline double MassScale(int method, double mass, const int & tau_type1, const int & tau_type2);
288
289
290 // factor out parameter space walking and probablity computing
291 inline int DitauMassCalculatorV9walk();
292
293
294 // Calculates mass of lep+tau system in LFV X->lep+tau decays
295 // It is based on DitauMassCalculatorV9, not optimized for speed yet, simple phase-space scan
296 inline int DitauMassCalculatorV9lfv(bool refit);
297
298
299
300 // only compute probability
301 inline int probCalculatorV9fast(
302 const double & phi1, const double & phi2,
303 const double & M_nu1, const double & M_nu2);
304
305 // initialize the walker
306 inline void SpaceWalkerInit();
307
308 //walk the walker
309 inline bool SpaceWalkerWalk();
310
311 inline bool precomputeCache();
312
313
314 inline bool checkMEtInRange () ;
315 inline bool checkAllParamInRange () ;
316
317 //----------------------------------------------
318 //
319 // ------------ Public methods ---------------
320 //
321 //______________________________________________
322
323public:
324
326
327 MissingMassCalculator(MMCCalibrationSet::e aset, std::string paramFilePath) ;
328
331
335
336 int RunMissingMassCalculator( const xAOD::IParticle* part1, const xAOD::IParticle* part2, const xAOD::MissingET* met, const int& njets );
337
338 //-------- Set Input Parameters
339 void FinalizeSettings(const xAOD::IParticle* part1, const xAOD::IParticle* part2, const xAOD::MissingET* met, const int& njets );
340 void SetNiterFit1(const int val) { m_niter_fit1=val; } // number of iterations per loop in dPhi loop
341 void SetNiterFit2(const int val) { m_niter_fit2=val; } // number of iterations per loop in MET loop
342 void SetNiterFit3(const int val) { m_niter_fit3=val; } // number of iterations per loop in Mnu loop
343 void SetNiterRandom(const int val) { m_NiterRandom=val; } // number of random iterations
344 void SetNsucStop(const int val) { m_NsucStop=val; } // Arrest criteria for Nsuccesses
345 void SetRMSStop(const int val) { m_RMSStop=val;}
346 void SetMeanbinStop(const double val) {m_meanbinStop=val;}
347 void SetRndmSeedAltering(const int val) { m_RndmSeedAltering=val; } // number of iterations per loop in Mnu loop
348 void SetEventNumber(const int eventNumber) { m_eventNumber = eventNumber; }
349
350 void SetMnuScanRange(const double val) { m_MnuScanRange=val; }
351
352 void SetProposalTryMEt(const double val) {m_proposalTryMEt=val; }
353 void SetProposalTryPhi(const double val) {m_ProposalTryPhi=val;}
354 void SetProposalTryMnu(const double val) {m_ProposalTryMnu=val;}
355
358
359 int GetNiterFit1() const { return m_niter_fit1; } // number of iterations per loop in dPhi loop
360 int GetNiterFit2() const { return m_niter_fit2; } // number of iterations per loop in MET loop
361 int GetNiterFit3() const { return m_niter_fit3; } // number of iterations per loop in Mnu loop
362 int GetNiterRandom() const { return m_niterRandomLocal; } // number of random iterations
363
364 int GetNsucStop() const { return m_NsucStop; } // Arrest criteria for NSuc
365 int GetRMSStop() const { return m_RMSStop; }
366 double GetMeanbinStop() const { return m_meanbinStop;}
367 int GetRndmSeedAltering() const { return m_RndmSeedAltering; } // number of iterations per loop in Mnu loop
368
372 int GetMarkovNAccept() const { return m_markovNAccept; }
374 double GetProposalTryMEt() const {return m_proposalTryMEt;}
375 double GetProposalTryPhi() const {return m_ProposalTryPhi;}
376 double GetProposalTryMnu() const {return m_ProposalTryMnu;}
377
378 void SetNsigmaMETscan_ll(const double val) { m_nsigma_METscan_ll=val; } // number of sigma's for MET-scan in ll events
379 void SetNsigmaMETscan_lh(const double val) { m_nsigma_METscan_lh=val; } // number of sigma's for MET-scan in lh events
380 void SetNsigmaMETscan_hh(const double val) { m_nsigma_METscan_hh=val; } // number of sigma's for MET-scan in hh events
381 void SetNsigmaMETscan(const double val) { m_nsigma_METscan=val; } // number of sigma's for MET-scan
382
383 void SetUseFloatStopping(const bool val); // switch for floating stopping criterion
386 void SetFloatStoppingComp(const double val) { m_fUseFloatStoppingComp = val;}
387 void SetBeamEnergy(const double val) { m_beamEnergy=val; }
388 void SetLFVLeplepRefit(const bool val) { m_lfvLeplepRefit=val; }
389 void SaveLlhHisto(const bool val);
390
391 double GetmMaxError() const {return m_PrintmMaxError;}
392 double GetmMeanError() const { return m_PrintmMeanError;}
394
395 int GetNNoSol() const {return m_markovNRejectNoSol;}
397 int GetNSol() const {return m_markovNAccept;}
398
399
400 //-------- Get results;
401 Double_t maxFitting(Double_t *x, Double_t *par);
402
403 // compute maximum from histo
404 double maxFromHist(TH1F *theHist, std::vector<double> & histInfo, const MaxHistStrategy::e maxHistStrategy=MaxHistStrategy::FIT,const int winHalfWidth=2,bool debug=false);
405 double maxFromHist(const std::shared_ptr<TH1F>& theHist, std::vector<double> & histInfo, const MaxHistStrategy::e maxHistStrategy=MaxHistStrategy::FIT,const int winHalfWidth=2,bool debug=false) {
406 return maxFromHist(theHist.get(), histInfo, maxHistStrategy, winHalfWidth, debug);
407 }
408
409 XYVector metvec_tmp;
410 inline double dTheta3DLimit(const int & tau_type, const int & limit_code,const double & P_tau);
411
412};
413} // namespace DiTauMassTools
414
415#endif
const bool debug
#define x
int TailCleanUp(const PtEtaPhiMVector &vis1, const PtEtaPhiMVector &nu1, const PtEtaPhiMVector &vis2, const PtEtaPhiMVector &nu2, const double &mmc_mass, const double &vis_mass, const double &eff_mass, const double &dphiTT)
double maxFromHist(TH1F *theHist, std::vector< double > &histInfo, const MaxHistStrategy::e maxHistStrategy=MaxHistStrategy::FIT, const int winHalfWidth=2, bool debug=false)
int probCalculatorV9fast(const double &phi1, const double &phi2, const double &M_nu1, const double &M_nu2)
void FinalizeSettings(const xAOD::IParticle *part1, const xAOD::IParticle *part2, const xAOD::MissingET *met, const int &njets)
std::vector< PtEtaPhiMVector > m_nu2FinalSolOldVec
int refineSolutions(const double &M_nu1, const double &M_nu2, const int nsol1, const int nsol2, const double &Mvis, const double &Meff)
MissingMassCalculator(MMCCalibrationSet::e aset, std::string paramFilePath)
MissingMassCalculator & operator=(const MissingMassCalculator &)=delete
int CheckSolutions(PtEtaPhiMVector nu_vec, PtEtaPhiMVector vis_vec, int decayType)
std::vector< PtEtaPhiMVector > m_nu1FinalSolOldVec
std::vector< PtEtaPhiMVector > m_tauvecsol1
std::vector< PtEtaPhiMVector > m_nuvec1_tmp
std::vector< PtEtaPhiMVector > m_nuvecsol1
std::vector< PtEtaPhiMVector > m_nuvec2_tmp
std::vector< PtEtaPhiMVector > m_nuvecsol2
std::vector< PtEtaPhiMVector > m_nu2FinalSolVec
Double_t maxFitting(Double_t *x, Double_t *par)
int NuPsolutionV3(const double &mNu1, const double &mNu2, const double &phi1, const double &phi2, int &nsol1, int &nsol2)
double MassScale(int method, double mass, const int &tau_type1, const int &tau_type2)
double dTheta3DLimit(const int &tau_type, const int &limit_code, const double &P_tau)
int NuPsolutionLFV(const XYVector &met_vec, const PtEtaPhiMVector &tau, const double &m_nu, std::vector< PtEtaPhiMVector > &nu_vec)
std::vector< PtEtaPhiMVector > m_tauvecsol2
std::vector< PtEtaPhiMVector > m_nu1FinalSolVec
MissingMassCalculator(const MissingMassCalculator &)=delete
int RunMissingMassCalculator(const xAOD::IParticle *part1, const xAOD::IParticle *part2, const xAOD::MissingET *met, const int &njets)
double maxFromHist(const std::shared_ptr< TH1F > &theHist, std::vector< double > &histInfo, const MaxHistStrategy::e maxHistStrategy=MaxHistStrategy::FIT, const int winHalfWidth=2, bool debug=false)
Class providing the definition of the 4-vector interface.
Definition part1.py:1
Definition part2.py:1
MissingET_v1 MissingET
Version control by type defintion.