BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page

BrModuleMatchTrack


class description - source file - inheritance tree

class BrModuleMatchTrack : public BrModule

    private:
Bool_t GetFiducialStatus(BrDetectorTrack* front, BrDetectorTrack* back) Bool_t MatchesCuts(Double_t dang, Double_t dy, Double_t daly, Double_t p) protected:
void CleanUpSharedBackTracks() void CleanUpSharedFrontTracks() void CleanUpTracks() void ClearClonesArrays() void CreateClonesArrays() virtual void DefineHistograms() void SetDefaultParameters() public:
BrModuleMatchTrack BrModuleMatchTrack() BrModuleMatchTrack BrModuleMatchTrack(const Char_t*, const Char_t*) BrModuleMatchTrack BrModuleMatchTrack(const Char_t*, const Char_t*, const Char_t*, const Char_t*, const Char_t*) BrModuleMatchTrack BrModuleMatchTrack(BrModuleMatchTrack&) virtual void ~BrModuleMatchTrack() static TClass* Class() void Clear() void CombineFrontBack(BrEventNode* InputEventNode) virtual void Event(BrEventNode* InputEventNode, BrEventNode* OutputEventNode) virtual void Finish() BrDetectorVolume* GetBackVolume() const Double_t GetFiducialCutDx() const Double_t GetFiducialCutDy() const BrDetectorVolume* GetFrontVolume() const virtual Char_t* GetHeaderVersion() const BrMagnetVolume* GetMagnetVolume() const Double_t GetMatchDalyOffset() const Double_t GetMatchDangOffset() const Double_t GetMatchDyOffset() const TClonesArray* GetMatchedTracks() const Double_t GetSigmaCut() const Double_t GetSigmaDaly() const Double_t GetSigmaDang() const Double_t GetSigmaDy() const Double_t GetSigmaRange() const virtual Char_t* GetSourceVersion() const virtual void Init() virtual TClass* IsA() const void ListMatchingParameters() const virtual void Print(Option_t* option = "B") const void SetFiducialCutDx(Double_t value) void SetFiducialCutDy(Double_t value) void SetHistMomentumRange(Double_t range = 20.0) void SetMatchDalyOffset(Double_t value) void SetMatchDangOffset(Double_t value) void SetMatchDyOffset(Double_t value) void SetNtuple(Bool_t n = kFALSE) void SetOldFiducialCuts() void SetRotationFactors(Double_t f = 1.0, Double_t b = 1.0) void SetSigmaCut(Double_t val) void SetSigmaDaly(Double_t val) void SetSigmaDang(Double_t val) void SetSigmaDy(Double_t val) void SetSigmaRange(Double_t val) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members

    private:
BrGeometryDbManager* fGeomDb ! Database manager descriptor BrDetectorVolume* fFrontVolume Volume Descriptor for front tr. chamber BrDetectorVolume* fBackVolume Volume Descriptor for back tr. chamber BrMagnetVolume* fMagnetVolume Volume Descriptor for Magnet Double_t fFiducialCutdx Fiducal volume cut X (from edge of magnet at entrance/exit) Double_t fFiducialCutdy Fiducal volume cut Y Double_t fMatchDangOffset Effective edge matching angle offset Double_t fMatchDalyOffset Effective edge matching cosine alphay offset Double_t fMatchDyOffset Offset for dy matching (cm) Double_t fSigmaMatchedDang Fitted sigma for match parameter (Dang) used for selection Double_t fSigmaMatchedDaly Fitted sigma for match parameter (Daly) used for selection Double_t fSigmaMatchedDy Fitted sigma for match parameter (Dy) used for selection Double_t fSigmaRange how many sigmas used for initial selection Double_t fSigmaCut how many sigmas for final selection. Double_t fFrontRot Rotation of front tracks (to make dAng independent of p) Double_t fBackRot Rotation of back tracks (to make dAng independent of p) TClonesArray* fTracksFrontBack Results of track matching TString fFrontName TString fBackName TString fMagnetName Bool_t fNewFiducialCuts Int_t fNoEvents Int_t fMatchedTrackSum Double_t fHistMomentumRange define range for several diagnostic histograms TH1F* hFrontTracks TH1F* hBackTracks TH1F* hMatchedTracks TH1F* hGoodMatchedTracks TH1F* hMatchDyAll TH1F* hMatchDxAll TH1F* hMatchDalxAll TH1F* hMatchDalyAll TH1F* hMatchDangAll TH1F* hMatchDyAcpt TH1F* hMatchDxAcpt TH1F* hMatchDangAcpt TH2F* hMatchDangDxAcpt TH2F* hMatchDangDxAll TH2F* hMatchDxDyAll TH2F* hMatchDxDyAcpt TH1F* hMatchDalxAcpt TH1F* hMatchDslAngxAcpt TH1F* hMatchDalyAcpt TH2F* hMatchDalyDyAcpt TH2F* hMatchDalyDyAll TH2F* hMatchDangTrack TH2F* hMatchDangVsP TH2F* hMatchDyVsP TH1F* hMatchP TH1F* hMatchPl TH1F* hMatchPGood TH1F* hMatchMidX TH1F* hMatchMidY TH1F* hAllMidX TH1F* hAllMidY TH1F* hNoFrontMatch TH1F* hMatchChisqAll TH1F* hMatchChisq TNtuple* fMatchInfo correlation between all variables Bool_t fNtuple switch on/off ntuple Bool_t fFieldParametersRead

