Fork me on GitHub

Ignore:
File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/ParticlePropagator.cc

    r764f12a0 r187fc41  
    6666{
    6767}
     68
    6869
    6970//------------------------------------------------------------------------------
     
    9192  fItInputArray = fInputArray->MakeIterator();
    9293
     94  // import beamspot
     95 
     96  fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
     97
    9398  // create output arrays
    9499
     
    111116{
    112117  Candidate *candidate, *mother;
    113   TLorentzVector candidatePosition, candidateMomentum;
     118  TLorentzVector candidatePosition, candidateMomentum, beamSpotPosition;
    114119  Double_t px, py, pz, pt, pt2, e, q;
    115120  Double_t x, y, z, t, r, phi;
     
    120125  Double_t tmp, discr, discr2;
    121126  Double_t delta, gammam, omega, asinrho;
    122   Double_t rcu, rc2, dxy, xd, yd, zd;
    123 
     127  Double_t rcu, rc2, xd, yd, zd;
     128  Double_t l, d0, dz, p, ctgTheta, phip, etap, alpha;
     129  Double_t bsx, bsy, bsz;
     130         
    124131  const Double_t c_light = 2.99792458E8;
    125 
     132 
     133  if (!fBeamSpotInputArray->GetSize () || !fBeamSpotInputArray->At(0))
     134    beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
     135  else
     136  {
     137    Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At(0));
     138    beamSpotPosition = beamSpotCandidate.Position;
     139  }
     140 
    126141  fItInputArray->Reset();
    127142  while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
     
    132147    y = candidatePosition.Y()*1.0E-3;
    133148    z = candidatePosition.Z()*1.0E-3;
     149   
     150    bsx = beamSpotPosition.X()*1.0E-3;
     151    bsy = beamSpotPosition.Y()*1.0E-3;
     152    bsz = beamSpotPosition.Z()*1.0E-3;
     153   
    134154    q = candidate->Charge;
    135155
     
    181201      y_t = y + py*t;
    182202      z_t = z + pz*t;
     203     
     204      l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
    183205
    184206      mother = candidate;
     
    186208
    187209      candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*e*1.0E3);
    188 
     210      candidate->L = l*1.0E3;
     211   
    189212      candidate->Momentum = candidateMomentum;
    190213      candidate->AddCandidate(mother);
     
    238261      yd = (rc2 > 0.0) ? yd / rc2 : -999;
    239262      zd = z + (TMath::Sqrt(xd*xd + yd*yd) - TMath::Sqrt(x*x + y*y))*pz/pt;
    240 
    241       // calculate impact paramater
    242       dxy = (xd*py - yd*px)/pt;
    243 
     263     
     264      // use perigee momentum rather than original particle
     265      // momentum, since the orignal particle momentum isn't known
     266     
     267      px = TMath::Sign(1.0,r) * pt * (-y_c / r_c);
     268      py = TMath::Sign(1.0,r) * pt * (x_c / r_c);
     269      etap = candidateMomentum.Eta();
     270      phip = TMath::ATan2(py, px);
     271       
     272      candidateMomentum.SetPtEtaPhiE(pt, etap, phip, candidateMomentum.E());
     273     
     274      // calculate additional track parameters (correct for beamspot position)
     275       
     276      d0        = (  (x - bsx) * py - (y - bsy) * px) / pt;
     277      dz        = z - ((x - bsx) * px + (py - bsy) * py) / pt * (pz / pt);
     278      p         = candidateMomentum.P();
     279      ctgTheta  = 1.0 / TMath::Tan (candidateMomentum.Theta ());
     280         
    244281      // 3. time evaluation t = TMath::Min(t_r, t_z)
    245282      //    t_r : time to exit from the sides
     
    287324      r_t = TMath::Hypot(x_t, y_t);
    288325
     326     
     327      // compute path length for an helix
     328     
     329      alpha = pz*1.0E9 / c_light / gammam;
     330      l = t * TMath::Sqrt(alpha*alpha + r*r*omega*omega);
     331         
    289332      if(r_t > 0.0)
    290333      {
     
    295338
    296339        candidate->Momentum = candidateMomentum;
    297         candidate->Dxy = dxy*1.0E3;
    298         candidate->Xd = xd*1.0E3;
     340       
     341        candidate->L =  l*1.0E3;
     342        candidate->D0 = d0*1.0E3;
     343        candidate->DZ = dz*1.0E3;
     344        candidate->P = p;
     345        candidate->PT = pt;
     346        candidate->CtgTheta = ctgTheta;
     347        candidate->Phi = phi;
     348
     349        candidate->Xd = xd*1.0E3;
    299350        candidate->Yd = yd*1.0E3;
    300351        candidate->Zd = zd*1.0E3;
Note: See TracChangeset for help on using the changeset viewer.