BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page

BrTpcTrackFollowModule


class description - source file - inheritance tree

class BrTpcTrackFollowModule : public BrTpcTrackModule

    private:
virtual void FindHitsNextRow(const Int_t rowNo, const Float_t* proj, const Float_t dx1, const Float_t dx2, const Float_t dy1, const Float_t dy2, TObjArray& FoundHits) virtual void FindTracks() virtual void FindTracksFromRow(const Int_t istart, const Int_t ilast) void TrackSingleHit(const Int_t istart, const Int_t ilast) public:
BrTpcTrackFollowModule BrTpcTrackFollowModule() BrTpcTrackFollowModule BrTpcTrackFollowModule(const Char_t* name, const Char_t* title) BrTpcTrackFollowModule BrTpcTrackFollowModule(BrTpcTrackFollowModule&) virtual void ~BrTpcTrackFollowModule() static TClass* Class() virtual void DefineHistograms() virtual void Event(BrEventNode* inNode, BrEventNode* outNode) virtual Int_t GetMaxRowsMissed() const virtual Double_t GetSearchWidthX() const virtual Double_t GetSearchWidthY() const virtual void Init() virtual TClass* IsA() const virtual void Print(Option_t* option = "B") const virtual void SetMaxRowsMissed(Int_t val = 3) virtual void SetSearchWidthX(Double_t val = 0.2) virtual void SetSearchWidthY(Double_t val = 0.4) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members

    private:
Float_t fSearchWidthX search width in x direction Float_t fSearchWidthY search width in y direction Int_t fMaxRowsMissed max number of rows allowed to miss BrDetectorVolume* fTpcVolume TPC detector volume pointer TH1F* hHitSearchDx The pad hit distribtuion at the proj TH2F* hHitRowSearchDx The pad hit dist vs row at the proj TH1F* hHitSearchDy The time hit distribtuion at the proj TH2F* hHitRowSearchDy The time hit dist vs row at the proj TH1F* hHitAcceptedDx The pad hit dist for accepted hits TH2F* hHitRowAcceptedDx The pad hit dist vs row for acc hits TH1F* hHitAcceptedDy The time hit dist for accepted hits TH2F* hHitRowAcceptedDy The time hit dist vs row for acc hits

Class Description




BrTpcTrackFollowModule() : BrTpcTrackModule()
 Default constructor. DO NOT USE

BrTpcTrackFollowModule(const Char_t* name, const Char_t* title) : BrTpcTrackModule(name, title)
 Named Constructor

void DefineHistograms()
 Define histograms. They are:
 <fill in here>

void Init()
 Job-level initialisation
 Init BrTpcTrackModule as well

void Event(BrEventNode* inNode, BrEventNode* outNode)
 Event method.
 First local tables and clonesarrays are cleared.
 Then hits are filled in the clonesarray for each row..
 The algorithm loops over hits in the last row. looking at the
 next active row in a search window that covers the front of the
 TPC the following windows is determined by the fSearchWidthX - and
 Y - which are angles.

 After the tracking the outputnode is filled with the best tracks

 If the tracking fails because of too many local track candidates
 (> 5000) then the Status is set to 1 and can be accessed with
 GetStatus. 0 is okay.

void FindTracks()
 Loop over rows to find tracks.
 Going from the back of the Tpc to the front

void FindTracksFromRow(const Int_t istart, const Int_t ilast)
 This is the method that loops over the hits in a single row
 It calls BrTrackSingleHit for each hit

void TrackSingleHit(const Int_t istart, const Int_t ilast)
 This routine tracks a single hit through all the remaining rows
 and forks the track candidate into more if more than 1 hit is
 found in the search window

void FindHitsNextRow(const Int_t rowNo, const Float_t *proj, const Float_t dx1, const Float_t dx2, const Float_t dy1, const Float_t dy2, TObjArray &FoundHits)
 Loop over the hits in row rownumber-1 and accepts the hits in the
 interval
 proj[0] - dx1 < x < proj[0] + dx2 &&
 proj[1] - dy1 < y < proj[1] + dy2
 The hits are then stored in the array FoundHits
 A lot of diagnosic histograms are filled!

void Print(Option_t* option) const
 Print module information
 See BrModule::Print for options.
 In addition this module defines the Option:
 <fill in here>



Inline Functions


                         Int_t GetMaxRowsMissed() const
                      Double_t GetSearchWidthX() const
                      Double_t GetSearchWidthY() const
                          void SetMaxRowsMissed(Int_t val = 3)
                          void SetSearchWidthX(Double_t val = 0.2)
                          void SetSearchWidthY(Double_t val = 0.4)
                       TClass* Class()
                       TClass* IsA() const
                          void ShowMembers(TMemberInspector& insp, char* parent)
                          void Streamer(TBuffer& b)
                          void StreamerNVirtual(TBuffer& b)
        BrTpcTrackFollowModule BrTpcTrackFollowModule(BrTpcTrackFollowModule&)
                          void ~BrTpcTrackFollowModule()

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on 2002/01/03 19:53:42 $ by cholm $

Validate HTML
Validate CSS