BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page

BrMrsTofMatchingModule


class description - source file - inheritance tree

class BrMrsTofMatchingModule : public BrModule

    private:
void CheckTofCal() Int_t GetPanelId(Int_t slat) void InitGeo() void InitParams() Bool_t IsValid(Double_t x) public:
BrMrsTofMatchingModule BrMrsTofMatchingModule() BrMrsTofMatchingModule BrMrsTofMatchingModule(const Char_t* name, const Char_t* title) BrMrsTofMatchingModule BrMrsTofMatchingModule(BrMrsTofMatchingModule&) virtual void ~BrMrsTofMatchingModule() virtual void Begin() static TClass* Class() virtual void DefineHistograms() virtual void Event(BrEventNode* inNode, BrEventNode* outNode) virtual void Finish() virtual void Init() virtual TClass* IsA() const virtual void Print(Option_t* option = "B") const void SetDefaultParameters() void SetMaxAdc(Double_t max = 100000) void SetMaxDX(Float_t dx = 1.5) void SetNoPedWidth(Int_t n = 50) void SetNtuple(Bool_t b = kFALSE) void SetTdcRange(Double_t min = 10, Double_t max = 4000) void SetTpm2PanXOffset(Int_t pn, Float_t x = 0) void SetXHitCorr(Float_t da = 0.0095) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members

    private:
Float_t fMaxDX Float_t fTpm2PanXOffset[6] BrTofCalibration* fTofwCalib TOFW calibration parameter element BrDetectorVolume* fTpm2Vol tracking chamber volumes BrDetectorVolume** fPanelVol tof volumes BrDetectorParamsTof** fPanelPar detector parameters BrDetectorParamsTof* fTofwPar detector parameters Double_t fNoPedWidth for matching between hit and tracks only Double_t fMaxAdc for hit selection Double_t fMaxTdc idem Double_t fMinTdc idem Double_t fXhitCorr hitcorrection for TOFW (2001/02 run) Bool_t fNtuple Bool_t* fValidSlat Int_t fNumMatchedSlats TNtuple* fNtMatch TH1F* fDxAll TH1F* fDxMatch TH1F** fPanelDx TH2F** fX2dAll TH2F** fX2dMatch TH1F* fPanelStat TH1F* fMatchStat TH1F** fMatchPanStat TH1F* fSlatDiff TH1F* fAccptSlatDist TH1F* fAccptMatchedSlatDist

Class Description

 Class BrMrsTofMatchingModule

 Module scanning the global track table and picking up those
 matching valid tof hits.
 Selected pairs saved into table of BrTofTrackMatch named "MRS Matching"

 Algorithm (Event method):
 -------------------------

 1- pick up tables of TOFW hits and MRS tracks
    if one of them missing, return

 2-
    a- Select tof hits by checking the TDCs and ADCs of both tubes

    b- loop over tracks and sub-loop over hits
          i - evaluate the dX between track projection and hit position
         ii - minimize this distance and store the "winner" hit
              slat number in an array that has the dimension of
              the track table. It might be that the minimal dX is
              beyond the matching criterium, then, the hit id is -1

    c- check "multiple hits" when more than 1 track share the same
       hit. In that case, ignore these tracks (it doesn't remove
       completely multiple hits since there might be one track
       only pointing to such hits - the other tracks having not
       been reconstructed)

    d- loop over remaining matchings and store result into table
       of BrTofTrackMatch object (cf class implementation file)

    Note: due to the problem in the TPM2 tracks (slopes) that causes tracks to
    be off at the slats at the edge, an effective correction to the hit pos is introduced.
    to correct the Hit - used in finding the 'best' hit.

____________________________________________________________________

 $Id: BrMrsTofMatchingModule.cxx,v 1.9 2002/08/30 16:27:09 hagel Exp $
 $Author: hagel $
 $Date: 2002/08/30 16:27:09 $
 $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>


BrMrsTofMatchingModule() : BrModule()
 Default constructor. DO NOT USE

BrMrsTofMatchingModule(const Char_t* name, const Char_t* title) : BrModule(name, title)

void SetDefaultParameters()
 default parameters

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

void Init()
 Job-level initialisation

void InitParams()

void InitGeo()
 private

void CheckTofCal()

void Begin()
 Run-level initialisation

void Event(BrEventNode* inNode, BrEventNode* outNode)
 Per event method

void Finish()

Int_t GetPanelId(Int_t slat)
 -----------------------------------------
 returns TOFW panel id given a slat number
 -----------------------------------------

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


                        Bool_t IsValid(Double_t x)
                          void SetTdcRange(Double_t min = 10, Double_t max = 4000)
                          void SetTpm2PanXOffset(Int_t pn, Float_t x = 0)
                          void SetNoPedWidth(Int_t n = 50)
                          void SetNtuple(Bool_t b = kFALSE)
                          void SetMaxAdc(Double_t max = 100000)
                          void SetMaxDX(Float_t dx = 1.5)
                          void SetXHitCorr(Float_t da = 0.0095)
                       TClass* Class()
                       TClass* IsA() const
                          void ShowMembers(TMemberInspector& insp, char* parent)
                          void Streamer(TBuffer& b)
                          void StreamerNVirtual(TBuffer& b)
        BrMrsTofMatchingModule BrMrsTofMatchingModule(BrMrsTofMatchingModule&)
                          void ~BrMrsTofMatchingModule()

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on 2002/08/30 16:27:09 $ by hagel $

Validate HTML
Validate CSS