TEveTrack.cxx

Go to the documentation of this file.
00001 // @(#)root/eve:$Id: TEveTrack.cxx 38366 2011-03-10 15:05:31Z 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 #include "TEveTrack.h"
00013 #include "TEveTrackPropagator.h"
00014 #include "TEvePointSet.h"
00015 #include "TEveVSDStructs.h"
00016 
00017 #include "TPolyLine3D.h"
00018 #include "TMarker.h"
00019 #include "TPolyMarker3D.h"
00020 #include "TColor.h"
00021 #include "TParticlePDG.h"
00022 
00023 // Updates
00024 #include "TEveManager.h"
00025 #include "TEveBrowser.h"
00026 #include "TEveTrackProjected.h"
00027 
00028 #include "Riostream.h"
00029 
00030 #include <vector>
00031 #include <algorithm>
00032 #include <functional>
00033 
00034 //==============================================================================
00035 //==============================================================================
00036 // TEveTrack
00037 //==============================================================================
00038 
00039 //______________________________________________________________________________
00040 //
00041 // Visual representation of a track.
00042 //
00043 
00044 ClassImp(TEveTrack);
00045 
00046 //______________________________________________________________________________
00047 TEveTrack::TEveTrack() :
00048    TEveLine(),
00049 
00050    fV(),
00051    fP(),
00052    fPEnd(),
00053    fBeta(0),
00054    fPdg(0),
00055    fCharge(0),
00056    fLabel(kMinInt),
00057    fIndex(kMinInt),
00058    fStatus(0),
00059    fLockPoints(kFALSE),
00060    fPathMarks(),
00061    fLastPMIdx(0),
00062    fPropagator(0)
00063 {
00064    // Default constructor.
00065 }
00066 
00067 //______________________________________________________________________________
00068 TEveTrack::TEveTrack(TParticle* t, Int_t label, TEveTrackPropagator* prop):
00069    TEveLine(),
00070 
00071    fV(t->Vx(), t->Vy(), t->Vz()),
00072    fP(t->Px(), t->Py(), t->Pz()),
00073    fPEnd(),
00074    fBeta(t->P()/t->Energy()),
00075    fPdg(0),
00076    fCharge(0),
00077    fLabel(label),
00078    fIndex(kMinInt),
00079    fStatus(t->GetStatusCode()),
00080    fLockPoints(kFALSE),
00081    fPathMarks(),
00082    fLastPMIdx(0),
00083    fPropagator(0)
00084 {
00085    // Constructor from TParticle.
00086 
00087    SetPropagator(prop);
00088    fMainColorPtr = &fLineColor;
00089 
00090    TParticlePDG* pdgp = t->GetPDG();
00091    if (pdgp) {
00092       fPdg    = pdgp->PdgCode();
00093       fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3);
00094    }
00095 
00096    SetName(t->GetName());
00097 }
00098 
00099 //______________________________________________________________________________
00100 TEveTrack::TEveTrack(TEveMCTrack* t, TEveTrackPropagator* prop):
00101    TEveLine(),
00102 
00103    fV(t->Vx(), t->Vy(), t->Vz()),
00104    fP(t->Px(), t->Py(), t->Pz()),
00105    fPEnd(),
00106    fBeta(t->P()/t->Energy()),
00107    fPdg(0),
00108    fCharge(0),
00109    fLabel(t->fLabel),
00110    fIndex(t->fIndex),
00111    fStatus(t->GetStatusCode()),
00112    fLockPoints(kFALSE),
00113    fPathMarks(),
00114    fLastPMIdx(0),
00115    fPropagator(0)
00116 {
00117    // Constructor from TEveUtil Monte Carlo track.
00118 
00119    SetPropagator(prop);
00120    fMainColorPtr = &fLineColor;
00121 
00122    TParticlePDG* pdgp = t->GetPDG();
00123    if (pdgp == 0) {
00124       t->ResetPdgCode(); pdgp = t->GetPDG();
00125    }
00126    fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3);
00127 
00128    SetName(t->GetName());
00129 }
00130 
00131 //______________________________________________________________________________
00132 TEveTrack::TEveTrack(TEveRecTrack* t, TEveTrackPropagator* prop) :
00133    TEveLine(),
00134 
00135    fV(t->fV),
00136    fP(t->fP),
00137    fPEnd(),
00138    fBeta(t->fBeta),
00139    fPdg(0),
00140    fCharge(t->fSign),
00141    fLabel(t->fLabel),
00142    fIndex(t->fIndex),
00143    fStatus(t->fStatus),
00144    fLockPoints(kFALSE),
00145    fPathMarks(),
00146    fLastPMIdx(0),
00147    fPropagator(0)
00148 {
00149    // Constructor from TEveUtil reconstructed track.
00150 
00151    SetPropagator(prop);
00152    fMainColorPtr = &fLineColor;
00153 
00154    SetName(t->GetName());
00155 }
00156 
00157 //______________________________________________________________________________
00158 TEveTrack::TEveTrack(const TEveTrack& t) :
00159    TEveLine(),
00160    fV(t.fV),
00161    fP(t.fP),
00162    fPEnd(),
00163    fBeta(t.fBeta),
00164    fPdg(t.fPdg),
00165    fCharge(t.fCharge),
00166    fLabel(t.fLabel),
00167    fIndex(t.fIndex),
00168    fStatus(t.fStatus),
00169    fLockPoints(t.fLockPoints),
00170    fPathMarks(),
00171    fLastPMIdx(t.fLastPMIdx),
00172    fPropagator(0)
00173 {
00174    // Copy constructor. Track paremeters are copied but the
00175    // extrapolation is not perfermed so you should still call
00176    // MakeTrack() to do that.
00177    // If points of 't' are locked, they are cloned.
00178 
00179    if (fLockPoints)
00180       ClonePoints(t);
00181 
00182    SetPathMarks(t);
00183    SetPropagator (t.fPropagator);
00184 
00185    CopyVizParams(&t);
00186 }
00187 
00188 //______________________________________________________________________________
00189 TEveTrack::~TEveTrack()
00190 {
00191    // Destructor.
00192 
00193    SetPropagator(0);
00194 }
00195 
00196 
00197 //______________________________________________________________________________
00198 const TGPicture* TEveTrack::GetListTreeIcon(Bool_t)
00199 {
00200    // Returns list-tree icon for TEveTrack.
00201 
00202    return fgListTreeIcons[4];
00203 }
00204 
00205 //______________________________________________________________________________
00206 void TEveTrack::ComputeBBox()
00207 {
00208    // Compute the bounding box of the track.
00209 
00210    if (Size() > 0 || ! fPathMarks.empty())
00211    {
00212       BBoxInit();
00213       Int_t    n = Size();
00214       Float_t* p = TPolyMarker3D::fP;
00215       for (Int_t i = 0; i < n; ++i, p += 3)
00216       {
00217          BBoxCheckPoint(p);
00218       }
00219       for (vPathMark_ci i = fPathMarks.begin(); i != fPathMarks.end(); ++i)
00220       {
00221          BBoxCheckPoint(i->fV);
00222       }
00223    }
00224    else
00225    {
00226       BBoxZero();
00227    }
00228 }
00229 
00230 //==============================================================================
00231 
00232 //______________________________________________________________________________
00233 void TEveTrack::SetStdTitle()
00234 {
00235    // Set standard track title based on most data-member values.
00236 
00237    TString idx(fIndex == kMinInt ? "<undef>" : Form("%d", fIndex));
00238    TString lbl(fLabel == kMinInt ? "<undef>" : Form("%d", fLabel));
00239    SetTitle(Form("Index=%s, Label=%s\nChg=%d, Pdg=%d\n"
00240                  "pT=%.3f, pZ=%.3f\nV=(%.3f, %.3f, %.3f)",
00241                  idx.Data(), lbl.Data(), fCharge, fPdg,
00242                  fP.Perp(), fP.fZ, fV.fX, fV.fY, fV.fZ));
00243 }
00244 
00245 //______________________________________________________________________________
00246 void TEveTrack::SetTrackParams(const TEveTrack& t)
00247 {
00248    // Copy track parameters from t. Track-propagator is set, too.
00249    // PathMarks are cleared - you can copy them via SetPathMarks(t).
00250    // If track 't' is locked, you should probably clone its points
00251    // over - use TEvePointSet::ClonePoints(t);
00252 
00253    fV          = t.fV;
00254    fP          = t.fP;
00255    fBeta       = t.fBeta;
00256    fPdg        = t.fPdg;
00257    fCharge     = t.fCharge;
00258    fLabel      = t.fLabel;
00259    fIndex      = t.fIndex;
00260 
00261    fPathMarks.clear();
00262    SetPropagator(t.fPropagator);
00263 }
00264 
00265 //______________________________________________________________________________
00266 void TEveTrack::SetPathMarks(const TEveTrack& t)
00267 {
00268    // Copy path-marks from t.
00269 
00270    std::copy(t.RefPathMarks().begin(), t.RefPathMarks().end(),
00271              std::back_insert_iterator<vPathMark_t>(fPathMarks));
00272 }
00273 
00274 //==============================================================================
00275 
00276 //______________________________________________________________________________
00277 void TEveTrack::SetPropagator(TEveTrackPropagator* prop)
00278 {
00279    // Set track's render style.
00280    // Reference counts of old and new propagator are updated.
00281 
00282    if (fPropagator == prop) return;
00283    if (fPropagator) fPropagator->DecRefCount(this);
00284    fPropagator = prop;
00285    if (fPropagator) prop->IncRefCount(this);
00286 }
00287 
00288 //==============================================================================
00289 
00290 //______________________________________________________________________________
00291 void TEveTrack::SetAttLineAttMarker(TEveTrackList* tl)
00292 {
00293    // Set line and marker attributes from TEveTrackList.
00294 
00295    SetRnrLine(tl->GetRnrLine());
00296    SetLineColor(tl->GetLineColor());
00297    SetLineStyle(tl->GetLineStyle());
00298    SetLineWidth(tl->GetLineWidth());
00299 
00300    SetRnrPoints(tl->GetRnrPoints());
00301    SetMarkerColor(tl->GetMarkerColor());
00302    SetMarkerStyle(tl->GetMarkerStyle());
00303    SetMarkerSize(tl->GetMarkerSize());
00304 }
00305 
00306 //==============================================================================
00307 
00308 //______________________________________________________________________________
00309 void TEveTrack::MakeTrack(Bool_t recurse)
00310 {
00311    // Calculate track representation based on track data and current
00312    // settings of the propagator.
00313    // If recurse is true, descend into children.
00314 
00315    if (!fLockPoints)
00316    {
00317       Reset(0);
00318       fLastPMIdx = 0;
00319 
00320       TEveTrackPropagator& rTP((fPropagator != 0) ? *fPropagator : TEveTrackPropagator::fgDefault);
00321 
00322       const Float_t maxRsq = rTP.GetMaxR() * rTP.GetMaxR();
00323       const Float_t maxZ   = rTP.GetMaxZ();
00324 
00325       if ( ! TEveTrackPropagator::IsOutsideBounds(fV, maxRsq, maxZ))
00326       {
00327          TEveVector currP = fP;
00328          Bool_t decay = kFALSE;
00329          fPropagator->InitTrack(fV, fCharge);
00330          for (vPathMark_i pm = fPathMarks.begin(); pm != fPathMarks.end(); ++pm, ++fLastPMIdx)
00331          {
00332             if (rTP.GetFitReferences() && pm->fType == TEvePathMark::kReference)
00333             {
00334                if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
00335                   break;
00336                // printf("%s fit reference  \n", fName.Data());
00337                if (fPropagator->GoToVertex(pm->fV, currP)) {
00338                   currP.fX = pm->fP.fX; currP.fY = pm->fP.fY; currP.fZ = pm->fP.fZ;
00339                }
00340                else
00341                {
00342                   break;
00343                }
00344             }
00345             else if (rTP.GetFitDaughters() && pm->fType == TEvePathMark::kDaughter)
00346             {
00347                if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
00348                   break;
00349                // printf("%s fit daughter  \n", fName.Data());
00350                if (fPropagator->GoToVertex(pm->fV, currP)) {
00351                   currP.fX -= pm->fP.fX; currP.fY -= pm->fP.fY; currP.fZ -= pm->fP.fZ;
00352                }
00353                else
00354                {
00355                   break;
00356                }
00357             }
00358             else if (rTP.GetFitDecay() && pm->fType == TEvePathMark::kDecay)
00359             {
00360                if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
00361                   break;
00362                // printf("%s fit decay \n", fName.Data());
00363                fPropagator->GoToVertex(pm->fV, currP);
00364                decay = kTRUE;
00365                ++fLastPMIdx;
00366                break;
00367             }
00368             else if (rTP.GetFitCluster2Ds() && pm->fType == TEvePathMark::kCluster2D)
00369             {
00370                TEveVector itsect;
00371                if (fPropagator->IntersectPlane(currP, pm->fV, pm->fP, itsect))
00372                {
00373                   TEveVector delta   = itsect - pm->fV;
00374                   TEveVector vtopass = pm->fV + pm->fE*(pm->fE.Dot(delta));
00375                   if (TEveTrackPropagator::IsOutsideBounds(vtopass, maxRsq, maxZ))
00376                      break;
00377                   if ( ! fPropagator->GoToVertex(vtopass, currP))
00378                      break;
00379                }
00380                else
00381                {
00382                   Warning("TEveTrack::MakeTrack", "Failed to intersect plane for Cluster2D. Ignoring path-mark.");
00383                }
00384             }
00385             else
00386             {
00387                if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
00388                   break;
00389             }            
00390          } // loop path-marks
00391 
00392          if (!decay)
00393          {
00394             // printf("%s loop to bounds  \n",fName.Data() );
00395             fPropagator->GoToBounds(currP);
00396          }
00397          fPEnd = currP;
00398          //  make_polyline:
00399          fPropagator->FillPointSet(this);
00400          fPropagator->ResetTrack();
00401       }
00402    }
00403 
00404    if (recurse)
00405    {
00406       for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
00407       {
00408          TEveTrack* t = dynamic_cast<TEveTrack*>(*i);
00409          if (t) t->MakeTrack(recurse);
00410       }
00411    }
00412 }
00413 
00414 //==============================================================================
00415 
00416 //______________________________________________________________________________
00417 void TEveTrack::CopyVizParams(const TEveElement* el)
00418 {
00419    // Copy visualization parameters from element el.
00420 
00421    // No local parameters.
00422 
00423    // const TEveTrack* t = dynamic_cast<const TEveTrack*>(el);
00424    // if (t)
00425    // {}
00426 
00427    TEveLine::CopyVizParams(el);
00428 }
00429 
00430 //______________________________________________________________________________
00431 void TEveTrack::WriteVizParams(ostream& out, const TString& var)
00432 {
00433    // Write visualization parameters.
00434 
00435    TEveLine::WriteVizParams(out, var);
00436 
00437    // TString t = "   " + var + "->";
00438 }
00439 
00440 //______________________________________________________________________________
00441 TClass* TEveTrack::ProjectedClass(const TEveProjection*) const
00442 {
00443    // Virtual from TEveProjectable, return TEveTrackProjected class.
00444 
00445    return TEveTrackProjected::Class();
00446 }
00447 
00448 //==============================================================================
00449 
00450 namespace
00451 {
00452    struct Cmp_pathmark_t
00453    {
00454       bool operator()(TEvePathMark const & a, TEvePathMark const & b)
00455       { return a.fTime < b.fTime; }
00456    };
00457 }
00458 
00459 //______________________________________________________________________________
00460 void TEveTrack::SortPathMarksByTime()
00461 {
00462    // Sort registerd pat-marks by time.
00463 
00464    std::sort(fPathMarks.begin(), fPathMarks.end(), Cmp_pathmark_t());
00465 }
00466 
00467 //______________________________________________________________________________
00468 void TEveTrack::PrintPathMarks()
00469 {
00470    // Print registered path-marks.
00471 
00472    static const TEveException eh("TEveTrack::PrintPathMarks ");
00473 
00474    printf("TEveTrack '%s', number of path marks %d, label %d\n",
00475           GetName(), (Int_t)fPathMarks.size(), fLabel);
00476 
00477    for (vPathMark_i pm = fPathMarks.begin(); pm != fPathMarks.end(); ++pm)
00478    {
00479       printf("  %-9s  p: %8f %8f %8f Vertex: %8e %8e %8e %g \n",
00480              pm->TypeName(),
00481              pm->fP.fX,  pm->fP.fY, pm->fP.fZ,
00482              pm->fV.fX,  pm->fV.fY, pm->fV.fZ,
00483              pm->fTime);
00484    }
00485 }
00486 
00487 //------------------------------------------------------------------------------
00488 
00489 //______________________________________________________________________________
00490 void TEveTrack::SecSelected(TEveTrack* track)
00491 {
00492    // Emits "SecSelected(TEveTrack*)" signal.
00493    // Called from TEveTrackGL on secondary-selection.
00494 
00495    Emit("SecSelected(TEveTrack*)", (Long_t)track);
00496 }
00497 
00498 //------------------------------------------------------------------------------
00499 
00500 //______________________________________________________________________________
00501 Bool_t TEveTrack::ShouldBreakTrack() const
00502 {
00503    // Should this track be broken in projections.
00504 
00505    Error("ShouldBreakTrack", "Deprected -- use TEveTrackPropagator functions.");
00506    return fPropagator->GetProjTrackBreaking() == TEveTrackPropagator::kPTB_Break;
00507 }
00508 
00509 //______________________________________________________________________________
00510 UChar_t TEveTrack::GetBreakProjectedTracks() const
00511 {
00512    // Deprected -- use TEveTrackPropagator functions.
00513    Error("GetBreakProjectedTracks", "Deprected -- use TEveTrackPropagator functions.");
00514    return 0;
00515 }
00516 
00517 //______________________________________________________________________________
00518 void TEveTrack::SetBreakProjectedTracks(UChar_t)
00519 {
00520    // Deprected -- use TEveTrackPropagator functions.
00521 
00522    Error("SetBreakProjectedTracks", "Deprected -- use TEveTrackPropagator functions.");
00523 }
00524 
00525 //______________________________________________________________________________
00526 Bool_t TEveTrack::GetDefaultBreakProjectedTracks()
00527 {
00528    // Deprected -- use TEveTrackPropagator functions.
00529    // Return true if tracks get broken into several segments when the
00530    // projected space consists of separate domains (like Rho-Z).
00531    // Static function.
00532 
00533    ::Error("TEveTrack::GetDefaultBreakProjectedTracks", "Deprected -- use TEveTrackPropagator functions.");
00534    return kTRUE;
00535 }
00536 
00537 //______________________________________________________________________________
00538 void TEveTrack::SetDefaultBreakProjectedTracks(Bool_t)
00539 {
00540    // Deprected -- use TEveTrackPropagator functions.
00541 
00542    ::Error("TEveTrack::SetDefaultBreakProjectedTracks", "Deprected -- use TEveTrackPropagator functions.");
00543 }
00544 
00545 
00546 //==============================================================================
00547 //==============================================================================
00548 // TEveTrackList
00549 //==============================================================================
00550 
00551 //______________________________________________________________________________
00552 //
00553 // A list of tracks supporting change of common attributes and
00554 // selection based on track parameters.
00555 
00556 ClassImp(TEveTrackList);
00557 
00558 //______________________________________________________________________________
00559 TEveTrackList::TEveTrackList(TEveTrackPropagator* prop) :
00560    TEveElementList(),
00561    TAttMarker(1, 20, 1),
00562    TAttLine(1,1,1),
00563 
00564    fPropagator(0),
00565    fRecurse(kTRUE),
00566    fRnrLine(kTRUE),
00567    fRnrPoints(kFALSE),
00568 
00569    fMinPt (0), fMaxPt (0), fLimPt (0),
00570    fMinP  (0), fMaxP  (0), fLimP  (0)
00571 {
00572    // Constructor. If track-propagator argument is 0, a new default
00573    // one is created.
00574 
00575    fChildClass = TEveTrack::Class(); // override member from base TEveElementList
00576 
00577    fMainColorPtr = &fLineColor;
00578 
00579    if (prop == 0) prop = new TEveTrackPropagator;
00580    SetPropagator(prop);
00581 }
00582 
00583 //______________________________________________________________________________
00584 TEveTrackList::TEveTrackList(const char* name, TEveTrackPropagator* prop) :
00585    TEveElementList(name),
00586    TAttMarker(1, 20, 1),
00587    TAttLine(1,1,1),
00588 
00589    fPropagator(0),
00590    fRecurse(kTRUE),
00591    fRnrLine(kTRUE),
00592    fRnrPoints(kFALSE),
00593 
00594    fMinPt (0), fMaxPt (0), fLimPt (0),
00595    fMinP  (0), fMaxP  (0), fLimP  (0)
00596 {
00597    // Constructor. If track-propagator argument is 0, a new default
00598    // one is created.
00599 
00600    fChildClass = TEveTrack::Class(); // override member from base TEveElementList
00601 
00602    fMainColorPtr = &fLineColor;
00603 
00604    if (prop == 0) prop = new TEveTrackPropagator;
00605    SetPropagator(prop);
00606 }
00607 
00608 //______________________________________________________________________________
00609 TEveTrackList::~TEveTrackList()
00610 {
00611    // Destructor.
00612 
00613    SetPropagator(0);
00614 }
00615 
00616 //==============================================================================
00617 
00618 //______________________________________________________________________________
00619 void TEveTrackList::SetPropagator(TEveTrackPropagator* prop)
00620 {
00621    // Set default propagator for tracks.
00622    // This is not enforced onto the tracks themselves but this is the
00623    // propagator that is shown in the TEveTrackListEditor.
00624 
00625    if (fPropagator == prop) return;
00626    if (fPropagator) fPropagator->DecRefCount();
00627    fPropagator = prop;
00628    if (fPropagator) prop->IncRefCount();
00629 }
00630 
00631 //==============================================================================
00632 
00633 //______________________________________________________________________________
00634 void TEveTrackList::MakeTracks(Bool_t recurse)
00635 {
00636    // Regenerate the visual representations of tracks.
00637    // The momentum limits are rescanned during the same traversal.
00638 
00639    fLimPt = fLimP = 0;
00640 
00641    for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
00642    {
00643       TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
00644       if (track)
00645       {
00646          track->MakeTrack(recurse);
00647 
00648          fLimPt = TMath::Max(fLimPt, track->fP.Perp());
00649          fLimP  = TMath::Max(fLimP,  track->fP.Mag());
00650       }
00651       if (recurse)
00652          FindMomentumLimits(*i, recurse);
00653    }
00654 
00655    fLimPt = RoundMomentumLimit(fLimPt);
00656    fLimP  = RoundMomentumLimit(fLimP);
00657 
00658    SanitizeMinMaxCuts();
00659 }
00660 
00661 //______________________________________________________________________________
00662 void TEveTrackList::FindMomentumLimits(Bool_t recurse)
00663 {
00664    // Loop over children and find highest pT and p of contained TEveTracks.
00665    // These are stored in members fLimPt and fLimP.
00666 
00667    fLimPt = fLimP = 0;
00668 
00669    if (HasChildren())
00670    {
00671       for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00672       {
00673          TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
00674          if (track)
00675          {
00676             fLimPt = TMath::Max(fLimPt, track->fP.Perp());
00677             fLimP  = TMath::Max(fLimP,  track->fP.Mag());
00678          }
00679          if (recurse)
00680             FindMomentumLimits(*i, recurse);
00681       }
00682 
00683       fLimPt = RoundMomentumLimit(fLimPt);
00684       fLimP  = RoundMomentumLimit(fLimP);
00685    }
00686 
00687    SanitizeMinMaxCuts();
00688 }
00689 
00690 //______________________________________________________________________________
00691 void TEveTrackList::FindMomentumLimits(TEveElement* el, Bool_t recurse)
00692 {
00693    // Loop over track elements of argument el and find highest pT and p.
00694    // These are stored in members fLimPt and fLimP.
00695 
00696    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00697    {
00698       TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
00699       if (track)
00700       {
00701          fLimPt = TMath::Max(fLimPt, track->fP.Perp());
00702          fLimP  = TMath::Max(fLimP,  track->fP.Mag());
00703       }
00704       if (recurse)
00705          FindMomentumLimits(*i, recurse);
00706    }
00707 }
00708 
00709 //______________________________________________________________________________
00710 Float_t TEveTrackList::RoundMomentumLimit(Float_t x)
00711 {
00712    // Round the momentum limit up to a nice value.
00713 
00714    using namespace TMath;
00715 
00716    if (x < 1e-3) return 1e-3;
00717 
00718    Double_t fac = Power(10, 1 - Floor(Log10(x)));
00719    return Ceil(fac*x) / fac;
00720 }
00721 
00722 //______________________________________________________________________________
00723 void TEveTrackList::SanitizeMinMaxCuts()
00724 {
00725    // Set Min/Max cuts so that they are within detected limits.
00726 
00727    using namespace TMath;
00728 
00729    fMinPt = Min(fMinPt, fLimPt);
00730    fMaxPt = fMaxPt == 0 ? fLimPt : Min(fMaxPt, fLimPt);
00731    fMinP  = Min(fMinP,  fLimP);
00732    fMaxP  = fMaxP  == 0 ? fLimP  : Min(fMaxP,  fLimP);
00733 }
00734 
00735 //==============================================================================
00736 
00737 //______________________________________________________________________________
00738 void TEveTrackList::SetRnrLine(Bool_t rnr)
00739 {
00740    // Set rendering of track as line for the list and the elements.
00741 
00742    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00743    {
00744       TEveTrack* track = (TEveTrack*)(*i);
00745       if (track->GetRnrLine() == fRnrLine)
00746          track->SetRnrLine(rnr);
00747       if (fRecurse)
00748          SetRnrLine(rnr, *i);
00749    }
00750    fRnrLine = rnr;
00751 }
00752 
00753 //______________________________________________________________________________
00754 void TEveTrackList::SetRnrLine(Bool_t rnr, TEveElement* el)
00755 {
00756    // Set rendering of track as line for children of el.
00757 
00758    TEveTrack* track;
00759    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00760    {
00761       track = dynamic_cast<TEveTrack*>(*i);
00762       if (track && (track->GetRnrLine() == fRnrLine))
00763          track->SetRnrLine(rnr);
00764       if (fRecurse)
00765          SetRnrLine(rnr, *i);
00766    }
00767 }
00768 
00769 //==============================================================================
00770 
00771 //______________________________________________________________________________
00772 void TEveTrackList::SetRnrPoints(Bool_t rnr)
00773 {
00774    // Set rendering of track as points for the list and the elements.
00775 
00776    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00777    {
00778       TEveTrack* track = (TEveTrack*)(*i);
00779       if (track->GetRnrPoints() == fRnrPoints)
00780          track->SetRnrPoints(rnr);
00781       if (fRecurse)
00782          SetRnrPoints(rnr, *i);
00783    }
00784    fRnrPoints = rnr;
00785 }
00786 
00787 //______________________________________________________________________________
00788 void TEveTrackList::SetRnrPoints(Bool_t rnr, TEveElement* el)
00789 {
00790    // Set rendering of track as points for children of el.
00791 
00792    TEveTrack* track;
00793    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00794    {
00795       track = dynamic_cast<TEveTrack*>(*i);
00796       if (track)
00797          if (track->GetRnrPoints() == fRnrPoints)
00798             track->SetRnrPoints(rnr);
00799       if (fRecurse)
00800          SetRnrPoints(rnr, *i);
00801    }
00802 }
00803 
00804 //==============================================================================
00805 
00806 //______________________________________________________________________________
00807 void TEveTrackList::SetMainColor(Color_t col)
00808 {
00809    // Set main (line) color for the list and the elements.
00810 
00811    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00812    {
00813       TEveTrack* track = (TEveTrack*)(*i);
00814       if (track->GetLineColor() == fLineColor)
00815          track->SetLineColor(col);
00816       if (fRecurse)
00817          SetLineColor(col, *i);
00818    }
00819    TEveElement::SetMainColor(col);
00820 }
00821 
00822 //______________________________________________________________________________
00823 void TEveTrackList::SetLineColor(Color_t col, TEveElement* el)
00824 {
00825    // Set line color for children of el.
00826 
00827    TEveTrack* track;
00828    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00829    {
00830       track = dynamic_cast<TEveTrack*>(*i);
00831       if (track && track->GetLineColor() == fLineColor)
00832          track->SetLineColor(col);
00833       if (fRecurse)
00834          SetLineColor(col, *i);
00835    }
00836 }
00837 
00838 //==============================================================================
00839 
00840 //______________________________________________________________________________
00841 void TEveTrackList::SetLineWidth(Width_t width)
00842 {
00843    // Set line width for the list and the elements.
00844 
00845    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00846    {
00847       TEveTrack* track = (TEveTrack*)(*i);
00848       if (track->GetLineWidth() == fLineWidth)
00849          track->SetLineWidth(width);
00850       if (fRecurse)
00851          SetLineWidth(width, *i);
00852    }
00853    fLineWidth = width;
00854 }
00855 
00856 //______________________________________________________________________________
00857 void TEveTrackList::SetLineWidth(Width_t width, TEveElement* el)
00858 {
00859    // Set line width for children of el.
00860 
00861    TEveTrack* track;
00862    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00863    {
00864       track = dynamic_cast<TEveTrack*>(*i);
00865       if (track && track->GetLineWidth() == fLineWidth)
00866          track->SetLineWidth(width);
00867       if (fRecurse)
00868          SetLineWidth(width, *i);
00869    }
00870 }
00871 
00872 //==============================================================================
00873 
00874 //______________________________________________________________________________
00875 void TEveTrackList::SetLineStyle(Style_t style)
00876 {
00877    // Set line style for the list and the elements.
00878 
00879    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00880    {
00881       TEveTrack* track = (TEveTrack*)(*i);
00882       if (track->GetLineStyle() == fLineStyle)
00883          track->SetLineStyle(style);
00884       if (fRecurse)
00885          SetLineStyle(style, *i);
00886    }
00887    fLineStyle = style;
00888 }
00889 
00890 //______________________________________________________________________________
00891 void TEveTrackList::SetLineStyle(Style_t style, TEveElement* el)
00892 {
00893    // Set line style for children of el.
00894 
00895    TEveTrack* track;
00896    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00897    {
00898       track = dynamic_cast<TEveTrack*>(*i);
00899       if (track && track->GetLineStyle() == fLineStyle)
00900          track->SetLineStyle(style);
00901       if (fRecurse)
00902          SetLineStyle(style, *i);
00903    }
00904 }
00905 
00906 //==============================================================================
00907 
00908 //______________________________________________________________________________
00909 void TEveTrackList::SetMarkerStyle(Style_t style)
00910 {
00911    // Set marker style for the list and the elements.
00912 
00913    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00914    {
00915       TEveTrack* track = (TEveTrack*)(*i);
00916       if (track->GetMarkerStyle() == fMarkerStyle)
00917          track->SetMarkerStyle(style);
00918       if (fRecurse)
00919          SetMarkerStyle(style, *i);
00920    }
00921    fMarkerStyle = style;
00922 }
00923 
00924 //______________________________________________________________________________
00925 void TEveTrackList::SetMarkerStyle(Style_t style, TEveElement* el)
00926 {
00927    // Set marker style for children of el.
00928 
00929    TEveTrack* track;
00930    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00931    {
00932       track = dynamic_cast<TEveTrack*>(*i);
00933       if (track && track->GetMarkerStyle() == fMarkerStyle)
00934          track->SetMarkerStyle(style);
00935       if (fRecurse)
00936          SetMarkerStyle(style, *i);
00937    }
00938 }
00939 
00940 //==============================================================================
00941 
00942 //______________________________________________________________________________
00943 void TEveTrackList::SetMarkerColor(Color_t col)
00944 {
00945    // Set marker color for the list and the elements.
00946 
00947    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00948    {
00949       TEveTrack* track = (TEveTrack*)(*i);
00950       if (track->GetMarkerColor() == fMarkerColor)
00951          track->SetMarkerColor(col);
00952       if (fRecurse)
00953          SetMarkerColor(col, *i);
00954    }
00955    fMarkerColor = col;
00956 }
00957 
00958 //______________________________________________________________________________
00959 void TEveTrackList::SetMarkerColor(Color_t col, TEveElement* el)
00960 {
00961    // Set marker color for children of el.
00962 
00963    TEveTrack* track;
00964    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00965    {
00966       track = dynamic_cast<TEveTrack*>(*i);
00967       if (track && track->GetMarkerColor() == fMarkerColor)
00968          track->SetMarkerColor(col);
00969       if (fRecurse)
00970          SetMarkerColor(col, *i);
00971    }
00972 }
00973 
00974 //==============================================================================
00975 
00976 //______________________________________________________________________________
00977 void TEveTrackList::SetMarkerSize(Size_t size)
00978 {
00979    // Set marker size for the list and the elements.
00980 
00981    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
00982    {
00983       TEveTrack* track = (TEveTrack*)(*i);
00984       if (track->GetMarkerSize() == fMarkerSize)
00985          track->SetMarkerSize(size);
00986       if (fRecurse)
00987          SetMarkerSize(size, *i);
00988    }
00989    fMarkerSize = size;
00990 }
00991 
00992 //______________________________________________________________________________
00993 void TEveTrackList::SetMarkerSize(Size_t size, TEveElement* el)
00994 {
00995    // Set marker size for children of el.
00996 
00997    TEveTrack* track;
00998    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
00999    {
01000       track = dynamic_cast<TEveTrack*>(*i);
01001       if (track && track->GetMarkerSize() == fMarkerSize)
01002          track->SetMarkerSize(size);
01003       if (fRecurse)
01004          SetMarkerSize(size, *i);
01005    }
01006 }
01007 
01008 //==============================================================================
01009 
01010 //______________________________________________________________________________
01011 void TEveTrackList::SelectByPt(Float_t min_pt, Float_t max_pt)
01012 {
01013    // Select visibility of tracks by transverse momentum.
01014    // If data-member fRecurse is set, the selection is applied
01015    // recursively to all children.
01016 
01017    fMinPt = min_pt;
01018    fMaxPt = max_pt;
01019 
01020    const Float_t minptsq = min_pt*min_pt;
01021    const Float_t maxptsq = max_pt*max_pt;
01022 
01023    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
01024    {
01025       const Float_t ptsq = ((TEveTrack*)(*i))->fP.Perp2();
01026       Bool_t on = ptsq >= minptsq && ptsq <= maxptsq;
01027       (*i)->SetRnrState(on);
01028       if (on && fRecurse)
01029          SelectByPt(min_pt, max_pt, *i);
01030    }
01031 }
01032 
01033 //______________________________________________________________________________
01034 void TEveTrackList::SelectByPt(Float_t min_pt, Float_t max_pt, TEveElement* el)
01035 {
01036    // Select visibility of el's children tracks by transverse momentum.
01037 
01038    const Float_t minptsq = min_pt*min_pt;
01039    const Float_t maxptsq = max_pt*max_pt;
01040 
01041    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
01042    {
01043       TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
01044       if (track)
01045       {
01046          const Float_t ptsq = track->fP.Perp2();
01047          Bool_t on = ptsq >= minptsq && ptsq <= maxptsq;
01048          track->SetRnrState(on);
01049          if (on && fRecurse)
01050             SelectByPt(min_pt, max_pt, *i);
01051       }
01052    }
01053 }
01054 
01055 //______________________________________________________________________________
01056 void TEveTrackList::SelectByP(Float_t min_p, Float_t max_p)
01057 {
01058    // Select visibility of tracks by momentum.
01059    // If data-member fRecurse is set, the selection is applied
01060    // recursively to all children.
01061 
01062    fMinP = min_p;
01063    fMaxP = max_p;
01064 
01065    const Float_t minpsq = min_p*min_p;
01066    const Float_t maxpsq = max_p*max_p;
01067 
01068    for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
01069    {
01070       const Float_t psq  = ((TEveTrack*)(*i))->fP.Mag2();
01071       Bool_t on = psq >= minpsq && psq <= maxpsq;
01072       (*i)->SetRnrState(psq >= minpsq && psq <= maxpsq);
01073       if (on && fRecurse)
01074          SelectByP(min_p, max_p, *i);
01075    }
01076 }
01077 
01078 //______________________________________________________________________________
01079 void TEveTrackList::SelectByP(Float_t min_p, Float_t max_p, TEveElement* el)
01080 {
01081    // Select visibility of el's children tracks by momentum.
01082 
01083    const Float_t minpsq = min_p*min_p;
01084    const Float_t maxpsq = max_p*max_p;
01085 
01086    for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
01087    {
01088       TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
01089       if (track)
01090       {
01091          const Float_t psq  = ((TEveTrack*)(*i))->fP.Mag2();
01092          Bool_t on = psq >= minpsq && psq <= maxpsq;
01093          track->SetRnrState(on);
01094          if (on && fRecurse)
01095             SelectByP(min_p, max_p, *i);
01096       }
01097    }
01098 }
01099 
01100 //==============================================================================
01101 
01102 //______________________________________________________________________________
01103 TEveTrack* TEveTrackList::FindTrackByLabel(Int_t label)
01104 {
01105    // Find track by label, select it and display it in the editor.
01106 
01107    for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
01108    {
01109       if (((TEveTrack*)(*i))->GetLabel() == label)
01110       {
01111          TGListTree     *lt   = gEve->GetLTEFrame()->GetListTree();
01112          TGListTreeItem *mlti = lt->GetSelected();
01113          if (mlti->GetUserData() != this)
01114             mlti = FindListTreeItem(lt);
01115          TGListTreeItem *tlti = (*i)->FindListTreeItem(lt, mlti);
01116          lt->HighlightItem(tlti);
01117          lt->SetSelected(tlti);
01118          gEve->EditElement(*i);
01119          return (TEveTrack*) *i;
01120       }
01121    }
01122    return 0;
01123 }
01124 
01125 //______________________________________________________________________________
01126 TEveTrack* TEveTrackList::FindTrackByIndex(Int_t index)
01127 {
01128    // Find track by index, select it and display it in the editor.
01129 
01130    for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
01131    {
01132       if (((TEveTrack*)(*i))->GetIndex() == index)
01133       {
01134          TGListTree     *lt   = gEve->GetLTEFrame()->GetListTree();
01135          TGListTreeItem *mlti = lt->GetSelected();
01136          if (mlti->GetUserData() != this)
01137             mlti = FindListTreeItem(lt);
01138          TGListTreeItem *tlti = (*i)->FindListTreeItem(lt, mlti);
01139          lt->HighlightItem(tlti);
01140          lt->SetSelected(tlti);
01141          gEve->EditElement(*i);
01142          return (TEveTrack*) *i;
01143       }
01144    }
01145    return 0;
01146 }
01147 
01148 //==============================================================================
01149 
01150 //______________________________________________________________________________
01151 void TEveTrackList::CopyVizParams(const TEveElement* el)
01152 {
01153    // Copy visualization parameters from element el.
01154 
01155    const TEveTrackList* m = dynamic_cast<const TEveTrackList*>(el);
01156    if (m)
01157    {
01158       TAttMarker::operator=(*m);
01159       TAttLine::operator=(*m);
01160       fRecurse = m->fRecurse;
01161       fRnrLine = m->fRnrLine;
01162       fRnrPoints = m->fRnrPoints;
01163       fMinPt   = m->fMinPt;
01164       fMaxPt   = m->fMaxPt;
01165       fLimPt   = m->fLimPt;
01166       fMinP    = m->fMinP;
01167       fMaxP    = m->fMaxP;
01168       fLimP    = m->fLimP;
01169    }
01170 
01171    TEveElement::CopyVizParams(el);
01172 }
01173 
01174 //______________________________________________________________________________
01175 void TEveTrackList::WriteVizParams(ostream& out, const TString& var)
01176 {
01177    // Write visualization parameters.
01178 
01179    TEveElement::WriteVizParams(out, var);
01180 
01181    TString t = "   " + var + "->";
01182    TAttMarker::SaveMarkerAttributes(out, var);
01183    TAttLine  ::SaveLineAttributes  (out, var);
01184    out << t << "SetRecurse("   << ToString(fRecurse)   << ");\n";
01185    out << t << "SetRnrLine("   << ToString(fRnrLine)   << ");\n";
01186    out << t << "SetRnrPoints(" << ToString(fRnrPoints) << ");\n";
01187    // These setters are not available -- need proper AND/OR mode.
01188    // out << t << "SetMinPt(" << fMinPt << ");\n";
01189    // out << t << "SetMaxPt(" << fMaxPt << ");\n";
01190    // out << t << "SetLimPt(" << fLimPt << ");\n";
01191    // out << t << "SetMinP("  << fMinP  << ");\n";
01192    // out << t << "SetMaxP("  << fMaxP  << ");\n";
01193    // out << t << "SetLimP("  << fLimP  << ");\n";
01194 }
01195 
01196 //______________________________________________________________________________
01197 TClass* TEveTrackList::ProjectedClass(const TEveProjection*) const
01198 {
01199    // Virtual from TEveProjectable, returns TEveTrackListProjected class.
01200 
01201    return TEveTrackListProjected::Class();
01202 }

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