BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page

BrRichPidModule


class description - source file - inheritance tree

class BrRichPidModule : public BrModule


    protected:
virtual Float_t CalculateMass(Float_t p, Float_t radius) virtual Float_t EstimateNumberOfPhotons(Float_t radius) virtual Float_t FindRing(BrRichPidModule::RICHHIT* richHits, BrRichPidModule::TRACKHIT track) virtual Float_t GuessRadius(BrRichPidModule::RICHHIT* richHits, Int_t trackNumber) virtual BrVector3D InImagePlaneCoordinates(const BrVector3D& hit) virtual BrVector3D IntersectionWithMirror(const BrLine3D& trackLine) virtual BrLine3D ReflectTrackLineOnMirror(const BrLine3D& trackLine) virtual void Snapshot(BrRichPidModule::RICHHIT* richHits, BrRichPidModule::TRACKHIT track) public:
BrRichPidModule BrRichPidModule() BrRichPidModule BrRichPidModule(const Char_t* name, const Char_t* title) BrRichPidModule BrRichPidModule(BrRichPidModule&) virtual void ~BrRichPidModule() virtual void Begin() static TClass* Class() virtual void DefineHistograms() virtual void End() 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 SetCheckGeometry(Bool_t b = kFALSE) void SetDRWindow(Float_t dr = 2.0) void SetImagePlaneOffset(Float_t x = 0.185, Float_t y = -0.083) void SetMaxNSnapshots(Int_t n = 0) void SetMaxRadius(Float_t r = 11) void SetMinHitsInRing(Int_t n = 4) void SetRefractiveIndex(Float_t n = 1.00202) void SetUseD3Momentum(Bool_t d3 = kFALSE) void SetUseD4Momentum(Bool_t d4 = kFALSE) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members

    private:
BrDetectorVolume* fRichVolume RICH volume BrDetectorVolume* fT5Volume T5 volume BrChkvParameters* fRichParams RICH parameters BrPlane3D fImagePlane image plane BrVector3D fImagePlaneCenter center of image plane BrVector3D fImagePlaneOffset BrVector3D fMirrorSphereCenter center of mirror Float_t fMirrorSphereRadius radius of mirror sphere Float_t fMirrorFocalLength focal length of mirror Float_t fMirrorFocalAngle focal angle Float_t fRefractiveIndex refractive index of gas Float_t fDRWindow Int_t fMinHitsInRing Minimum n hits in ring Float_t fMaxRadius Maximum ring radius Int_t fNumberOfHits Number of tubes with signal Int_t fNumberOfTracks Number of tracks Bool_t fUseD3Momentum Bool_t fUseD4Momentum Bool_t fCheckGeometry if on, geo check hist is created TF1* fGausFit TDirectory* fHistDir Int_t fMaxNSnapshots max number of snapshots Int_t fNSnapshots snapshots taken TList* fSnapshots list of snapshot hist Int_t fNDEDR TList* fDEDR ring fitting hists TH2F* fRadVsP radius vs momentum TH2F* fRadVsP1 radius vs momentum, ntracks=1 TH2F* fRadVsP2 radius vs momentum, ntracks>1 TH2F* fMassVsP mass vs momnetum TH2F* fMassVsRad mass vs radius TH1F* fMass TH2F* fRingCenters x,y on image plane TH2F* fDRVsPhi TH2F* fDRVsPhiAcc TH2F* fRing TH2F* fRingAcc TH1F* fHitsInRing number of tubes hit in ring TH2F* fHitsVsRadius TH2F* fHitsVsPhotons TH2F* fHitsVsEstPhotons TF1* fElectronRvsP TF1* fMuonRvsP TF1* fPionRvsP TF1* fKaonRvsP TF1* fProtonRvsP

Class Description

 RICH pid module - makes particle identification from RICH data.

 Find the Bfs track table and the RichRdo object in the inNode.
 Calculate the track lines reflection on the RICH mirror and find
 the intersection with the image plane - this should be the center
 of a ring generated by the particle.

 Find a ring by looking at the number of hits in delta radius areas
 around the center.


BrRichPidModule()
 Default constructor. DO NOT USE

BrRichPidModule(const Char_t* name, const Char_t* title) : BrModule(name, title)
 Named Constructor

void DefineHistograms()
 Define histograms.

