LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
MinMassesSqInvJigsaw.cc
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 #include <TMatrixD.h>
31 #include <TDecompSVD.h>
32 
33 #include "RestFrames/MinMassesSqInvJigsaw.hh"
35 
36 namespace RestFrames {
37 
38  MinMassesSqInvJigsaw::MinMassesSqInvJigsaw(const std::string& sname,
39  const std::string& stitle,
40  int N_vis_inv_pair) :
41  InvisibleJigsaw(sname, stitle, N_vis_inv_pair, N_vis_inv_pair),
42  m_Npair(N_vis_inv_pair)
43  {
44  m_InvMassDependancy = true;
45  }
46 
48  : InvisibleJigsaw(), m_Npair(0) {}
49 
50  MinMassesSqInvJigsaw::~MinMassesSqInvJigsaw() {}
51 
53  return MinMassesSqInvJigsaw::m_Empty;
54  }
55 
57  if(!IsSoundMind())
58  return 0.;
59 
60  if(m_Npair <= 0)
61  return 0.;
62  if(m_Npair == 1)
63  return std::max(0.,GetChildState(0).GetMinimumMass());
64 
65  TLorentzVector PvisTOT(0.,0.,0.,0.);
66  m_Pvis.clear();
67  for(int i = 0; i < m_Npair; i++){
68  m_Pvis.push_back(GetDependancyStates(i+m_Npair).GetFourVector());
69  PvisTOT += m_Pvis[i];
70  }
71 
72  double ECM = 0.;
73  double Minv;
74  double Pvis;
75  TVector3 Boost = PvisTOT.BoostVector();
76 
77  if(Boost.Mag() >= 1){
78  for(int i = 0; i < m_Npair; i++)
79  ECM += std::max(0.,GetChildState(i).GetMinimumMass());
80 
81  return ECM;
82  }
83 
84  for(int i = 0; i < m_Npair; i++){
85  m_Pvis[i].Boost(-Boost);
86  Pvis = m_Pvis[i].P();
87 
88  Minv = std::max(0.,GetChildState(i).GetMinimumMass());
89 
90  ECM += sqrt(Pvis*Pvis + Minv*Minv);
91  }
92 
93  return ECM;
94  }
95 
97  if(!IsSoundMind())
98  return SetSpirit(false);
99 
100  if(m_Npair <= 1)
101  return false;
102 
103  TLorentzVector INV = GetParentState().GetFourVector();
104  double Minv = INV.M();
105 
106  TLorentzVector VIS(0.,0.,0.,0.);
107  m_Pvis.clear();
108  int Nnonzero = 0;
109  for(int i = 0; i < m_Npair; i++){
110  m_Pvis.push_back(GetDependancyStates(i).GetFourVector());
111  VIS += m_Pvis[i];
112  if(m_Pvis[i].E() > 1e-6)
113  Nnonzero++;
114  }
115 
116  TLorentzVector VISTOT = VIS;
117 
118  // go to INV+VIS CM frame
119  TVector3 BoostCM = (VIS+INV).BoostVector();
120 
121  if(BoostCM.Mag() >= 1 || Nnonzero < m_Npair){
122  TLorentzVector Pinv;
123  Pinv.SetVectM(INV.Vect(), std::max(0.,GetChildState(0).GetMinimumMass()));
124  GetChildState(0).SetFourVector(Pinv);
125  for(int i = 1; i < m_Npair; i++){
126  Pinv.SetVectM(TVector3(), std::max(0.,GetChildState(i).GetMinimumMass()));
127  GetChildState(i).SetFourVector(Pinv);
128  }
129  return SetSpirit(true);
130  }
131 
132  INV.Boost(-BoostCM);
133  VIS.Boost(-BoostCM);
134  for(int i = 0; i < m_Npair; i++)
135  m_Pvis[i].Boost(-BoostCM);
136 
137  TVector3 BoostINV = INV.BoostVector();
138  TVector3 BoostVIS = VIS.BoostVector();
139 
140  if(BoostINV.Mag() >= 1 || BoostVIS.Mag() >= 1){
141  TLorentzVector Pinv;
142  Pinv.SetVectM(INV.Vect(), std::max(0.,GetChildState(0).GetMinimumMass()));
143  Pinv.Boost(BoostCM);
144  GetChildState(0).SetFourVector(Pinv);
145  for(int i = 1; i < m_Npair; i++){
146  Pinv.SetVectM(TVector3(), std::max(0.,GetChildState(i).GetMinimumMass()));
147  GetChildState(i).SetFourVector(Pinv);
148  }
149  return SetSpirit(true);
150  }
151 
152  // INV states defined in the VIS frame
153  m_Pinv.clear();
154  m_Minv.clear();
155  for(int i = 0; i < m_Npair; i++){
156  m_Minv.push_back(std::max(0.,GetChildState(i).GetMinimumMass()));
157  m_Pinv.push_back(m_Pvis[i]);
158  m_Pinv[i].Boost(-BoostVIS);
159  m_Pinv[i].SetVectM(m_Pinv[i].Vect(), m_Minv[i]);
160  }
161 
162  // VIS states in INV frame
163  for(int i = 0; i < m_Npair; i++)
164  m_Pvis[i].Boost(-BoostINV);
165 
166  if(m_Npair == 2){
167  TVector3 Vdiff = (m_Pvis[0].Vect()-m_Pvis[1].Vect()).Unit();
168  double pinv = GetP(Minv, m_Minv[0], m_Minv[1]);
169  m_Pinv[0].SetVectM( pinv*Vdiff, m_Minv[0]);
170  m_Pinv[1].SetVectM(-pinv*Vdiff, m_Minv[1]);
171  } else {
172  ApplyOptimalRotation();
173 
174  double k = GetPScale(Minv);
175  for(int i = 0; i < m_Npair; i++)
176  m_Pinv[i].SetVectM(k*m_Pinv[i].Vect(), m_Minv[i]);
177  }
178 
179  // return to original frame
180  for(int i = 0; i < m_Npair; i++){
181  m_Pinv[i].Boost(BoostINV);
182  m_Pinv[i].Boost(BoostCM);
183  GetChildState(i).SetFourVector(m_Pinv[i]);
184  }
185 
186  return SetSpirit(true);
187  }
188 
189  // Based on Newton-Raphson root finding
190  double MinMassesSqInvJigsaw::GetPScale(double Minv){
191  std::vector<double> Pinv2;
192  double Ek = 0.;
193  double fk = 0.;
194  double dfk = 0.;
195  for(int i = 0; i < m_Npair; i++){
196  Pinv2.push_back(m_Pinv[i].P()*m_Pinv[i].P());
197  Ek += sqrt(m_Minv[i]*m_Minv[i]+Pinv2[i]);
198  fk += m_Minv[i];
199  }
200  if(fk > Minv || Ek <= 0.) return 0.;
201 
202  double k2 = Minv/Ek;
203  k2 *= k2;
204  double dk2 = k2;
205  int count = 0;
206  while(fabs(dk2) >= 1e-10 && count < 100){
207  fk = -Minv;
208  dfk = 0.;
209  for(int i = 0; i < m_Npair; i++){
210  Ek = sqrt(m_Minv[i]*m_Minv[i]+k2*Pinv2[i]);
211  fk += Ek;
212  dfk += Ek > 0 ? Pinv2[i]/Ek : 0.;
213  }
214  dk2 = 2.*fk/dfk;
215  k2 -= dk2;
216  count++;
217  }
218 
219  return sqrt(std::max(0.,k2));
220  }
221 
222  // Optimal rotation found using Singular Value Decomposition
223  void MinMassesSqInvJigsaw::ApplyOptimalRotation(){
224  // first, check dimensionality of points
225  TVector3 Z(0.,0.,0.);
226  int index = 0;
227 
228  while(Z.Mag() <= 0. && index < m_Npair){
229  Z = m_Pinv[index].Vect().Cross(m_Pvis[index].Vect());
230  index++;
231  }
232  if(Z.Mag() <= 0.) return; // already aligned
233  Z = Z.Unit();
234 
235  bool b_2D = true;
236  for(int i = 0; i < m_Npair; i++){
237  if(Z.Dot(m_Pvis[i].Vect().Unit()) > 1e-8){
238  b_2D = false;
239  break;
240  }
241  if(Z.Dot(m_Pinv[i].Vect().Unit()) > 1e-8){
242  b_2D = false;
243  break;
244  }
245  }
246 
247  // two dimensional problem, one rotation
248  if(b_2D){
249  TVector3 X(0.,0.,0.);
250  for(int i = 0; i < m_Npair; i++)
251  X += m_Pvis[i].Vect() - m_Pinv[i].Vect();
252  if(X.Mag() <= 0) return; // can't improve
253 
254  X = X.Unit();
255  TVector3 Y = Z.Cross(X).Unit();
256 
257  double num = 0.;
258  double den = 0.;
259  for(int i = 0; i < m_Npair; i++){
260  num += m_Pvis[i].Vect().Dot(Y)*
261  (m_Pvis[i].Vect().Dot(X)-m_Pinv[i].Vect().Dot(X));
262  den += m_Pvis[i].Vect().Dot(m_Pinv[i].Vect());
263  }
264 
265  double theta = atan2(num,den);
266  for(int i = 0; i < m_Npair; i++)
267  m_Pinv[i].Rotate(theta, -Z);
268 
269  return;
270  }
271 
272  // three dimensional problem - R from SVD
273  TMatrixD H(3,3);
274  double val;
275  for(int i = 0; i < 3; i++)
276  for(int j = 0; j < 3; j++)
277  H(i,j) = 0.;
278 
279  for(int p = 0; p < m_Npair; p++){
280  H(0,0) += m_Pinv[p].Px()*m_Pvis[p].Px();
281  H(1,0) += m_Pinv[p].Px()*m_Pvis[p].Py();
282  H(2,0) += m_Pinv[p].Px()*m_Pvis[p].Pz();
283  H(0,1) += m_Pinv[p].Py()*m_Pvis[p].Px();
284  H(1,1) += m_Pinv[p].Py()*m_Pvis[p].Py();
285  H(2,1) += m_Pinv[p].Py()*m_Pvis[p].Pz();
286  H(0,2) += m_Pinv[p].Pz()*m_Pvis[p].Px();
287  H(1,2) += m_Pinv[p].Pz()*m_Pvis[p].Py();
288  H(2,2) += m_Pinv[p].Pz()*m_Pvis[p].Pz();
289  }
290 
291  TDecompSVD SVD(H);
292  SVD.Decompose();
293  TMatrixD UT(TMatrixD::kTransposed,SVD.GetU());
294  TMatrixD R(SVD.GetV(),TMatrixD::kMult,UT);
295 
296  TVector3 V;
297  for(int p = 0; p < m_Npair; p++){
298  V = m_Pinv[p].Vect();
299  V.SetXYZ(R(0,0)*V.Px()+R(1,0)*V.Py()+R(2,0)*V.Pz(),
300  R(0,1)*V.Px()+R(1,1)*V.Py()+R(2,1)*V.Pz(),
301  R(0,2)*V.Px()+R(1,2)*V.Py()+R(2,2)*V.Pz());
302  m_Pinv[p].SetVectM(V, m_Pinv[p].M());
303  }
304 
305  return;
306  }
307 
308  MinMassesSqInvJigsaw MinMassesSqInvJigsaw::m_Empty;
309 
310 }
RestFrames::MinMassesSqInvJigsaw::Empty
static MinMassesSqInvJigsaw & Empty()
Returns empty MinMassesSqInvJigsaw.
Definition: MinMassesSqInvJigsaw.cc:52
InvisibleState.hh
RestFrames::State::SetFourVector
virtual void SetFourVector(const TLorentzVector &V)
Sets four-vector of this frame.
Definition: State.cc:161
RestFrames::State::GetFourVector
virtual TLorentzVector GetFourVector() const
Returns four vector of this frame.
Definition: State.cc:165
RestFrames::MinMassesSqInvJigsaw
Definition: MinMassesSqInvJigsaw.hh:37
RestFrames::MinMassesSqInvJigsaw::GetMinimumMass
virtual double GetMinimumMass() const
Returns minimum Lorentz invariant mass.
Definition: MinMassesSqInvJigsaw.cc:56
RestFrames::MinMassesSqInvJigsaw::MinMassesSqInvJigsaw
MinMassesSqInvJigsaw()
Empty constructor.
Definition: MinMassesSqInvJigsaw.cc:47
RestFrames::InvisibleJigsaw
Definition: InvisibleJigsaw.hh:39
RestFrames::MinMassesSqInvJigsaw::AnalyzeEvent
virtual bool AnalyzeEvent()
Analyzes event for this jigsaw.
Definition: MinMassesSqInvJigsaw.cc:96