Class Description

 The BrModuleMatchTrack class provides an interface to the BRAHMS
 Tracking Combining of two detectors with a magnet in between.
 It is a general module combining any two tracking detectors
 surrounding a single dipole. eg. TPM1-D5-TPM2 or T1-D2-T2

 At present the module assumes the effective edge approximation
 for the magnetic field.

  Initial track matching is done base on 3 cut parameters, namely
  All the following cuts are sigmas in the variable. The actual cut
  in done in MatchesCuts(dang, dy, daly, p)
     SigmaMatchedDang : cut for difference in signed angle between tracks
                 alpah1-alpha2. Good tracks have DAlpha=0. This cut
                 is momentum and particle dependent.
     SigmaMatchedDy   : cut for vertical position difference on effective midplane.
     SigmaMatchedDaly : cut for vertical slope difference (take as dy/dz for track
                 front and back of magnet.
     SigmaRange: The actual cut is SigmaRange*Sigmaparameter.

     MatchDyOffset : Quite often small offset are present in the
                     tracking matching. This is introduced in order to make
                     symmetric cuts. The parameter is only used when
                     checking for the cuts limits.

  The calculations of tracks are now done in relative coordinates, and not
  in global as before. This implies changes for the members of the content of
  BrMatchedTrack.


Char_t* GetSourceVersion() const

BrModuleMatchTrack() : BrModule()
 Default Constructor. Do not use - Is only here for the purpose
 of various ROOT functionality.

BrModuleMatchTrack(const Char_t *Name, const Char_t *Title) : BrModule(Name,Title)
Set volume parameters to null,
Create clones arrays for internal use
Read simulation parameters from default file (simpar.dat) in current directory

BrModuleMatchTrack(const Char_t *Name, const Char_t *Title, const Char_t *tfront, const Char_t *tback, const Char_t* tmagnet) : BrModule(Name,Title)
 Constructor; Setup module.

 Get volume parameters from BrGeometryDbManager. This removes any
 need for the application to deal with the geometry. This module
 do not need to know where the Db information comes from. Please
 note that the geometry db objects are owned by the
 Geometry Manager and MUST not be deleted by the module.

 Create clones arrays for internal use
 Set default matching parameters.


void Init()

void DefineHistograms()
 define diagnostics histograms for track combination module   
 Book histograms (called from BrModule::Book)


void SetDefaultParameters()
 Set/Reset internal tracking parameters, cuts
 to default values.


void CreateClonesArrays()
Create the clones arrays used for internal structures.

void Clear()
 Clear event structures.

void ClearClonesArrays()
 Clear the clones arrays used for internal structures.


~BrModuleMatchTrack()
Destructor


void Event(BrEventNode* InputNode, BrEventNode* OutputNode)
 Normal Event Method for BrModule. The local detectortracks for the front and back
 chamber must be found in the InputNode. If matching is succesfull a table of
 BrMatchedTrack are added are added to Outputnode as well as a table for

void CombineFrontBack(BrEventNode *InputNode)
 Loop over tracks in TFront / TBack and see how well they match and
 create candidates for tracks
 This algorithm is the effective edge.
 For tracks matched they are added to the
 Clones array for matched BrTracksFrontBack (i.e. front back)


Bool_t GetFiducialStatus(BrDetectorTrack* front,BrDetectorTrack* back)
 This method provides fiducial cuts base on BrMagnetVolume::GetSwimStatus

void Finish()
 Finish
 Summarize performance

void CleanUpTracks()
 This cleanup method calls the two methods:
 CleanUpSharedBackTracks() and CleanUpSharedFrontTracks()

 CleanUpSharedBackTracks() finds mathced tracks that shares back
 tracks, find the best by looking at overall chisq in matching
 parameter, and set the status of the tracks (kOk for the best
 one, ...) - this is the old CleanUpTracks() method.

 CleanUpSharedFrontTracks() does the same as CleanUpSharedBackTracks()
 just for the front tracks. The only difference is that this method
 checks that the status is kOk before it compares the front tracks.

void CleanUpSharedBackTracks()
 The cleanup process implemented here is reasonably simple
 and incorporates the ideas used in the first analysis of
 the 2000 MRS data. The initial matching applies rather
 wide cuts in matching tracks. It can therefor well happen
 that e.g. the same backtrack is matched up with more than one
 front track. This happened in about 2-5% of the cases. The present
 algorithm sorts those out, or rather set the status after having selected
 the better one.

 The check for vertex cuts are not done here. It was deemed better to
 do so in a specific filter module.
 a) Identify tracks tracks that have multiple usage
 b) Identify best tracks by looking at overall chisq in matching
    parameter.
 c) at this time an overall estimat based on cut parameters and
    fidicual cuts are done.

