BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page

BrMultRdoModule


class description - source file - inheritance tree

class BrMultRdoModule : public BrModule, public BrMultUtil


    protected:
virtual Double_t CalibrateAdc(Int_t adc, Bool_t& below) virtual void CalibrateArray(BrMultRdo* rdo) virtual Double_t CalibrateArrayMultiplicity() virtual Double_t CalibrateAvgEta(Double_t& ringScale) virtual Double_t CalibrateEnergy(Double_t calAdc) virtual Double_t CalibrateEta() virtual Double_t CalibrateMultiplicity(Double_t ringEta) virtual void CalibrateRings(BrMultRdo* rdo) virtual void CalibrateRingsCloseNeighbors(BrMultRdo* rdo) virtual void CalibrateRingsFractionArray(BrMultRdo* rdo) public:
virtual void ~BrMultRdoModule() virtual void Begin() static TClass* Class() virtual void DefineHistograms() virtual void Finish() virtual TClass* IsA() const virtual void Print(Option_t* option = "B") const virtual void SetAdcGapLimit(Int_t gapLimit = 0) virtual void SetOutlierMethod(BrMultRdoModule::EOutlierMethod m = kNoCorrection) virtual void SetSaturationChannel(Int_t channel = 1000000) virtual void SetSingleAdcComp(Float_t comp) virtual void SetSingleAdcHigh(Float_t high) virtual void SetSingleAdcLow(Float_t low) virtual void SetThresholdFactor(Float_t factor = 5) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
BrMultTmpCalibration* fCalibration BrMultParameters* fParameters const BrRunInfo* fCurrentRun TH2F* fSingleAdcHisto Single vs Adc histogram TH2F* fSingleCalAdcHisto Single vs Calibrated ADC histogram TH2F* fSingleEnergyHisto Single vs deposited energy histogram TH2F* fSingleMultHisto Single vs Multiplicity histogram TH2F* fRingCalAdcHisto Ring vs Calibraed ADC histogram TH2F* fRingEnergyHisto Ring vs deposited energy histogram TH2F* fRingMultHisto Ring vs Multiplicity histogram TH2F* fRingHitsHisto Ring vs # Hits histogram TH1F* fArrayCalAdcHisto Total Calibrated ADC histogram TH1F* fArrayEnergyHisto Total energy histogram TH1F* fArrayMultHisto Total Multiplicity histogram TH1F* fArrayHitsHisto Total # Hits histogram TH2F* fVtxVsCalAdcHisto TH2F* fVtxVsMultHisto TH1F* fEtaDistHisto psuedo dN/deta TH1F* fdNdEtaSumHisto dN/deta Sum TH1F* fdNdEtaFillHisto dN/deta Fills TH1F* fdNdEtaHisto dN/deta Normalised UInt_t fEtaFills Count how many times dN/deta was filled Int_t fAdcGapLimit Lower bound of ADCs gap Float_t fThresholdFactor Thresholdfactor for ADC's Int_t fSaturationChannel Saturation value for ADC's BrMultRdoModule::EOutlierMethod fOutlierMethod Method to correct for outliers Short_t fCurrentSingle Short_t fCurrentRing Double_t fCurrentRingZ Double_t fCurrentRingR Double_t fCurrentRingEta Double_t fCurrentVtxZ Float_t fSingleAdcLow Float_t fSingleAdcHigh Float_t fSingleAdcComp Float_t fBaseCalAdc Float_t fBaseMult Float_t fBaseEnergy public:
static const BrMultRdoModule::EOutlierMethod kCloseNeighbors static const BrMultRdoModule::EOutlierMethod kFractionArray static const BrMultRdoModule::EOutlierMethod kHitCorrection static const BrMultRdoModule::EOutlierMethod kNoCorrection


See also

BrSiRdoModule, BrTileRdoModule

Class Description




void DefineHistograms()
 Define some histograms:
   Dimension       Contents
 --------------------------------------------------
     2D            ADC vs Single
     2D            Energy vs Single
     2D            <Energy> vs Ring
     1D            Total <Energy>
     2D            Hits vs Ring
     1D            Hits

void Begin()

void Finish()
 Job-level finalisation
 Normalise the dN/deta histogram to number of events and bin size

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

void CalibrateRings(BrMultRdo* rdo)
 PRIVATE METHOD:
 Make the calibrated data for the individual detector groups i.e.,
 the rings.
 This method does the outlier correction as well. There are three
 ways of doing that:

   "Close Neighbors"
   See the method BrMultRdoModule::CalibrateRingsCloseNeighbors()

   "Fraction Array"
   See the method BrMultRdoModule::CalibrateRingsFractionArray()

   "None"
   No correction for outliers, just sum it.

 As a side effect, this method also calculates the sum energy and
 multiplicity over the full array, but a further correction need
 to be applied.

