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
96
99 double m_nsigma_METscan_lfv_lh{}, m_beamEnergy{}; // number of sigmas for MET-scan
100
103
104 double m_prob_tmp{};
105
107 double m_mtautauSum{};
108
110 int m_seed{};
111 // data member for the spaceWalker approach
112
113 int m_iter0{};
119 int m_nosol1{};
120 int m_nosol2{};
122 bool m_switch1{};
123 bool m_switch2{};
124
126
132
136
141
142
143 double m_mTau{},m_mTau2{};
145 double m_eTau1{}, m_eTau2{};
146 double m_eTau10{}, m_eTau20{};
155
161 std::vector<double> m_probFinalSolOldVec;
162 std::vector<double> m_mtautauFinalSolOldVec;
163 std::vector<PtEtaPhiMVector> m_nu1FinalSolOldVec;
164 std::vector<PtEtaPhiMVector> m_nu2FinalSolOldVec;
165
166 int m_nsol{};
167 std::vector<double> m_probFinalSolVec;
168 std::vector<double> m_mtautauFinalSolVec;
169 std::vector<PtEtaPhiMVector> m_nu1FinalSolVec;
170 std::vector<PtEtaPhiMVector> m_nu2FinalSolVec;
171
172
175 double m_walkWeight{};
177
179
180 PtEtaPhiMVector m_tauVec1,m_tauVec2;
186 double m_tauVec1E{};
187 double m_tauVec2E{};
188 double m_m2Nu1{};
189 double m_m2Nu2{};
190 double m_ET2v1{};
191 double m_ET2v2{};
192 double m_E2v1{};
193 double m_E2v2{};
194 double m_Ev2{};
195 double m_Ev1{};
196 double m_Mvis{},m_Meff{};
198
199 //--- define histograms for histogram method
200 //--- upper limits need to be revisied in the future!!! It may be not enough for some analyses
201 std::shared_ptr<TH1F> m_fMfit_all;
202 std::shared_ptr<TH1F> m_fMEtP_all;
203 std::shared_ptr<TH1F> m_fMEtL_all;
204 std::shared_ptr<TH1F> m_fMnu1_all;
205 std::shared_ptr<TH1F> m_fMnu2_all;
206 std::shared_ptr<TH1F> m_fPhi1_all;
207 std::shared_ptr<TH1F> m_fPhi2_all;
208 std::shared_ptr<TGraph> m_fMfit_allGraph;
209 std::shared_ptr<TH1F> m_fMfit_allNoWeight;
210
211 std::shared_ptr<TH1F> m_fPXfit1;
212 std::shared_ptr<TH1F> m_fPYfit1;
213 std::shared_ptr<TH1F> m_fPZfit1;
214 std::shared_ptr<TH1F> m_fPXfit2;
215 std::shared_ptr<TH1F> m_fPYfit2;
216 std::shared_ptr<TH1F> m_fPZfit2;
217
218 // these histograms are used for the floating stopping criterion
219 std::shared_ptr<TH1F> m_fMmass_split1;
220 std::shared_ptr<TH1F> m_fMEtP_split1;
221 std::shared_ptr<TH1F> m_fMEtL_split1;
222 std::shared_ptr<TH1F> m_fMnu1_split1;
223 std::shared_ptr<TH1F> m_fMnu2_split1;
224 std::shared_ptr<TH1F> m_fPhi1_split1;
225 std::shared_ptr<TH1F> m_fPhi2_split1;
226 std::shared_ptr<TH1F> m_fMmass_split2;
227 std::shared_ptr<TH1F> m_fMEtP_split2;
228 std::shared_ptr<TH1F> m_fMEtL_split2;
229 std::shared_ptr<TH1F> m_fMnu1_split2;
230 std::shared_ptr<TH1F> m_fMnu2_split2;
231 std::shared_ptr<TH1F> m_fPhi1_split2;
232 std::shared_ptr<TH1F> m_fPhi2_split2;
233
235
236 TH1F* m_fPhi1{};
237 TH1F* m_fPhi2{};
238 TH1F* m_fMnu1{};
239 TH1F* m_fMnu2{};
240 TH1F* m_fMetx{};
241 TH1F* m_fMety{};
242 TH1F* m_fTheta3D{};
243 TH1F* m_fTauProb{};
244
245 // for intermediate calc
246 PtEtaPhiMVector m_TLVdummy;
247
248 //---------------- protected variables
249 DitauStuff m_fDitauStuffFit; // results based on fit method
250 DitauStuff m_fDitauStuffHisto; // results based on histo method
251
252 int m_niter_fit1{}; // number of iterations for dR-dPhi scan
253 int m_niter_fit2{}; // number of iterations for MET-scan
254 int m_niter_fit3{}; // number of iterations for Mnu-scan
255 int m_NiterRandom{}; // number of random iterations (for lh, multiply or divide by 10 for ll and hh)
258 int m_RndmSeedAltering{}; // reset seed (not necessary by default)
259
260 double m_dTheta3d_binMin{}; // minimal step size for dTheta3D
261 double m_dTheta3d_binMax{}; // maximum step size for dTheta3D
262 double m_dRmax_tau{}; // maximum dR(nu-visTau)
263
264 double m_MnuScanRange{}; // range of M(nunu) scan; M(nunu) range can be affected by selection cuts
265
266
267 //---------------- protected functions
268 void ClearDitauStuff(DitauStuff &fStuff);
269 void DoOutputInfo();
270 void PrintOtherInput();
271 void PrintResults();
272
273 inline int NuPsolutionV3(const double & mNu1, const double & mNu2, const double & phi1, const double & phi2,
274 int & nsol1, int & nsol2);
275
276 inline int NuPsolutionLFV(const XYVector & met_vec, const PtEtaPhiMVector & tau,
277 const double & m_nu, std::vector<PtEtaPhiMVector> &nu_vec);
278
279
280 protected:
281 inline int CheckSolutions(PtEtaPhiMVector nu_vec, PtEtaPhiMVector vis_vec, int decayType);
282 inline int TailCleanUp(const PtEtaPhiMVector & vis1, const PtEtaPhiMVector & nu1,
283 const PtEtaPhiMVector & vis2, const PtEtaPhiMVector & nu2,
284 const double & mmc_mass, const double & vis_mass, const double & eff_mass, const double & dphiTT);
285
286
287 inline int refineSolutions ( const double & M_nu1, const double & M_nu2,
288 const int nsol1, const int nsol2,
289 const double & Mvis, const double & Meff);
290
291
292
293 inline void handleSolutions();
294
295 inline double MassScale(int method, double mass, const int & tau_type1, const int & tau_type2);
296
297
298 // factor out parameter space walking and probablity computing
299 inline int DitauMassCalculatorV9walk();
300
301
302 // Calculates mass of lep+tau system in LFV X->lep+tau decays
303 // It is based on DitauMassCalculatorV9, not optimized for speed yet, simple phase-space scan
304 inline int DitauMassCalculatorV9lfv(bool refit);
305
306
307
308 // only compute probability
309 inline int probCalculatorV9fast(
310 const double & phi1, const double & phi2,
311 const double & M_nu1, const double & M_nu2);
312
313 // initialize the walker
314 inline void SpaceWalkerInit();
315
316 //walk the walker
317 inline bool SpaceWalkerWalk();
318
319 inline bool precomputeCache();
320
321
322 inline bool checkMEtInRange () ;
323 inline bool checkAllParamInRange () ;
324
325 //----------------------------------------------
326 //
327 // ------------ Public methods ---------------
328 //
329 //______________________________________________
330
331public:
332
334
335 MissingMassCalculator(MMCCalibrationSet::e aset, std::string paramFilePath) ;
336
339
343
344 int RunMissingMassCalculator( const xAOD::IParticle* part1, const xAOD::IParticle* part2, const xAOD::MissingET* met, const int& njets );
345
346 //-------- Set Input Parameters
347 void FinalizeSettings(const xAOD::IParticle* part1, const xAOD::IParticle* part2, const xAOD::MissingET* met, const int& njets );
348 void SetNiterFit1(const int val) { m_niter_fit1=val; } // number of iterations per loop in dPhi loop
349 void SetNiterFit2(const int val) { m_niter_fit2=val; } // number of iterations per loop in MET loop
350 void SetNiterFit3(const int val) { m_niter_fit3=val; } // number of iterations per loop in Mnu loop
351 void SetNiterRandom(const int val) { m_NiterRandom=val; } // number of random iterations
352 void SetNsucStop(const int val) { m_NsucStop=val; } // Arrest criteria for Nsuccesses
353 void SetRMSStop(const int val) { m_RMSStop=val;}
354 void SetMeanbinStop(const double val) {m_meanbinStop=val;}
355 void SetRndmSeedAltering(const int val) { m_RndmSeedAltering=val; } // number of iterations per loop in Mnu loop
356 void SetdTheta3d_binMax(const double val) { m_dTheta3d_binMax=val; } // maximum step size for dTheta3D
357 void SetdTheta3d_binMin(const double val) { m_dTheta3d_binMin=val; } // minimal step size for dTheta3D
358 void SetEventNumber(const int eventNumber) { m_eventNumber = eventNumber; }
359
360 void SetMnuScanRange(const double val) { m_MnuScanRange=val; }
361
362 void SetProposalTryMEt(const double val) {m_proposalTryMEt=val; }
363 void SetProposalTryPhi(const double val) {m_ProposalTryPhi=val;}
364 void SetProposalTryMnu(const double val) {m_ProposalTryMnu=val;}
365 void SetProposalTryEtau(const double val) {m_ProposalTryEtau=val;}
366
369
370 int GetNiterFit1() const { return m_niter_fit1; } // number of iterations per loop in dPhi loop
371 int GetNiterFit2() const { return m_niter_fit2; } // number of iterations per loop in MET loop
372 int GetNiterFit3() const { return m_niter_fit3; } // number of iterations per loop in Mnu loop
373 int GetNiterRandom() const { return m_niterRandomLocal; } // number of random iterations
374
375 int GetNsucStop() const { return m_NsucStop; } // Arrest criteria for NSuc
376 int GetRMSStop() const { return m_RMSStop; }
377 double GetMeanbinStop() const { return m_meanbinStop;}
378 int GetRndmSeedAltering() const { return m_RndmSeedAltering; } // number of iterations per loop in Mnu loop
379
383 int GetMarkovNAccept() const { return m_markovNAccept; }
385 double GetProposalTryMEt() const {return m_proposalTryMEt;}
386 double GetProposalTryPhi() const {return m_ProposalTryPhi;}
387 double GetProposalTryMnu() const {return m_ProposalTryMnu;}
388 double GetProposalTryEtau() const {return m_ProposalTryEtau;}
389
390 void SetNsigmaMETscan_ll(const double val) { m_nsigma_METscan_ll=val; } // number of sigma's for MET-scan in ll events
391 void SetNsigmaMETscan_lh(const double val) { m_nsigma_METscan_lh=val; } // number of sigma's for MET-scan in lh events
392 void SetNsigmaMETscan_hh(const double val) { m_nsigma_METscan_hh=val; } // number of sigma's for MET-scan in hh events
393 void SetNsigmaMETscan(const double val) { m_nsigma_METscan=val; } // number of sigma's for MET-scan
394
395 void SetUseFloatStopping(const bool val); // switch for floating stopping criterion
398 void SetFloatStoppingComp(const double val) { m_fUseFloatStoppingComp = val;}
399 void SetBeamEnergy(const double val) { m_beamEnergy=val; }
400 void SetLFVLeplepRefit(const bool val) { m_lfvLeplepRefit=val; }
401 void SaveLlhHisto(const bool val);
402
403 double GetmMaxError() const {return m_PrintmMaxError;}
404 double GetmMeanError() const { return m_PrintmMeanError;}
406
407 int GetNNoSol() const {return m_markovNRejectNoSol;}
409 int GetNSol() const {return m_markovNAccept;}
410
411
412 //-------- Get results;
413 Double_t maxFitting(Double_t *x, Double_t *par);
414
415 // compute maximum from histo
416 double maxFromHist(TH1F *theHist, std::vector<double> & histInfo, const MaxHistStrategy::e maxHistStrategy=MaxHistStrategy::FIT,const int winHalfWidth=2,bool debug=false);
417 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) {
418 return maxFromHist(theHist.get(), histInfo, maxHistStrategy, winHalfWidth, debug);
419 }
420
421 XYVector metvec_tmp;
422 inline double dTheta3DLimit(const int & tau_type, const int & limit_code,const double & P_tau);
423
424};
425} // namespace DiTauMassTools
426
427#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.