LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
MaxProbBreitWignerInvJigsaw.hh
Go to the documentation of this file.
1 // RestFrames: particle physics event analysis library
3 // --------------------------------------------------------------------
4 // Copyright (c) 2014-2016, Christopher Rogan
14 // This file is part of RestFrames.
15 //
16 // RestFrames is free software; you can redistribute it and/or modify
17 // it under the terms of the GNU General Public License as published by
18 // the Free Software Foundation; either version 2 of the License, or
19 // (at your option) any later version.
20 //
21 // RestFrames is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with RestFrames. If not, see <http://www.gnu.org/licenses/>.
29 
30 #ifndef MaxProbBreitWignerInvJigsaw_HH
31 #define MaxProbBreitWignerInvJigsaw_HH
32 
33 #include "Math/Minimizer.h"
34 #include "Math/Functor.h"
35 
37 
38 namespace RestFrames {
39 
41  public:
42 
48  MaxProbBreitWignerInvJigsaw(const std::string& sname,
49  const std::string& stitle,
50  int N_vis_inv_pair);
52 
56  virtual std::string GetLabel() const { return "Max Prob Breit-Wigner"; }
57 
62  virtual void SetPoleMass(double mass, int i = 0);
63 
68  virtual void SetWidth(double width, int i = 0);
69 
76  virtual double GetMinimumMass() const;
77 
83  virtual bool AnalyzeEvent();
84 
85  private:
86  const int m_Npair;
87  mutable std::vector<TLorentzVector> m_Pvis;
88  mutable std::vector<TLorentzVector> m_Pinv;
89  std::vector<double> m_Minv;
90 
91  std::vector<double> m_Mass;
92  std::vector<double> m_Width;
93 
94  bool m_2D;
95  TVector3 m_X;
96  TVector3 m_Y;
97  TVector3 m_Z;
98 
99  double GetPScale(double Minv);
100  void ApplyOptimalRotation();
101 
102  ROOT::Math::Minimizer* m_minimizer;
103  ROOT::Math::Functor* m_functor;
104 
105  double EvaluateMetric(const double* param);
106 
107  };
108 
109 }
110 
111 #endif
RestFrames::MaxProbBreitWignerInvJigsaw::GetMinimumMass
virtual double GetMinimumMass() const
Returns minimum Lorentz invariant mass.
Definition: MaxProbBreitWignerInvJigsaw.cc:66
InvisibleJigsaw.hh
RestFrames::MaxProbBreitWignerInvJigsaw::GetLabel
virtual std::string GetLabel() const
Returns name of this Jigsaw.
Definition: MaxProbBreitWignerInvJigsaw.hh:56
RestFrames::MaxProbBreitWignerInvJigsaw::AnalyzeEvent
virtual bool AnalyzeEvent()
Analyzes event for this jigsaw.
Definition: MaxProbBreitWignerInvJigsaw.cc:94
RestFrames::MaxProbBreitWignerInvJigsaw::SetWidth
virtual void SetWidth(double width, int i=0)
Sets width of a frame.
Definition: MaxProbBreitWignerInvJigsaw.cc:308
RestFrames::MaxProbBreitWignerInvJigsaw::MaxProbBreitWignerInvJigsaw
MaxProbBreitWignerInvJigsaw(const std::string &sname, const std::string &stitle, int N_vis_inv_pair)
Standard constructor.
Definition: MaxProbBreitWignerInvJigsaw.cc:36
RestFrames::MaxProbBreitWignerInvJigsaw
Definition: MaxProbBreitWignerInvJigsaw.hh:40
RestFrames::MaxProbBreitWignerInvJigsaw::SetPoleMass
virtual void SetPoleMass(double mass, int i=0)
Sets pole mass of a frame.
Definition: MaxProbBreitWignerInvJigsaw.cc:301
RestFrames::InvisibleJigsaw
Definition: InvisibleJigsaw.hh:39