void CalibrateRingsCloseNeighbors(BrMultRdo* rdo)
 PRIVATE METHOD:
 Correct for outliers.
 The "Close Neighbors" method calculates the avarage and variance
 over a given ring and it's immediate neighbors. Then for each
 single in the ring and it's immediate neighbors it tags a single as
 an outlier, if it's energy signal is higher then the mean plus
 three variances. For those singles tag as outliers, the method
 substitute the energy with avarage over the the non-tagged singles
 in the ring and it's immediate neighbors.

void CalibrateRingsFractionArray(BrMultRdo* rdo)
 PRIVATE METHOD:
 Correct for outliers.
 The "Fraction Array" method tags single that have a calibrated ADC
 signal greater then 10% of the total array calibrated ADC signal
 as an outlier. The outliers energy signal is substituted with the
 mean over the rest of the array.

void CalibrateArray(BrMultRdo* rdo)
 PRIVATE METHOD:
 Make the calibrated data for the full detector.

Double_t CalibrateEta()
 PRIVATE METHOD:
 Calculate the pseudo-rapidity for the current ring.

Double_t CalibrateAvgEta(Double_t & ringScale)
PRIVATE METHOD
 Calculate average pseudorapidity accounting for depth of detectors.
 The method returns the geometric scaling factor (ringScale), where
 dN/deta  = ringScale * (# of primary hits).

         |--- width -------|

          -----------------                ---
         |                 |               |
         |    si or tile   |             Depth
         |                 |               |
          -----------------               ---  ---
                                               |
               --     delta -eta                |
                                               |
                         --                   front
                                               |
                                 --             |
                                               |
                    ------------------XX (vertex) --------------------

 Assuming the azmulthal acceptance of a given detector is delta-phi,
 the geometric factor to obtain dN/dEta for that detector is:

             ringScale = ( 2*pi/delta-phi ) * 1/delta-eta


Double_t CalibrateAdc(Int_t adc, Bool_t& below)
 PRIVATE METHOD:
 Calibrate the ADC from DAQ. That is, do the pedestal subtraction
 if above threshold.

Double_t CalibrateEnergy(Double_t calAdc)
 PRIVATE METHOD:
 The energy in a given single is calculated as

                       sin(theta + tilt)
  E = calADC * gain *  -----------------   (1)
                           sin(theta)

 The third factor in (1) is a correction for the tilt of the
 detector.
 Given that

   sin(theta + tilt) = sin(theta)*cos(tilt) + sin(tilt)cos(theta)

 and that

                          ringR
   sin(theta)   = -----------------------
                  sqrt(distZ^2 + ringR^2)

                          distZ
   cos(theta)   = -----------------------
                  sqrt(distZ^2 + ringR^2)

                  distZ
   cotan(theta) = -----
                  ringR


        ------- ---+--- ------- -------
                  /|
                 / | ringR
                /  |
        =======*===+==================
               distZ

 where distZ = ringZ - vtxZ, we get

   sin(theta - tilt)    sin(theta)*cos(tilt) - sin(tilt)cos(theta)
   ----------------- =  ------------------------------------------
       sin(theta)                         sin(theta)

                                   cos(theta)
                     = cos(tilt) - ---------- * sin(tilt)
                                   sin(theta)

                     = cos(tilt) - cotan(theta) * sin(tilt)

                                   distZ
                     = cos(tilt) - ----- * sin(tilt)
                                   ringR

 so we get


                                    distZ
  E = calADC * gain * ( cos(tilt) - ----- * sin(tilt) )
                                    ringR

 We still need to correct for the vertex position, but that is
 done in the Event method.


Double_t CalibrateMultiplicity(Double_t ringEta)
 PRIVATE METHOD:
 Calculate the Multiplicity from the energy deposited in a single



Inline Functions


           Double_t CalibrateArrayMultiplicity()
               void SetAdcGapLimit(Int_t gapLimit = 0)
               void SetThresholdFactor(Float_t factor = 5)
               void SetSaturationChannel(Int_t channel = 1000000)
               void SetOutlierMethod(BrMultRdoModule::EOutlierMethod m = kNoCorrection)
               void SetSingleAdcLow(Float_t low)
               void SetSingleAdcHigh(Float_t high)
               void SetSingleAdcComp(Float_t comp)
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void Streamer(TBuffer& b)
               void StreamerNVirtual(TBuffer& b)
               void ~BrMultRdoModule()

This page automatically generated by script docBrat by Christian Holm

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

Validate HTML
Validate CSS