void CleanUpSharedFrontTracks()

void ListMatchingParameters() const

Bool_t MatchesCuts(Double_t dang, Double_t dy, Double_t daly, Double_t p)
 Function to check if a matching fit is found between
 calculated parameters for deviate in angle, y and slopes.
 The (estimated) momentum is also a parameter to allow for better
 criteri with time.
 The cuts are detector, magnet and momentum dependent due to multiple
 scattering and detector resolution.
 The offset evaluation is applied here, as the only place in the code.

void Print(Option_t* option) const
 Print info on module



Inline Functions


                      void SetNtuple(Bool_t n = kFALSE)
             TClonesArray* GetMatchedTracks() const
         BrDetectorVolume* GetFrontVolume() const
         BrDetectorVolume* GetBackVolume() const
           BrMagnetVolume* GetMagnetVolume() const
                      void SetMatchDangOffset(Double_t value)
                      void SetMatchDalyOffset(Double_t value)
                      void SetMatchDyOffset(Double_t value)
                      void SetFiducialCutDx(Double_t value)
                      void SetFiducialCutDy(Double_t value)
                  Double_t GetSigmaDang() const
                  Double_t GetSigmaDaly() const
                  Double_t GetSigmaDy() const
                  Double_t GetMatchDangOffset() const
                  Double_t GetMatchDalyOffset() const
                  Double_t GetMatchDyOffset() const
                  Double_t GetSigmaRange() const
                  Double_t GetSigmaCut() const
                  Double_t GetFiducialCutDx() const
                  Double_t GetFiducialCutDy() const
                      void SetSigmaDang(Double_t val)
                      void SetSigmaDaly(Double_t val)
                      void SetSigmaDy(Double_t val)
                      void SetSigmaRange(Double_t val)
                      void SetSigmaCut(Double_t val)
                      void SetRotationFactors(Double_t f = 1.0, Double_t b = 1.0)
                      void SetHistMomentumRange(Double_t range = 20.0)
                      void SetOldFiducialCuts()
                   Char_t* GetHeaderVersion() const
                   TClass* Class()
                   TClass* IsA() const
                      void ShowMembers(TMemberInspector& insp, char* parent)
                      void Streamer(TBuffer& b)
                      void StreamerNVirtual(TBuffer& b)
        BrModuleMatchTrack BrModuleMatchTrack(BrModuleMatchTrack&)

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on 2002/08/08 15:36:13 $ by ufstasze $

Validate HTML
Validate CSS