31 #include <TDecompSVD.h>
33 #include "RestFrames/MinMassesSqInvJigsaw.hh"
36 namespace RestFrames {
39 const std::string& stitle,
42 m_Npair(N_vis_inv_pair)
44 m_InvMassDependancy =
true;
50 MinMassesSqInvJigsaw::~MinMassesSqInvJigsaw() {}
53 return MinMassesSqInvJigsaw::m_Empty;
65 TLorentzVector PvisTOT(0.,0.,0.,0.);
67 for(
int i = 0; i < m_Npair; i++){
68 m_Pvis.push_back(GetDependancyStates(i+m_Npair).GetFourVector());
75 TVector3 Boost = PvisTOT.BoostVector();
78 for(
int i = 0; i < m_Npair; i++)
84 for(
int i = 0; i < m_Npair; i++){
85 m_Pvis[i].Boost(-Boost);
90 ECM += sqrt(Pvis*Pvis + Minv*Minv);
98 return SetSpirit(
false);
104 double Minv = INV.M();
106 TLorentzVector VIS(0.,0.,0.,0.);
109 for(
int i = 0; i < m_Npair; i++){
110 m_Pvis.push_back(GetDependancyStates(i).GetFourVector());
112 if(m_Pvis[i].E() > 1e-6)
116 TLorentzVector VISTOT = VIS;
119 TVector3 BoostCM = (VIS+INV).BoostVector();
121 if(BoostCM.Mag() >= 1 || Nnonzero < m_Npair){
123 Pinv.SetVectM(INV.Vect(), std::max(0.,GetChildState(0).
GetMinimumMass()));
125 for(
int i = 1; i < m_Npair; i++){
126 Pinv.SetVectM(TVector3(), std::max(0.,GetChildState(i).
GetMinimumMass()));
129 return SetSpirit(
true);
134 for(
int i = 0; i < m_Npair; i++)
135 m_Pvis[i].Boost(-BoostCM);
137 TVector3 BoostINV = INV.BoostVector();
138 TVector3 BoostVIS = VIS.BoostVector();
140 if(BoostINV.Mag() >= 1 || BoostVIS.Mag() >= 1){
142 Pinv.SetVectM(INV.Vect(), std::max(0.,GetChildState(0).
GetMinimumMass()));
145 for(
int i = 1; i < m_Npair; i++){
146 Pinv.SetVectM(TVector3(), std::max(0.,GetChildState(i).
GetMinimumMass()));
149 return SetSpirit(
true);
155 for(
int i = 0; i < m_Npair; i++){
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]);
163 for(
int i = 0; i < m_Npair; i++)
164 m_Pvis[i].Boost(-BoostINV);
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]);
172 ApplyOptimalRotation();
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]);
180 for(
int i = 0; i < m_Npair; i++){
181 m_Pinv[i].Boost(BoostINV);
182 m_Pinv[i].Boost(BoostCM);
186 return SetSpirit(
true);
190 double MinMassesSqInvJigsaw::GetPScale(
double Minv){
191 std::vector<double> Pinv2;
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]);
200 if(fk > Minv || Ek <= 0.)
return 0.;
206 while(fabs(dk2) >= 1e-10 && count < 100){
209 for(
int i = 0; i < m_Npair; i++){
210 Ek = sqrt(m_Minv[i]*m_Minv[i]+k2*Pinv2[i]);
212 dfk += Ek > 0 ? Pinv2[i]/Ek : 0.;
219 return sqrt(std::max(0.,k2));
223 void MinMassesSqInvJigsaw::ApplyOptimalRotation(){
225 TVector3 Z(0.,0.,0.);
228 while(Z.Mag() <= 0. && index < m_Npair){
229 Z = m_Pinv[index].Vect().Cross(m_Pvis[index].Vect());
232 if(Z.Mag() <= 0.)
return;
236 for(
int i = 0; i < m_Npair; i++){
237 if(Z.Dot(m_Pvis[i].Vect().Unit()) > 1e-8){
241 if(Z.Dot(m_Pinv[i].Vect().Unit()) > 1e-8){
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;
255 TVector3 Y = Z.Cross(X).Unit();
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());
265 double theta = atan2(num,den);
266 for(
int i = 0; i < m_Npair; i++)
267 m_Pinv[i].Rotate(theta, -Z);
275 for(
int i = 0; i < 3; i++)
276 for(
int j = 0; j < 3; j++)
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();
293 TMatrixD UT(TMatrixD::kTransposed,SVD.GetU());
294 TMatrixD R(SVD.GetV(),TMatrixD::kMult,UT);
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());
308 MinMassesSqInvJigsaw MinMassesSqInvJigsaw::m_Empty;