TEveTrackPropagator.cxx

Go to the documentation of this file.
00001 // @(#)root/eve:$Id: TEveTrackPropagator.cxx 36373 2010-10-19 17:43:35Z matevz $
00002 // Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 
00013 #include "TEveTrackPropagator.h"
00014 #include "TEveTrack.h"
00015 #include "TEveTrans.h"
00016 
00017 #include "TMath.h"
00018 
00019 #include <cassert>
00020 
00021 namespace
00022 {
00023    const Float_t kBMin     = 1e-6;
00024    const Float_t kPtMinSqr = 1e-20;
00025    const Float_t kAMin     = 1e-10;
00026    const Float_t kStepEps  = 1e-3;
00027 }
00028 
00029 //______________________________________________________________________________
00030 TEveTrackPropagator::Helix_t::Helix_t() :
00031    fCharge(0),
00032    fMaxAng(45), fMaxStep(20.f), fDelta(0.1),
00033    fPhi(0), fValid(kFALSE),
00034    fLam(-1), fR(-1), fPhiStep(-1), fSin(-1), fCos(-1),
00035    fRKStep(20.0f),
00036    fPtMag(-1), fPlMag(-1), fLStep(-1)
00037 {
00038    // Default constructor.
00039 }
00040 
00041 //______________________________________________________________________________
00042 void TEveTrackPropagator::Helix_t::UpdateCommon(const TEveVector& p, const TEveVector& b)
00043 {
00044    // Common update code for helix and RK propagation.
00045 
00046    fB = b;
00047 
00048    // base vectors
00049    fE1 = b;
00050    fE1.Normalize();
00051    fPlMag = p.Dot(fE1);
00052    fPl    = fE1*fPlMag;
00053 
00054    fPt    = p - fPl;
00055    fPtMag = fPt.Mag();
00056    fE2    = fPt;
00057    fE2.Normalize();
00058 }
00059 
00060 //______________________________________________________________________________
00061 void TEveTrackPropagator::Helix_t::UpdateHelix(const TEveVector& p, const TEveVector& b,
00062                                                Bool_t full_update, Bool_t enforce_max_step)
00063 {
00064    // Update helix parameters.
00065 
00066    UpdateCommon(p, b);
00067 
00068    // helix parameters
00069    TMath::Cross(fE1.Arr(), fE2.Arr(), fE3.Arr());
00070    if (fCharge < 0) fE3.NegateXYZ();
00071 
00072    if (full_update)
00073    {
00074       using namespace TMath;
00075       Float_t a = fgkB2C * b.Mag() * Abs(fCharge);
00076       if (a > kAMin && fPtMag*fPtMag > kPtMinSqr)
00077       {
00078          fValid = kTRUE;
00079 
00080          fR   = Abs(fPtMag / a);
00081          fLam = fPlMag / fPtMag;
00082 
00083          // get phi step, compare fMaxAng with fDelta
00084          fPhiStep = fMaxAng * DegToRad();
00085          if (fR > fDelta)
00086          {
00087             Float_t ang  = 2.0 * ACos(1.0f - fDelta/fR);
00088             if (ang < fPhiStep)
00089                fPhiStep = ang;
00090          }
00091 
00092          // check max step size
00093          Float_t curr_step = fR * fPhiStep * Sqrt(1.0f + fLam*fLam);
00094          if (curr_step > fMaxStep || enforce_max_step)
00095             fPhiStep *= fMaxStep / curr_step;
00096 
00097          fLStep = fR * fPhiStep * fLam;
00098          fSin   = Sin(fPhiStep);
00099          fCos   = Cos(fPhiStep);
00100       }
00101       else
00102       {
00103          fValid = kFALSE;
00104       }
00105    }
00106 }
00107 
00108 //______________________________________________________________________________
00109 void TEveTrackPropagator::Helix_t::UpdateRK(const TEveVector& p, const TEveVector& b)
00110 {
00111    // Update helix for stepper RungeKutta.
00112 
00113    UpdateCommon(p, b);
00114 
00115    if (fCharge)
00116    {
00117       fValid = kTRUE;
00118 
00119       // cached values for propagator
00120       fB = b;
00121       fPlMag = p.Dot(fB);
00122    }
00123    else
00124    {
00125       fValid = kFALSE;
00126    }
00127 }
00128 
00129 //______________________________________________________________________________
00130 void TEveTrackPropagator::Helix_t::Step(const TEveVector4& v, const TEveVector& p,
00131                                         TEveVector4& vOut, TEveVector& pOut)
00132 {
00133    // Step helix for given momentum p from vertex v.
00134 
00135    vOut = v;
00136 
00137    if (fValid)
00138    {
00139       TEveVector d = fE2*(fR*fSin) + fE3*(fR*(1-fCos)) + fE1*fLStep;
00140       vOut    += d;
00141       vOut.fT += TMath::Abs(fLStep);
00142 
00143       pOut = fPl + fE2*(fPtMag*fCos) + fE3*(fPtMag*fSin);
00144 
00145       fPhi += fPhiStep;
00146    }
00147    else
00148    {
00149       // case: pT < kPtMinSqr or B < kBMin
00150       // might happen if field directon changes pT ~ 0 or B becomes zero
00151       vOut    += p * (fMaxStep / p.Mag());
00152       vOut.fT += fMaxStep;
00153       pOut  = p;
00154    }
00155 }
00156 
00157 
00158 //==============================================================================
00159 // TEveTrackPropagator
00160 //==============================================================================
00161 
00162 //______________________________________________________________________________
00163 //
00164 // Holding structure for a number of track rendering parameters.
00165 // Calculates path taking into account the parameters.
00166 //
00167 // This is decoupled from TEveTrack/TEveTrackList to allow sharing of the
00168 // Propagator among several instances. Back references are kept so the
00169 // tracks can be recreated when the parameters change.
00170 //
00171 // TEveTrackList has Get/Set methods for RnrStlye. TEveTrackEditor and
00172 // TEveTrackListEditor provide editor access.
00173 //
00174 // Specify whether 2D projected tracks get broken into several
00175 // segments when the projected space consists of separate domains
00176 // (like Rho-Z). This is true by default.
00177 
00178 ClassImp(TEveTrackPropagator);
00179 
00180 Float_t             TEveTrackPropagator::fgDefMagField = 0.5;
00181 const Float_t       TEveTrackPropagator::fgkB2C        = 0.299792458e-2;
00182 TEveTrackPropagator TEveTrackPropagator::fgDefault;
00183 
00184 Float_t             TEveTrackPropagator::fgEditorMaxR  = 2000;
00185 Float_t             TEveTrackPropagator::fgEditorMaxZ  = 4000;
00186 
00187 //______________________________________________________________________________
00188 TEveTrackPropagator::TEveTrackPropagator(const char* n, const char* t,
00189                                          TEveMagField *field, Bool_t own_field) :
00190    TEveElementList(n, t),
00191    TEveRefBackPtr(),
00192 
00193    fStepper(kHelix),
00194    fMagFieldObj(field),
00195    fOwnMagFiledObj(own_field),
00196 
00197    fMaxR    (350),   fMaxZ    (450),
00198    fNMax    (4096),  fMaxOrbs (0.5),
00199 
00200    fEditPathMarks (kTRUE),
00201    fFitDaughters  (kTRUE),   fFitReferences (kTRUE),
00202    fFitDecay      (kTRUE),   fFitCluster2Ds (kTRUE),
00203    fRnrDaughters  (kFALSE),  fRnrReferences (kFALSE),
00204    fRnrDecay      (kFALSE),  fRnrCluster2Ds (kFALSE),
00205    fRnrFV         (kFALSE),
00206    fPMAtt(), fFVAtt(),
00207 
00208    fProjTrackBreaking(kPTB_Break), fRnrPTBMarkers(kFALSE), fPTBAtt(),
00209 
00210    fV()
00211 {
00212    // Default constructor.
00213 
00214    fPMAtt.SetMarkerColor(kYellow);
00215    fPMAtt.SetMarkerStyle(2);
00216    fPMAtt.SetMarkerSize(2);
00217 
00218    fFVAtt.SetMarkerColor(kRed);
00219    fFVAtt.SetMarkerStyle(4);
00220    fFVAtt.SetMarkerSize(1.5);
00221 
00222    fPTBAtt.SetMarkerColor(kBlue);
00223    fPTBAtt.SetMarkerStyle(4);
00224    fPTBAtt.SetMarkerSize(0.8);
00225 
00226    if (fMagFieldObj == 0) {
00227       fMagFieldObj = new TEveMagFieldConst(0., 0., fgDefMagField);
00228       fOwnMagFiledObj = kTRUE;
00229    }
00230 }
00231 
00232 //______________________________________________________________________________
00233 TEveTrackPropagator::~TEveTrackPropagator()
00234 {
00235    // Destructor.
00236 
00237    if (fOwnMagFiledObj)
00238    {
00239       delete fMagFieldObj;
00240    }
00241 }
00242 
00243 //______________________________________________________________________________
00244 void TEveTrackPropagator::OnZeroRefCount()
00245 {
00246    // Virtual from TEveRefBackPtr - track reference count has reached zero.
00247 
00248    CheckReferenceCount("TEveTrackPropagator::OnZeroRefCount ");
00249 }
00250 
00251 //______________________________________________________________________________
00252 void TEveTrackPropagator::CheckReferenceCount(const TEveException& eh)
00253 {
00254    // Check reference count - virtual from TEveElement.
00255    // Must also take into account references from TEveRefBackPtr.
00256 
00257    if (fRefCount <= 0)
00258    {
00259       TEveElementList::CheckReferenceCount(eh);
00260    }
00261 }
00262 
00263 //______________________________________________________________________________
00264 void TEveTrackPropagator::ElementChanged(Bool_t update_scenes, Bool_t redraw)
00265 {
00266    // Element-change notification.
00267    // Stamp all tracks as requiring display-list regeneration.
00268    // Virtual from TEveElement.
00269 
00270    TEveTrack* track;
00271    RefMap_i i = fBackRefs.begin();
00272    while (i != fBackRefs.end())
00273    {
00274       track = dynamic_cast<TEveTrack*>(i->first);
00275       track->StampObjProps();
00276       ++i;
00277    }
00278    TEveElementList::ElementChanged(update_scenes, redraw);
00279 }
00280 
00281 //==============================================================================
00282 
00283 //______________________________________________________________________________
00284 void TEveTrackPropagator::InitTrack(TEveVector &v,  Int_t charge)
00285 {
00286    // Initialize internal data-members for given particle parameters.
00287 
00288    fV.fX = v.fX;
00289    fV.fY = v.fY;
00290    fV.fZ = v.fZ;
00291 
00292    fPoints.push_back(fV);
00293 
00294    // init helix
00295    fH.fPhi    = 0;
00296    fH.fCharge = charge;
00297 }
00298 
00299 //______________________________________________________________________________
00300 void TEveTrackPropagator::ResetTrack()
00301 {
00302    // Reset cache holding particle trajectory.
00303 
00304    fPoints.clear();
00305 
00306    // reset helix
00307    fH.fPhi = 0;
00308 }
00309 
00310 //______________________________________________________________________________
00311 Bool_t TEveTrackPropagator::GoToVertex(TEveVector& v, TEveVector& p)
00312 {
00313    // Propagate particle with momentum p to vertex v.
00314 
00315    Update(fV, p, kTRUE);
00316 
00317    if ((v-fV).Mag() < kStepEps)
00318    {
00319       fPoints.push_back(v);
00320       return kTRUE;
00321    }
00322 
00323    return fH.fValid ? LoopToVertex(v, p) : LineToVertex(v);
00324 }
00325 
00326 //______________________________________________________________________________
00327 void TEveTrackPropagator::GoToBounds(TEveVector& p)
00328 {
00329    // Propagate particle to bounds.
00330    // Return TRUE if hit bounds.
00331 
00332    Update(fV, p, kTRUE);
00333 
00334    fH.fValid ? LoopToBounds(p): LineToBounds(p);
00335 }
00336 
00337 //______________________________________________________________________________
00338 void TEveTrackPropagator::Update(const TEveVector4& v, const TEveVector& p,
00339                                  Bool_t full_update, Bool_t enforce_max_step)
00340 {
00341    // Update helix / B-field projection state.
00342 
00343    if (fStepper == kHelix)
00344    {
00345       fH.UpdateHelix(p, fMagFieldObj->GetField(v), !fMagFieldObj->IsConst() || full_update, enforce_max_step);
00346    }
00347    else
00348    {
00349       fH.UpdateRK(p, fMagFieldObj->GetField(v));
00350 
00351       if (full_update)
00352       {
00353          using namespace TMath;
00354 
00355          Float_t a = fgkB2C * fMagFieldObj->GetMaxFieldMag() * Abs(fH.fCharge);
00356          if (a > kAMin)
00357          {
00358             fH.fR = p.Mag() / a;
00359 
00360             // get phi step, compare fDelta with MaxAng
00361             fH.fPhiStep = fH.fMaxAng * DegToRad();
00362             if (fH.fR > fH.fDelta )
00363             {
00364                Float_t ang  = 2.0 * ACos(1.0f - fH.fDelta/fH.fR);
00365                if (ang < fH.fPhiStep)
00366                   fH.fPhiStep = ang;
00367             }
00368 
00369             // check against maximum step-size
00370             fH.fRKStep = fH.fR * fH.fPhiStep * Sqrt(1 + fH.fLam*fH.fLam);
00371             if (fH.fRKStep > fH.fMaxStep || enforce_max_step)
00372             {
00373                fH.fPhiStep *= fH.fMaxStep / fH.fRKStep;
00374                fH.fRKStep   = fH.fMaxStep;
00375             }
00376          }
00377          else
00378          {
00379             fH.fRKStep = fH.fMaxStep; 
00380          }
00381       }
00382    }
00383 }
00384 
00385 //______________________________________________________________________________
00386 void TEveTrackPropagator::Step(const TEveVector4 &v, const TEveVector &p, TEveVector4 &vOut, TEveVector &pOut)
00387 {
00388    // Wrapper to step helix.
00389 
00390    if (fStepper == kHelix)
00391    {
00392       fH.Step(v, p, vOut, pOut);
00393    }
00394    else
00395    {
00396       Double_t vecRKIn[7];
00397       vecRKIn[0] = v.fX;
00398       vecRKIn[1] = v.fY;
00399       vecRKIn[2] = v.fZ;
00400       Double_t pm = p.Mag();
00401       Double_t nm = 1.0 / pm;
00402       vecRKIn[3] = p.fX*nm;
00403       vecRKIn[4] = p.fY*nm;
00404       vecRKIn[5] = p.fZ*nm;
00405       vecRKIn[6] = p.Mag();
00406 
00407       Double_t vecRKOut[7];
00408       StepRungeKutta(fH.fRKStep, vecRKIn, vecRKOut);
00409 
00410       vOut.fX = vecRKOut[0];
00411       vOut.fY = vecRKOut[1];
00412       vOut.fZ = vecRKOut[2];
00413       vOut.fT = v.fT + fH.fRKStep;
00414       pm = vecRKOut[6];
00415       pOut.fX = vecRKOut[3]*pm;
00416       pOut.fY = vecRKOut[4]*pm;
00417       pOut.fZ = vecRKOut[5]*pm;
00418    }
00419 }
00420 
00421 //______________________________________________________________________________
00422 void TEveTrackPropagator::LoopToBounds(TEveVector& p)
00423 {
00424    // Propagate charged particle with momentum p to bounds.
00425    // It is expected that Update() with full-update was called before.
00426 
00427    const Float_t maxRsq = fMaxR*fMaxR;
00428 
00429    TEveVector4 currV(fV);
00430    TEveVector4 forwV(fV);
00431    TEveVector  forwP (p);
00432 
00433    Int_t np = fPoints.size();
00434    Float_t maxPhi = fMaxOrbs*TMath::TwoPi();
00435 
00436    while (fH.fPhi < maxPhi && np<fNMax)
00437    {
00438       Step(currV, p, forwV, forwP);
00439 
00440       // cross R
00441       if (forwV.Perp2() > maxRsq)
00442       {
00443          Float_t t = (fMaxR - currV.R()) / (forwV.R() - currV.R());
00444          if (t < 0 || t > 1)
00445          {
00446             Warning("HelixToBounds", "In MaxR crossing expected t>=0 && t<=1: t=%f, r1=%f, r2=%f, MaxR=%f.",
00447                     t, currV.R(), forwV.R(), fMaxR);
00448             return;
00449          }
00450          TEveVector d(forwV);
00451          d -= currV;
00452          d *= t;
00453          d += currV;
00454          fPoints.push_back(d);
00455          return;
00456       }
00457 
00458       // cross Z
00459       else if (TMath::Abs(forwV.fZ) > fMaxZ)
00460       {
00461          Float_t t = (fMaxZ - TMath::Abs(currV.fZ)) / TMath::Abs((forwV.fZ - currV.fZ));
00462          if (t < 0 || t > 1)
00463          {
00464             Warning("HelixToBounds", "In MaxZ crossing expected t>=0 && t<=1: t=%f, z1=%f, z2=%f, MaxZ=%f.",
00465                     t, currV.fZ, forwV.fZ, fMaxZ);
00466             return;
00467          }
00468          TEveVector d(forwV -currV);
00469          d *= t;
00470          d += currV;
00471          fPoints.push_back(d);
00472          return;
00473       }
00474 
00475       currV = forwV;
00476       p     = forwP;
00477       Update(currV, p);
00478 
00479       fPoints.push_back(currV);
00480       ++np;
00481    }
00482 }
00483 
00484 //______________________________________________________________________________
00485 Bool_t TEveTrackPropagator::LoopToVertex(TEveVector& v, TEveVector& p)
00486 {
00487    // Propagate charged particle with momentum p to vertex v.
00488    // It is expected that Update() with full-update was called before.
00489 
00490    const Float_t maxRsq = fMaxR * fMaxR;
00491 
00492    TEveVector4 currV(fV);
00493    TEveVector4 forwV(fV);
00494    TEveVector  forwP(p);
00495 
00496    Int_t first_point = fPoints.size();
00497    Int_t np          = first_point;
00498 
00499    Float_t prod0=0, prod1;
00500 
00501    do
00502    {
00503       Step(currV, p, forwV, forwP);
00504       Update(forwV, forwP);
00505 
00506       if (PointOverVertex(v, forwV, &prod1))
00507       {
00508          break;
00509       }
00510 
00511       if (IsOutsideBounds(forwV, maxRsq, fMaxZ))
00512       {
00513          fV = currV;
00514          return kFALSE;
00515       }
00516 
00517       fPoints.push_back(forwV);
00518       currV = forwV;
00519       p     = forwP;
00520       prod0 = prod1;
00521       ++np;
00522    } while (np < fNMax);
00523 
00524    // make the remaining fractional step
00525    if (np > first_point)
00526    {
00527       if ((v - currV).Mag() > kStepEps)
00528       {
00529          Float_t step_frac = prod0 / (prod0 - prod1);
00530          if (step_frac > 0)
00531          {
00532             // Step for fraction of previous step size.
00533             // We pass 'enforce_max_step' flag to Update().
00534             Float_t orig_max_step = fH.fMaxStep;
00535             fH.fMaxStep = step_frac * (forwV - currV).Mag();
00536             Update(currV, p, kTRUE, kTRUE);
00537             Step(currV, p, forwV, forwP);
00538             p     = forwP;
00539             currV = forwV;
00540             fPoints.push_back(currV);
00541             ++np;
00542             fH.fMaxStep = orig_max_step;
00543          }
00544 
00545          // Distribute offset to desired crossing point over all segment.
00546 
00547          TEveVector off(v - currV);
00548          off *= 1.0f / currV.fT;
00549 
00550          // Calculate the required momentum rotation.
00551          // lpd - last-points-delta
00552          TEveVector lpd0(fPoints[np-1]);
00553          lpd0 -= fPoints[np-2];
00554          lpd0.Normalize();
00555 
00556          for (Int_t i = first_point; i < np; ++i)
00557          {
00558             fPoints[i] += off * fPoints[i].fT;
00559          }
00560 
00561          TEveVector lpd1(fPoints[np-1]);
00562          lpd1 -= fPoints[np-2];
00563          lpd1.Normalize();
00564 
00565          TEveTrans tt;
00566          tt.SetupFromToVec(lpd0, lpd1);
00567 
00568          // TEveVector pb4(p);
00569          // printf("Rotating momentum: p0 = "); p.Dump();
00570          tt.RotateIP(p);
00571          // printf("                   p1 = "); p.Dump();
00572          // printf("  n1=%f, n2=%f, dp = %f deg\n", pb4.Mag(), p.Mag(),
00573          //        TMath::RadToDeg()*TMath::ACos(p.Dot(pb4)/(pb4.Mag()*p.Mag())));
00574 
00575          fV = v;
00576          return kTRUE;
00577       }
00578    }
00579 
00580    fPoints.push_back(v);
00581    fV = v;
00582    return kTRUE;
00583 }
00584 
00585 //______________________________________________________________________________
00586 Bool_t TEveTrackPropagator::LineToVertex(TEveVector& v)
00587 {
00588    // Propagate neutral particle to vertex v.
00589 
00590    TEveVector4 currV = v;
00591 
00592    currV.fX = v.fX;
00593    currV.fY = v.fY;
00594    currV.fZ = v.fZ;
00595    fPoints.push_back(currV);
00596 
00597    fV = v;
00598    return kTRUE;
00599 }
00600 
00601 //______________________________________________________________________________
00602 void TEveTrackPropagator::LineToBounds(TEveVector& p)
00603 {
00604    // Propagatate neutral particle with momentum p to bounds.
00605 
00606    Float_t tZ = 0, tR = 0, tB = 0;
00607 
00608    // time where particle intersect +/- fMaxZ
00609    if (p.fZ > 0)
00610       tZ = (fMaxZ - fV.fZ) / p.fZ;
00611    else if (p.fZ < 0)
00612       tZ = - (fMaxZ + fV.fZ) / p.fZ;
00613 
00614    // time where particle intersects cylinder
00615    Double_t a = p.fX*p.fX + p.fY*p.fY;
00616    Double_t b = 2.0 * (fV.fX*p.fX + fV.fY*p.fY);
00617    Double_t c = fV.fX*fV.fX + fV.fY*fV.fY - fMaxR*fMaxR;
00618    Double_t d = b*b - 4.0*a*c;
00619    if (d >= 0) {
00620       Double_t sqrtD = TMath::Sqrt(d);
00621       tR = (-b - sqrtD) / (2.0 * a);
00622       if (tR < 0) {
00623          tR = (-b + sqrtD) / (2.0 * a);
00624       }
00625       tB = tR < tZ ? tR : tZ; // compare the two times
00626    } else {
00627       tB = tZ;
00628    }
00629    TEveVector nv(fV.fX + p.fX*tB, fV.fY + p.fY*tB, fV.fZ + p.fZ*tB);
00630    LineToVertex(nv);
00631 }
00632 
00633 //______________________________________________________________________________
00634 Bool_t TEveTrackPropagator::HelixIntersectPlane(const TEveVector& p,
00635                                                 const TEveVector& point,
00636                                                 const TEveVector& normal,
00637                                                 TEveVector& itsect)
00638 {
00639    // Intersect helix with a plane. Current position and argument p define
00640    // the helix.
00641 
00642    TEveVector pos(fV);
00643    TEveVector mom(p);
00644    if (fMagFieldObj->IsConst())
00645       fH.UpdateHelix(mom, fMagFieldObj->GetField(pos), kFALSE, kFALSE);
00646 
00647    TEveVector n(normal);
00648    TEveVector delta = pos - point;
00649    Float_t d = delta.Dot(n);
00650    if (d > 0) {
00651       n.NegateXYZ(); // Turn normal around so that we approach from negative side of the plane
00652       d = -d;
00653    }
00654 
00655    TEveVector4 forwV;
00656    TEveVector  forwP;
00657    TEveVector4 pos4(pos);
00658    while (kTRUE)
00659    {
00660       Update(pos4, mom);
00661       Step(pos4, mom, forwV , forwP);
00662       Float_t new_d = (forwV - point).Dot(n);
00663       if (new_d < d)
00664       {
00665          // We are going further away ... fail intersect.
00666          Warning("HelixIntersectPlane", "going away from the plane.");
00667          return kFALSE;
00668       }
00669       if (new_d > 0)
00670       {
00671          delta = forwV - pos;
00672          itsect = pos + delta * (d / (d - new_d));
00673          return kTRUE;
00674       }
00675       pos4 = forwV;
00676       mom  = forwP;
00677    }
00678 }
00679 
00680 //______________________________________________________________________________
00681 Bool_t TEveTrackPropagator::LineIntersectPlane(const TEveVector& p,
00682                                                const TEveVector& point,
00683                                                const TEveVector& normal,
00684                                                      TEveVector& itsect)
00685 {
00686    // Intersect line with a plane. Current position and argument p define
00687    // the line.
00688 
00689    TEveVector pos(fV.fX, fV.fY, fV.fZ);
00690    TEveVector delta = pos - point;
00691 
00692    Float_t d = delta.Dot(normal);
00693    if (d == 0) {
00694       itsect = pos;
00695       return kTRUE;
00696    }
00697 
00698    Float_t t = (p.Dot(normal)) / d;
00699    if (t < 0) {
00700       return kFALSE;
00701    } else {
00702       itsect = pos + p*t;
00703       return kTRUE;
00704    }
00705 }
00706 
00707 //______________________________________________________________________________
00708 Bool_t TEveTrackPropagator::IntersectPlane(const TEveVector& p,
00709                                            const TEveVector& point,
00710                                            const TEveVector& normal,
00711                                                  TEveVector& itsect)
00712 {
00713    // Find intersection of currently propagated track with a plane.
00714    // Current track position is used as starting point.
00715    //
00716    // Args:
00717    //  p        - track momentum to use for extrapolation
00718    //  point    - a point on a plane
00719    //  normal   - normal of the plane
00720    //  itsect   - output, point of intersection
00721    // Returns:
00722    //  kFALSE if intersection can not be found, kTRUE otherwise.
00723 
00724    if (fH.fCharge && fMagFieldObj && p.Perp2() > kPtMinSqr)
00725       return HelixIntersectPlane(p, point, normal, itsect);
00726    else
00727       return LineIntersectPlane(p, point, normal, itsect);
00728 }
00729 
00730 //______________________________________________________________________________
00731 void TEveTrackPropagator::FillPointSet(TEvePointSet* ps) const
00732 {
00733    // Reset ps and populate it with points in propagation cache.
00734 
00735    Int_t size = TMath::Min(fNMax, (Int_t)fPoints.size());
00736    ps->Reset(size);
00737    for (Int_t i = 0; i < size; ++i)
00738    {
00739       const TEveVector4& v = fPoints[i];
00740       ps->SetNextPoint(v.fX, v.fY, v.fZ);
00741    }
00742 }
00743 
00744 /******************************************************************************/
00745 
00746 //______________________________________________________________________________
00747 void TEveTrackPropagator::RebuildTracks()
00748 {
00749    // Rebuild all tracks using this render-style.
00750 
00751    TEveTrack* track;
00752    RefMap_i i = fBackRefs.begin();
00753    while (i != fBackRefs.end())
00754    {
00755       track = dynamic_cast<TEveTrack*>(i->first);
00756       track->MakeTrack();
00757       track->StampObjProps();
00758       ++i;
00759    }
00760 }
00761 
00762 //______________________________________________________________________________
00763 void TEveTrackPropagator::SetMagField(Float_t bX, Float_t bY, Float_t bZ)
00764 {
00765    // Set constant magnetic field and rebuild tracks.
00766 
00767    SetMagFieldObj(new TEveMagFieldConst(bX, bY, bZ));
00768 }
00769 
00770 //______________________________________________________________________________
00771 void TEveTrackPropagator::SetMagFieldObj(TEveMagField* field, Bool_t own_field)
00772 {
00773    // Set constant magnetic field and rebuild tracks.
00774 
00775    if (fMagFieldObj && fOwnMagFiledObj) delete fMagFieldObj;
00776 
00777    fMagFieldObj    = field;
00778    fOwnMagFiledObj = own_field;
00779 
00780    RebuildTracks();
00781 }
00782 
00783 //______________________________________________________________________________
00784 void TEveTrackPropagator::PrintMagField(Float_t x, Float_t y, Float_t z) const
00785 {
00786    if (fMagFieldObj) fMagFieldObj->PrintField(x, y, z);
00787 }
00788 
00789 //______________________________________________________________________________
00790 void TEveTrackPropagator::SetMaxR(Float_t x)
00791 {
00792    // Set maximum radius and rebuild tracks.
00793 
00794    fMaxR = x;
00795    RebuildTracks();
00796 }
00797 
00798 //______________________________________________________________________________
00799 void TEveTrackPropagator::SetMaxZ(Float_t x)
00800 {
00801    // Set maximum z and rebuild tracks.
00802 
00803    fMaxZ = x;
00804    RebuildTracks();
00805 }
00806 
00807 //______________________________________________________________________________
00808 void TEveTrackPropagator::SetMaxOrbs(Float_t x)
00809 {
00810    // Set maximum number of orbits and rebuild tracks.
00811 
00812    fMaxOrbs = x;
00813    RebuildTracks();
00814 }
00815 
00816 //______________________________________________________________________________
00817 void TEveTrackPropagator::SetMinAng(Float_t x)
00818 {
00819    // Set maximum step angle and rebuild tracks.
00820    // WARNING -- this method / variable was mis-named.
00821 
00822    Warning("SetMinAng", "This method was mis-named, use SetMaxAng() instead!");
00823    SetMaxAng(x);
00824 }
00825 //______________________________________________________________________________
00826 Float_t TEveTrackPropagator::GetMinAng() const
00827 {
00828    // Get maximum step angle.
00829    // WARNING -- this method / variable was mis-named.
00830 
00831    Warning("GetMinAng", "This method was mis-named, use GetMaxAng() instead!");
00832    return GetMaxAng();
00833 }
00834 
00835 //______________________________________________________________________________
00836 void TEveTrackPropagator::SetMaxAng(Float_t x)
00837 {
00838    // Set maximum step angle and rebuild tracks.
00839 
00840    fH.fMaxAng = x;
00841    RebuildTracks();
00842 }
00843 
00844 //______________________________________________________________________________
00845 void TEveTrackPropagator::SetMaxStep(Float_t x)
00846 {
00847    // Set maximum step-size and rebuild tracks.
00848 
00849    fH.fMaxStep = x;
00850    RebuildTracks();
00851 }
00852 
00853 //______________________________________________________________________________
00854 void TEveTrackPropagator::SetDelta(Float_t x)
00855 {
00856    // Set maximum error and rebuild tracks.
00857 
00858    fH.fDelta = x;
00859    RebuildTracks();
00860 }
00861 
00862 //______________________________________________________________________________
00863 void TEveTrackPropagator::SetFitDaughters(Bool_t x)
00864 {
00865    // Set daughter creation point fitting and rebuild tracks.
00866 
00867    fFitDaughters = x;
00868    RebuildTracks();
00869 }
00870 
00871 //______________________________________________________________________________
00872 void TEveTrackPropagator::SetFitReferences(Bool_t x)
00873 {
00874    // Set track-reference fitting and rebuild tracks.
00875 
00876    fFitReferences = x;
00877    RebuildTracks();
00878 }
00879 
00880 //______________________________________________________________________________
00881 void TEveTrackPropagator::SetFitDecay(Bool_t x)
00882 {
00883    // Set decay fitting and rebuild tracks.
00884 
00885    fFitDecay = x;
00886    RebuildTracks();
00887 }
00888 
00889 //______________________________________________________________________________
00890 void TEveTrackPropagator::SetFitCluster2Ds(Bool_t x)
00891 {
00892    // Set 2D-cluster fitting and rebuild tracks.
00893 
00894    fFitCluster2Ds = x;
00895    RebuildTracks();
00896 }
00897 
00898 //______________________________________________________________________________
00899 void TEveTrackPropagator::SetRnrDecay(Bool_t rnr)
00900 {
00901    // Set decay rendering and rebuild tracks.
00902 
00903    fRnrDecay = rnr;
00904    RebuildTracks();
00905 }
00906 
00907 //______________________________________________________________________________
00908 void TEveTrackPropagator::SetRnrCluster2Ds(Bool_t rnr)
00909 {
00910    // Set rendering of 2D-clusters and rebuild tracks.
00911 
00912    fRnrCluster2Ds = rnr;
00913    RebuildTracks();
00914 }
00915 
00916 //______________________________________________________________________________
00917 void TEveTrackPropagator::SetRnrDaughters(Bool_t rnr)
00918 {
00919    // Set daughter rendering and rebuild tracks.
00920 
00921    fRnrDaughters = rnr;
00922    RebuildTracks();
00923 }
00924 
00925 //______________________________________________________________________________
00926 void TEveTrackPropagator::SetRnrReferences(Bool_t rnr)
00927 {
00928    // Set track-reference rendering and rebuild tracks.
00929 
00930    fRnrReferences = rnr;
00931    RebuildTracks();
00932 }
00933 
00934 //______________________________________________________________________________
00935 void TEveTrackPropagator::SetRnrFV(Bool_t x)
00936 {
00937    // Set first-vertex rendering and rebuild tracks.
00938 
00939    fRnrFV = x;
00940    RebuildTracks();
00941 }
00942 
00943 //______________________________________________________________________________
00944 void TEveTrackPropagator::SetProjTrackBreaking(UChar_t x)
00945 {
00946    // Set projection break-point mode and rebuild tracks.
00947 
00948    fProjTrackBreaking = x;
00949    RebuildTracks();
00950 }
00951 
00952 //______________________________________________________________________________
00953 void TEveTrackPropagator::SetRnrPTBMarkers(Bool_t x)
00954 {
00955    // Set projection break-point rendering and rebuild tracks.
00956 
00957    fRnrPTBMarkers = x;
00958    RebuildTracks();
00959 }
00960 
00961 //______________________________________________________________________________
00962 void TEveTrackPropagator::StepRungeKutta(Double_t step,
00963                                          Double_t* vect, Double_t* vout)
00964 {
00965   // Wrapper to step with method RungeKutta.
00966 
00967   ///   ******************************************************************
00968   ///   *                                                                *
00969   ///   *  Runge-Kutta method for tracking a particle through a magnetic *
00970   ///   *  field. Uses Nystroem algorithm (See Handbook Nat. Bur. of     *
00971   ///   *  Standards, procedure 25.5.20)                                 *
00972   ///   *                                                                *
00973   ///   *  Input parameters                                              *
00974   ///   *       CHARGE    Particle charge                                *
00975   ///   *       STEP      Step size                                      *
00976   ///   *       VECT      Initial co-ords,direction cosines,momentum     *
00977   ///   *  Output parameters                                             *
00978   ///   *       VOUT      Output co-ords,direction cosines,momentum      *
00979   ///   *  User routine called                                           *
00980   ///   *       CALL GUFLD(X,F)                                          *
00981   ///   *                                                                *
00982   ///   *    ==>Called by : <USER>, GUSWIM                               *
00983   ///   *       Authors    R.Brun, M.Hansroul  *********                 *
00984   ///   *                  V.Perevoztchikov (CUT STEP implementation)    *
00985   ///   *                                                                *
00986   ///   *                                                                *
00987   ///   ******************************************************************
00988 
00989   Double_t h2, h4, f[4];
00990   Double_t xyzt[3], a, b, c, ph,ph2;
00991   Double_t secxs[4],secys[4],seczs[4],hxp[3];
00992   Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
00993   Double_t est, at, bt, ct, cba;
00994   Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
00995 
00996   Double_t x;
00997   Double_t y;
00998   Double_t z;
00999 
01000   Double_t xt;
01001   Double_t yt;
01002   Double_t zt;
01003 
01004   // const Int_t maxit = 1992;
01005   const Int_t maxit  = 500;
01006   const Int_t maxcut = 11;
01007 
01008   const Double_t hmin   = 1e-4; // !!! MT ADD,  should be member
01009   const Double_t kdlt   = 1e-3; // !!! MT CHANGE from 1e-4, should be member
01010   const Double_t kdlt32 = kdlt/32.;
01011   const Double_t kthird = 1./3.;
01012   const Double_t khalf  = 0.5;
01013   const Double_t kec    = 2.9979251e-3;
01014 
01015   const Double_t kpisqua = 9.86960440109;
01016   const Int_t kix  = 0;
01017   const Int_t kiy  = 1;
01018   const Int_t kiz  = 2;
01019   const Int_t kipx = 3;
01020   const Int_t kipy = 4;
01021   const Int_t kipz = 5;
01022 
01023   // *.
01024   // *.    ------------------------------------------------------------------
01025   // *.
01026   // *             this constant is for units cm,gev/c and kgauss
01027   // *
01028   Int_t iter = 0;
01029   Int_t ncut = 0;
01030   for(Int_t j = 0; j < 7; j++)
01031     vout[j] = vect[j];
01032 
01033   Double_t pinv   = kec * fH.fCharge / vect[6];
01034   Double_t tl     = 0.;
01035   Double_t h      = step;
01036   Double_t rest;
01037 
01038   do {
01039     rest  = step - tl;
01040     if (TMath::Abs(h) > TMath::Abs(rest))
01041        h = rest;
01042 
01043     f[0] = -fH.fB.fX;
01044     f[1] = -fH.fB.fY;
01045     f[2] = -fH.fB.fZ;
01046 
01047     // * start of integration
01048     x      = vout[0];
01049     y      = vout[1];
01050     z      = vout[2];
01051     a      = vout[3];
01052     b      = vout[4];
01053     c      = vout[5];
01054 
01055     h2     = khalf * h;
01056     h4     = khalf * h2;
01057     ph     = pinv * h;
01058     ph2    = khalf * ph;
01059     secxs[0] = (b * f[2] - c * f[1]) * ph2;
01060     secys[0] = (c * f[0] - a * f[2]) * ph2;
01061     seczs[0] = (a * f[1] - b * f[0]) * ph2;
01062     ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
01063     if (ang2 > kpisqua) break;
01064 
01065     dxt    = h2 * a + h4 * secxs[0];
01066     dyt    = h2 * b + h4 * secys[0];
01067     dzt    = h2 * c + h4 * seczs[0];
01068     xt     = x + dxt;
01069     yt     = y + dyt;
01070     zt     = z + dzt;
01071 
01072     // * second intermediate point
01073     est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
01074     if (est > h) {
01075       if (ncut++ > maxcut) break;
01076       h *= khalf;
01077       continue;
01078     }
01079 
01080     xyzt[0] = xt;
01081     xyzt[1] = yt;
01082     xyzt[2] = zt;
01083 
01084     fH.fB = fMagFieldObj->GetField(xt, yt, zt);
01085     f[0] = -fH.fB.fX;
01086     f[1] = -fH.fB.fY;
01087     f[2] = -fH.fB.fZ;
01088 
01089     at     = a + secxs[0];
01090     bt     = b + secys[0];
01091     ct     = c + seczs[0];
01092 
01093     secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
01094     secys[1] = (ct * f[0] - at * f[2]) * ph2;
01095     seczs[1] = (at * f[1] - bt * f[0]) * ph2;
01096     at     = a + secxs[1];
01097     bt     = b + secys[1];
01098     ct     = c + seczs[1];
01099     secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
01100     secys[2] = (ct * f[0] - at * f[2]) * ph2;
01101     seczs[2] = (at * f[1] - bt * f[0]) * ph2;
01102     dxt    = h * (a + secxs[2]);
01103     dyt    = h * (b + secys[2]);
01104     dzt    = h * (c + seczs[2]);
01105     xt     = x + dxt;
01106     yt     = y + dyt;
01107     zt     = z + dzt;
01108     at     = a + 2.*secxs[2];
01109     bt     = b + 2.*secys[2];
01110     ct     = c + 2.*seczs[2];
01111 
01112     est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
01113     if (est > 2.*TMath::Abs(h)) {
01114       if (ncut++ > maxcut) break;
01115       h *= khalf;
01116       continue;
01117     }
01118 
01119     xyzt[0] = xt;
01120     xyzt[1] = yt;
01121     xyzt[2] = zt;
01122 
01123     fH.fB = fMagFieldObj->GetField(xt, yt, zt);
01124     f[0] = -fH.fB.fX;
01125     f[1] = -fH.fB.fY;
01126     f[2] = -fH.fB.fZ;
01127 
01128     z      = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
01129     y      = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
01130     x      = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
01131 
01132     secxs[3] = (bt*f[2] - ct*f[1])* ph2;
01133     secys[3] = (ct*f[0] - at*f[2])* ph2;
01134     seczs[3] = (at*f[1] - bt*f[0])* ph2;
01135     a      = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
01136     b      = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
01137     c      = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
01138 
01139     est    = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
01140       + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
01141       + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
01142 
01143     if (est > kdlt && TMath::Abs(h) > hmin) {
01144       if (ncut++ > maxcut) break;
01145       h *= khalf;
01146       continue;
01147     }
01148 
01149     ncut = 0;
01150     // * if too many iterations, go to helix
01151     if (iter++ > maxit) break;
01152 
01153     tl += h;
01154     if (est < kdlt32)
01155       h *= 2.;
01156     cba    = 1./ TMath::Sqrt(a*a + b*b + c*c);
01157     vout[0] = x;
01158     vout[1] = y;
01159     vout[2] = z;
01160     vout[3] = cba*a;
01161     vout[4] = cba*b;
01162     vout[5] = cba*c;
01163     rest = step - tl;
01164     if (step < 0.) rest = -rest;
01165     if (rest < 1.e-5*TMath::Abs(step))
01166     {
01167        Float_t dot = (vout[3]*vect[3] + vout[4]*vect[4] + vout[5]*vect[5]);
01168        fH.fPhi += TMath::ACos(dot);
01169        return;
01170     }
01171 
01172   } while(1);
01173 
01174   // angle too big, use helix
01175 
01176   f1  = f[0];
01177   f2  = f[1];
01178   f3  = f[2];
01179   f4  = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
01180   rho = -f4*pinv;
01181   tet = rho * step;
01182 
01183   hnorm = 1./f4;
01184   f1 = f1*hnorm;
01185   f2 = f2*hnorm;
01186   f3 = f3*hnorm;
01187 
01188   hxp[0] = f2*vect[kipz] - f3*vect[kipy];
01189   hxp[1] = f3*vect[kipx] - f1*vect[kipz];
01190   hxp[2] = f1*vect[kipy] - f2*vect[kipx];
01191 
01192   hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
01193 
01194   rho1 = 1./rho;
01195   sint = TMath::Sin(tet);
01196   cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
01197 
01198   g1 = sint*rho1;
01199   g2 = cost*rho1;
01200   g3 = (tet-sint) * hp*rho1;
01201   g4 = -cost;
01202   g5 = sint;
01203   g6 = cost * hp;
01204 
01205   vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
01206   vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
01207   vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
01208 
01209   vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
01210   vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
01211   vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
01212 
01213   fH.fPhi += tet;
01214 }

Generated on Tue Jul 5 14:16:14 2011 for ROOT_528-00b_version by  doxygen 1.5.1