void Init()
 In this method the following is initialized:

 1) Geometry and Parameter db managers
 2) Rich parameters (tube numbering and positions)
 3) The Geometry of the RICH and T5
 4) Rich mirror information (focal length, mirror sphere radius...)
 5) Image plane

void Begin()
 Run-level initialisation

void Event(BrEventNode* inNode, BrEventNode* outNode)
 Following is done in the event method:


BrLine3D ReflectTrackLineOnMirror(const BrLine3D& trackLine)
 Reflects the track line on the focal mirror.

 1) find the reflection axis (the axis defined by the center of
    the mirror sphere and the track-mirror intersection point).

 2) find angle between track line and reflection axis

 3) find the rotation axis (perpendicular to trackLine and
    reflection axis)

 4) rotate trackline with 2*angle around reflection axis

BrVector3D IntersectionWithMirror(const BrLine3D& trackLine)
 Calculates the intersection between the trackLine and the focal
 mirror of the RICH (Djamels calculations).

 find intersection point with sphere:

 line equations:     sphere:
   x = ax*t + xo      (x - xs)^2 + (y - ys)^2 + (z - zs)^2 = r^2
   y = ay*t + yo
   z = az*t + zo

   Equation to solve:

   (axt + xo - xs)^2 + (ayt + yo - ys)^2 + (azt + zo - zs)^2 - r^2 = 0

   t^2(ax^2 + ay^2 + az^2)
        + 2t(ax*(xo - xs) + ay*(yo - ys) + az*(zo - zs))
              + ((xo - xs)^2 + (yo - ys)^2 + (zo - zs)^2 - r^2) = 0

   c1 t^2 + c2 t + c3 = 0;

BrVector3D InImagePlaneCoordinates(const BrVector3D& hit)
 Hit must be on image plane. Return hit coordinates (x,y,0) on
 Rich image plane.

Float_t GuessRadius(RICHHIT* richHits, Int_t trackNumber)
 Tries to find the radius of the ring.

 Seach for hits with distance to the ring center in a dR window and
 fills the number of hits into a histogram.

 Fits this histogram with a gaus to find the approximate maximum
 (radius) , if the maximum of the histogram is too far from the
 mean of the fit, simply use the max as radius.

Float_t FindRing(RICHHIT* richHits, TRACKHIT track)
 Rind a ring with center at (cX, cY) from hits.

 1) Requires nHits >= fMinHitsInRing
 2) Guesses a radius
 3) Checks that there is a reasonable number of hits in the ring
 4) If not, remove hits in ring from the hits array and try again
 5) Does this 5 times. If no rings found returns 0.
 6) If ring is found, refit the radius.

Float_t EstimateNumberOfPhotons(Float_t radius)
 The number of detected cherenkov photons is given by

 N = No * L sin^2(theta)

 where No is a detector quality parameter, L is the path length
 and theta is the cherenkov angle

Float_t CalculateMass(Float_t p, Float_t radius)
 Calculates the theoretical mass from the momentum and the ring
 radius.
                      ____________
 we have :  1/beta = V 1 + (m/p)^2   and   cos(theta) = 1/(n*beta)
                   _____________________
         => m = p V (n*cos(theta))^2 - 1


void Snapshot(RICHHIT* richHits, TRACKHIT track)

void End()
 Run-level finalisation

void Finish()
 Job-level finalisation

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


                   void SetDRWindow(Float_t dr = 2.0)
                   void SetMaxNSnapshots(Int_t n = 0)
                   void SetMinHitsInRing(Int_t n = 4)
                   void SetMaxRadius(Float_t r = 11)
                   void SetRefractiveIndex(Float_t n = 1.00202)
                   void SetUseD3Momentum(Bool_t d3 = kFALSE)
                   void SetUseD4Momentum(Bool_t d4 = kFALSE)
                   void SetCheckGeometry(Bool_t b = kFALSE)
                   void SetImagePlaneOffset(Float_t x = 0.185, Float_t y = -0.083)
                TClass* Class()
                TClass* IsA() const
                   void ShowMembers(TMemberInspector& insp, char* parent)
                   void Streamer(TBuffer& b)
                   void StreamerNVirtual(TBuffer& b)
        BrRichPidModule BrRichPidModule(BrRichPidModule&)
                   void ~BrRichPidModule()

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on 2002/06/10 13:30:57 $ by ekman $

Validate HTML
Validate CSS