CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TauA1NuConstrainedFitter.cc
Go to the documentation of this file.
1 /* From SimpleFits Package
2  * Designed an written by
3  * author: Ian M. Nugent
4  * Humboldt Foundations
5  */
6 #include <functional>
10 #include <iostream>
11 
12 using namespace tauImpactParameter;
13 
14 TauA1NuConstrainedFitter::TauA1NuConstrainedFitter(unsigned int ambiguity,const LorentzVectorParticle& A1,const TVector3& PVertex, const TMatrixTSym<double>& VertexCov):
16  ambiguity_(ambiguity)
17 {
18  TLorentzVector Tau(0,0,0,0);
19  //dummy substitution not used later
20  TVectorT<double> nu_par(LorentzVectorParticle::NLorentzandVertexPar,1);
21  TMatrixTSym<double> nu_cov(LorentzVectorParticle::NLorentzandVertexPar);
22  LorentzVectorParticle Nu(nu_par,nu_cov,PDGInfo::nu_tau,0.0,A1.bField());
23  particles_.push_back(A1);
24  particles_.push_back(Nu);
25 
26  // setup 13 by 13 matrix
28  TVectorT<double> inpar(size);
29  TMatrixTSym<double> incov(size);
30 
31  // Get primary vertex information
32  if(VertexCov.GetNrows()!=LorentzVectorParticle::NVertex)return;
33  inpar(LorentzVectorParticle::vx)=PVertex.X();
34  inpar(LorentzVectorParticle::vy)=PVertex.Y();
35  inpar(LorentzVectorParticle::vz)=PVertex.Z();
36  for(int i=0; i<LorentzVectorParticle::NVertex;i++){
37  for(int j=0; j<LorentzVectorParticle::NVertex;j++)incov(i,j)=VertexCov(i,j);
38  }
39  int A1offset=LorentzVectorParticle::NVertex;
42  inpar(i+A1offset)=A1.parameter(i);
43  inpar(i+Nuoffset)=Nu.parameter(i)+1.0;// offset by 1 GeV to prevent convergence on first iteration
45  incov(i+A1offset,j+A1offset)=A1.covariance(i,j);
46  incov(i+Nuoffset,j+Nuoffset)=Nu.covariance(i,j);
47  }
48  }
49 
50  exppar.ResizeTo(nexpandedpar,1);
54 
55  TVectorT<double> PAR_0(npar);
56  par_0.ResizeTo(npar);
57  cov_0.ResizeTo(npar,npar);
59  for(int i=0; i<npar;i++)par_0(i)=PAR_0(i);
61 
62  for(int i=0; i<npar;i++){
63  for(int j=0;j<npar;j++){cov_0(i,j)=expcov(i,j);}
64  }
65 
66  par.ResizeTo(npar);
67  par=par_0;
68  cov.ResizeTo(npar,npar);
69  cov=cov_0;
70 }
71 
72 TVectorT<double> TauA1NuConstrainedFitter::ComputeInitalExpPar(const TVectorT<double>& inpar){
73  TVectorT<double> outpar(nexpandedpar);
76  TVector3 sv(inpar(LorentzVectorParticle::vx+offset),inpar(LorentzVectorParticle::vy+offset),inpar(LorentzVectorParticle::vz+offset));
77  TVector3 TauDir=sv-pv;
78  outpar(tau_phi)=TauDir.Phi();
79  outpar(tau_theta)=TauDir.Theta();
80  outpar(a1_px)=inpar(LorentzVectorParticle::px+offset);
81  outpar(a1_py)=inpar(LorentzVectorParticle::py+offset);
82  outpar(a1_pz)=inpar(LorentzVectorParticle::pz+offset);
83  outpar(a1_m)=inpar(LorentzVectorParticle::m+offset);
84  outpar(a1_vx)=inpar(LorentzVectorParticle::vx+offset);
85  outpar(a1_vy)=inpar(LorentzVectorParticle::vy+offset);
86  outpar(a1_vz)=inpar(LorentzVectorParticle::vz+offset);
88  outpar(nu_px)=inpar(LorentzVectorParticle::px+offset);
89  outpar(nu_py)=inpar(LorentzVectorParticle::py+offset);
90  outpar(nu_pz)=inpar(LorentzVectorParticle::pz+offset);
91  return outpar;
92 }
93 
94 TVectorT<double> TauA1NuConstrainedFitter::ComputeExpParToPar(const TVectorT<double>& inpar){
95  TVectorT<double> outpar(npar);
96  for(int i=0;i<npar;i++){outpar(i)=inpar(i);}
97  return outpar;
98 }
99 
100 TVectorT<double> TauA1NuConstrainedFitter::ComputeNuLorentzVectorPar(const TVectorT<double>& inpar){
101  TVectorT<double> outpar(LorentzVectorParticle::NLorentzandVertexPar);
102  outpar(LorentzVectorParticle::vx)=inpar(a1_vx);
103  outpar(LorentzVectorParticle::vy)=inpar(a1_vy);
104  outpar(LorentzVectorParticle::vz)=inpar(a1_vz);
105  outpar(LorentzVectorParticle::px)=inpar(nu_px);
106  outpar(LorentzVectorParticle::py)=inpar(nu_py);
107  outpar(LorentzVectorParticle::pz)=inpar(nu_pz);
108  outpar(LorentzVectorParticle::m)=0;
109  return outpar;
110 }
111 
112 TVectorT<double> TauA1NuConstrainedFitter::ComputeA1LorentzVectorPar(const TVectorT<double>& inpar){
113  TVectorT<double> outpar(LorentzVectorParticle::NLorentzandVertexPar);
114  outpar(LorentzVectorParticle::vx)=inpar(a1_vx);
115  outpar(LorentzVectorParticle::vy)=inpar(a1_vy);
116  outpar(LorentzVectorParticle::vz)=inpar(a1_vz);
117  outpar(LorentzVectorParticle::px)=inpar(a1_px);
118  outpar(LorentzVectorParticle::py)=inpar(a1_py);
119  outpar(LorentzVectorParticle::pz)=inpar(a1_pz);
120  outpar(LorentzVectorParticle::m)=inpar(a1_m);
121  return outpar;
122 }
123 
124 TVectorT<double> TauA1NuConstrainedFitter::ComputeMotherLorentzVectorPar(const TVectorT<double>& inpar){
125  TVectorT<double> outpar(LorentzVectorParticle::NLorentzandVertexPar);
126  TVectorT<double> nupar=ComputeNuLorentzVectorPar(inpar);
127  TVectorT<double> a1par=ComputeA1LorentzVectorPar(inpar);
129  if(i<LorentzVectorParticle::NVertex){outpar(i)=a1par(i);}
130  else{outpar(i)=nupar(i)+a1par(i);}
131  //if(i==LorentzVectorParticle::m) outpar(i,0)=PDGInfo::tau_mass();
132  }
133  double nu_px = nupar(LorentzVectorParticle::px);
134  double nu_py = nupar(LorentzVectorParticle::py);
135  double nu_pz = nupar(LorentzVectorParticle::pz);
136  double Enu2 = nu_px*nu_px + nu_py*nu_py + nu_pz*nu_pz;
137  double a1_px = a1par(LorentzVectorParticle::px);
138  double a1_py = a1par(LorentzVectorParticle::py);
139  double a1_pz = a1par(LorentzVectorParticle::pz);
140  double a1_m = a1par(LorentzVectorParticle::m);
141  double Ea12 = a1_px*a1_px + a1_py*a1_py + a1_pz*a1_pz + a1_m*a1_m;
142  double outpar_px = outpar(LorentzVectorParticle::px);
143  double outpar_py = outpar(LorentzVectorParticle::py);
144  double outpar_pz = outpar(LorentzVectorParticle::pz);
145  double P2=outpar_px*outpar_px + outpar_py*outpar_py + outpar_pz*outpar_pz;
146  outpar(LorentzVectorParticle::m)=sqrt(fabs(Enu2 + Ea12 + 2*sqrt(Enu2*Ea12)-P2));
147  return outpar;
148 }
149 
151  // assumes changes to a1 correlation to vertex is small
152  if(par.GetNrows()==npar && cov.GetNrows()==npar && exppar.GetNrows()==npar && expcov.GetNrows()==npar) return;
153  for(int i=0; i<npar;i++){
154  exppar(i)=par(i);
155  for(int j=0; j<npar;j++){expcov(i,j)=cov(i,j);}
156  }
157 }
158 
159 std::vector<LorentzVectorParticle> TauA1NuConstrainedFitter::getRefitDaughters(){
160  std::vector<LorentzVectorParticle> refitParticles;
162  double c(0),b(0);
163  for(unsigned int i=0;i<particles_.size();i++){c+=particles_[i].charge();b=particles_[i].bField();}
164  TVectorT<double> a1=ComputeA1LorentzVectorPar(exppar);
166  refitParticles.push_back(LorentzVectorParticle(a1,a1cov,particles_[0].pdgId(),c,b));
167  TVectorT<double> nu=ComputeNuLorentzVectorPar(exppar);
169  refitParticles.push_back(LorentzVectorParticle(nu,nucov,PDGInfo::nu_tau,0.0,b));
170  return refitParticles;
171 }
172 
175  double c(0),b(0);
176  for(unsigned int i=0;i<particles_.size();i++){c+=particles_[i].charge();b=particles_[i].bField();}
177  TVectorT<double> m=ComputeMotherLorentzVectorPar(exppar);
179  LorentzVectorParticle mymother= LorentzVectorParticle(m,mcov,(int)(-1.0*fabs(PDGInfo::tau_minus)*c),c,b);
180  return mymother;
181 }
182 
183 void TauA1NuConstrainedFitter::CovertParToObjects(const TVectorD &v,TLorentzVector &a1,TLorentzVector &nu,double &phi,double &theta,TVector3 &TauDir){
184  a1=TLorentzVector(v(a1_px),v(a1_py),v(a1_pz),sqrt(v(a1_m)*v(a1_m)+v(a1_px)*v(a1_px)+v(a1_py)*v(a1_py)+v(a1_pz)*v(a1_pz)));
185  nu=TLorentzVector(v(nu_px),v(nu_py),v(nu_pz),sqrt(v(nu_px)*v(nu_px)+v(nu_py)*v(nu_py)+v(nu_pz)*v(nu_pz)));
186  phi=v(tau_phi);
187  theta=v(tau_theta);
188  TauDir.SetMagThetaPhi(1.0,theta,phi);
189 }
190 
193  // Check if Tau Direction is unphysical and if nessicary set the starting point to Theta_{GJ-Max}
194  TLorentzVector a1(par(a1_px),par(a1_py),par(a1_pz),sqrt(par(a1_m)*par(a1_m)+par(a1_px)*par(a1_px)+par(a1_py)*par(a1_py)+par(a1_pz)*par(a1_pz)));
195  double phi(par(tau_phi)),theta(par(tau_theta));
196  TLorentzVector Tau_plus,Tau_minus,nu_plus,nu_minus;
197  TVector3 TauDir(cos(phi)*sin(theta),sin(phi)*sin(theta),cos(theta));
198  bool isReal;
199  solveByRotation(TauDir,a1,Tau_plus,Tau_minus,nu_plus,nu_minus,isReal);
200 
201  //check that the do product of the a1 and tau is positive, otherwise there is no information for tau direction -> use zero solution
202  if(TauDir.Dot(a1.Vect())<0){
203  isReal=false;
204  }
205 
206  //case 1: is real then solve analytically
207  if(isReal && (ambiguity_==plus || ambiguity_==minus)){
208  // popogate errors
211  for(int i=0; i<npar;i++) par(i)=par_tmp(i);
212  return true;
213  }
214  // case 2 is in unphsyical region - rotate and substitue \theta_{GJ} with \theta_{GJ}^{Max} and then solve analytically
215  else if(!isReal && ambiguity_==zero){
218  for(int i=0; i<npar;i++) par(i)=par_tmp(i);
219  return true;
220  }
221  return false;
222 }
223 
224 TVectorT<double> TauA1NuConstrainedFitter::SolveAmbiguityAnalytically(const TVectorT<double>& inpar, unsigned int amb){
225  // Solve equation quadratic equation
226  TVectorT<double> outpar(inpar.GetNrows());
227  TLorentzVector a1,nu;
228  double phi(0),theta(0);
229  TVector3 TauDir;
230  CovertParToObjects(inpar,a1,nu,phi,theta,TauDir);
231  TLorentzVector a1_d=a1;
232  TLorentzVector nu_d=nu;
233  TLorentzVector Tau_plus,Tau_minus,nu_plus,nu_minus;
234  bool isReal;
235  solveByRotation(TauDir,a1_d,Tau_plus,Tau_minus,nu_plus,nu_minus,isReal,true);
236  if(amb==plus)nu=nu_plus;
237  else nu=nu_minus;
238  for(int i=0; i<outpar.GetNrows();i++){ outpar(i)=inpar(i);}
239  outpar(nu_px)=nu.Px();
240  outpar(nu_py)=nu.Py();
241  outpar(nu_pz)=nu.Pz();
242  return outpar;
243 }
244 
245 TVectorT<double> TauA1NuConstrainedFitter::SolveAmbiguityAnalyticallywithRot(const TVectorT<double>& inpar, unsigned int ambiguity){
246  // Rotate and subsitute \theta_{GJ} with \theta_{GJ}^{Max} - assumes uncertianty on thata and phi of the a1 or small compared to the tau direction.
247  TVectorT<double> outpar(inpar.GetNrows());
248  TLorentzVector a1,nu;
249  double phi(0),theta(0);
250  TVector3 TauDir;
251  CovertParToObjects(inpar,a1,nu,phi,theta,TauDir);
252  double theta_a1(a1.Theta()),phi_a1(a1.Phi()),theta_GJMax(thetaGJMax(a1));
253  TauDir.RotateZ(-phi_a1);
254  TauDir.RotateY(-theta_a1);
255  double phiprime(TauDir.Phi());
256  TauDir=TVector3(sin(theta_GJMax)*cos(phiprime),sin(theta_GJMax)*sin(phiprime),cos(theta_GJMax));
257  TauDir.RotateY(theta_a1);
258  TauDir.RotateZ(phi_a1);
259  for(int i=0; i<outpar.GetNrows();i++) outpar(i)=inpar(i);
260  outpar(tau_phi)=TauDir.Phi();
261  outpar(tau_theta)=TauDir.Theta();
262  return SolveAmbiguityAnalytically(outpar,ambiguity);
263 }
264 
265 // Return the significance of the rotation when the tau direction is in the unphysical region
267  TVectorT<double> par_tmp=TauA1NuConstrainedFitter::TauRot(par);
269  if(!(cov_tmp(0,0)>0)) return -999; // return invalid value if the covariance is unphysical (safety flag)
270  if(par_tmp(0)>0) return par_tmp(0)/sqrt(cov_tmp(0,0)); // return the significance if the value is in the unphysical region
271  return 0; // reutrn 0 for the rotation significance if the tau is in the physical region
272 }
273 
274 
275 TVectorT<double> TauA1NuConstrainedFitter::TauRot(const TVectorT<double>& inpar){
276  TVectorT<double> outpar(1);
277  TLorentzVector a1,nu;
278  double phi(0),theta(0);
279  TVector3 TauDir;
280  CovertParToObjects(inpar,a1,nu,phi,theta,TauDir);
281  double theta_a1(a1.Theta()),phi_a1(a1.Phi()),theta_GJMax(thetaGJMax(a1));
282  TauDir.RotateZ(-phi_a1);
283  TauDir.RotateY(-theta_a1);
284  outpar(0)=(TauDir.Theta()-theta_GJMax);
285  return outpar;
286 }
virtual double covariance(int i, int j) const
Definition: Particle.h:25
int i
Definition: DBlmapReader.cc:9
static const std::string Nu
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
static TVectorT< double > ComputeInitalExpPar(const TVectorT< double > &inpar)
static double thetaGJMax(const TLorentzVector &a1)
static TVectorT< double > ComputeA1LorentzVectorPar(const TVectorT< double > &inpar)
static TVectorT< double > ComputeMotherLorentzVectorPar(const TVectorT< double > &inpar)
std::vector< LorentzVectorParticle > getRefitDaughters()
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TauA1NuConstrainedFitter(unsigned int ambiguity, const LorentzVectorParticle &A1, const TVector3 &PVertex, const TMatrixTSym< double > &VertexCov)
int j
Definition: DBlmapReader.cc:9
static TVectorT< double > TauRot(const TVectorT< double > &inpar)
static TVectorT< double > SolveAmbiguityAnalyticallywithRot(const TVectorT< double > &inpar, unsigned int ambiguity)
std::vector< LorentzVectorParticle > particles_
double b
Definition: hdecay.h:120
Geom::Phi< T > phi() const
static void CovertParToObjects(const TVectorD &v, TLorentzVector &a1, TLorentzVector &nu, double &phi, double &theta, TVector3 &TauDir)
static TVectorT< double > ComputeNuLorentzVectorPar(const TVectorT< double > &inpar)
double amb
Definition: hdecay.h:21
virtual double bField() const
Definition: Particle.h:26
static TVectorT< double > ComputeExpParToPar(const TVectorT< double > &inpar)
static TVectorT< double > SolveAmbiguityAnalytically(const TVectorT< double > &inpar, unsigned int ambiguity)
tuple size
Write out results.
static void solveByRotation(const TVector3 &TauDir, const TLorentzVector &A1, TLorentzVector &Tau_plus, TLorentzVector &Tau_minus, TLorentzVector &nu_plus, TLorentzVector &nu_minus, bool &isReal, bool rotateback=true)
static TMatrixTSym< double > propagateError(std::function< TVectorT< double >(const TVectorT< double > &)> f, const TVectorT< double > &inPar, TMatrixTSym< double > &inCov, double epsilon=0.001, double errorEpsilonRatio=1000)