#include "hmdcdedx2.h"

HMdcDeDx2


class description - source file - inheritance tree (.pdf)

class HMdcDeDx2 : public HParCond

Inheritance Chart:
TObject
<-
TNamed
<-
HParSet
<-
HParCond
<-
HMdcDeDx2

    protected:
void calcSegPoints(HMdcSeg*, Double_t&, Double_t&, Double_t&, Double_t&, Double_t&, Double_t&) UChar_t fillingInput(HMdcSeg** seg) UChar_t select(Float_t, Float_t, UChar_t, Float_t wind = -99.) void sort(Int_t) public:
HMdcDeDx2(const char* name = "MdcDeDx2", const char* title = "Mdc lookup for MDC dEdX calculation", const char* context = "MdcDeDx2Production") HMdcDeDx2(const HMdcDeDx2&) ~HMdcDeDx2() Float_t calcDeDx(HMdcSeg** seg, Float_t*, Float_t*, UChar_t*, Float_t*, UChar_t*, Int_t vers = 2, Int_t mod = 2, Bool_t useTruncMean = kTRUE, Float_t truncMeanWindow = -99.) static TClass* Class() virtual void clear() Bool_t createMaxPar(Bool_t print = kFALSE) Double_t dEdXToToT(Int_t s, Int_t m, Double_t angle, Double_t mindist, Double_t dEdX) static Double_t energyLoss(Int_t id, Double_t p, Double_t hefr = 0.6) static TGraph* energyLossGraph(Int_t id, Double_t hefr = 0.6, TString opt = p, Bool_t exchange = kFALSE, Int_t markerstyle = 8, Int_t markercolor = 2, Float_t markersize = 0.7) void findBin(Int_t m, Double_t* angle, Double_t* mindist, Int_t* abin, Int_t* dbin) Double_t* getArray(Int_t& size) Double_t getFuncMaxPar(Int_t s, Int_t m, Int_t abin, Int_t dbin) void getFuncMaxPar(Double_t* p) Double_t* getFuncPar(Int_t s, Int_t m, Int_t abin, Int_t dbin) void getFuncPar(Double_t* p) Double_t* getFuncWidthPar(Int_t s, Int_t m, Int_t abin, Int_t dbin) void getFuncWidthPar(Double_t* p) Double_t getHeFraction() Int_t getMinimumWires() Int_t getN_Angle() Int_t getN_Dist() Int_t getN_Param() Int_t getN_Shift_Param() virtual Bool_t getParams(HParamList*) Float_t getWindow() virtual Bool_t init(HParIo*, Int_t*) Bool_t initContainer() virtual TClass* IsA() const Double_t normalize(Int_t s, Int_t m, Double_t angle, Double_t mindist, Double_t ToT) HMdcDeDx2& operator=(const HMdcDeDx2&) void printParam(TString opt = all) virtual void putParams(HParamList*) Double_t scaledTimeAboveThreshold(HGeantKine* kine = 0, Double_t p = -1, Float_t t1 = -999, Float_t t2 = -999, Float_t* t2err = 0, Int_t s = 0, Int_t m = 0, Double_t alpha = 0, Double_t mindist = 0) void setDebugMode(Bool_t dodebug) void setFuncMaxPar(Int_t s, Int_t m, Int_t abin, Int_t dbin, Double_t val) void setFuncMaxPar(Double_t* p) void setFuncPar(Int_t s, Int_t m, Int_t abin, Int_t dbin, Double_t* p, Int_t size) void setFuncPar(Double_t* p) void setFuncWidthPar(Int_t s, Int_t m, Int_t abin, Int_t dbin, Double_t* p, Int_t size) void setFuncWidthPar(Double_t* p) void setHeFraction(Double_t fr) void setMinimumWires(Int_t minwires) void setUseCalibration(Bool_t ok) void setWindow(Float_t win) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) Double_t toTSigma(Int_t s, Int_t m, Double_t angle, Double_t mindist, Int_t shift = 0) Double_t toTTodEdX(Int_t s, Int_t m, Double_t angle, Double_t mindist, Double_t ToT) virtual Int_t write(HParIo*)

Data Members


    protected:
Double_t par[6][4][18][40][4] fitparams for tot->dedx Double_t shiftpar[6][4][18][40][2] width of tot distribution Double_t parMax[6][4][18][40] max val for tot in tot->dEdx static Float_t MaxMdcMinDist[4] max value for mindist Float_t MinDistStep[4] step size for mindist Float_t AngleStep step size for angle Double_t hefr fraction of helium of the gas mix Int_t minimumWires minimum required wires Float_t window window for truncating arround mean (units of sigma) HMdcSizesCells* sizescells ! pointer to HMdcSizesCells container HCategory* catcal ! pointer to mdc cal1 HLocation loccal ! location object of cal1 Bool_t isInitialized ! flag after init Int_t method ! method switch for filling input Int_t module ! method switch for filling input for module 1/2 of segment TArrayD measurements ! array of measurements Bool_t useCalibration ! use/don't use normalization table Int_t ctskipmod0 ! counter for wires skipped with t2<=-998 in mod0 of seg Int_t ctskipmod1 ! counter for wires skipped with t2<=-998 in mod1 of seg Bool_t debug ! kTRUE print infos of trun mean

Class Description

*-- AUTHOR : J. Markert

 HMdcDeDx

 This transformation class calculates the "pseudo dEdx" from t2-t1 (time above threshold)
 of all fired drift cells included in one HMdcSeg. The transformation is performed doing a
 normalization of the measured t2-t1 with impact angle and minimum distance to wire
 (MDCCAL2 parametrization) and the algorithm to normalize the single measurements and
 average over all cells included in the segment.
-------------------------------------------------------------------------
 Calculation of dEdx:
 Float_t calcDeDx(HMdcSeg*,Float_t*,Float_t*,UChar_t*,Float_t*,UChar_t*,
		    Int_t vers=2,Int_t mod=2, Bool_t useTruncMean=kTRUE, Float_t truncMeanWindow=-99.)
 calculates dedx of a MDC (or 2) segment by doing normalization for
 impact angle and distance from wire for the fired cells of
 the segment. The measurements are sorted in decreasing order,
 the mean and sigma is calculated. Afterwards the truncated mean
 is calculated and returned as dedx. Mean,sigma,number of wires before
 truncating,truncated mean,truncated mean sigma and number of wires
 after truncating can be obtained by the returned pointers.
 Different methods can be selected:
 vers==0 : Input filling from inner HMdcSeg
 vers==1 : Input filling from outer HMdcSeg
 vers==2 : Input filling from inner+outer HMdcSeg (default)
 mod==0  : calc dedx from 1st module in segment
 mod==1  : calc dedx from 2nd module in segment
 mod==2  : calc dedx from whole segment (default)
 useTruncMean: kTRUE (default) apply truncated Mean
 truncMeanWindow (unit: SIGMA RMS of mean TOT): -99 (default) use standard window
 returns -99 if nothing was calculated
-------------------------------------------------------------------------
 Settings of the truncated Mean method:
 The truncated mean is calculated according to the set window (in sigma units)
 (setWindow(Float_t window)) arround the mean. For the calculation of the truncated
 mean only drift cells with a dEdx inside the window arround the mean (calculated from
 all drift cells of the segment) will be taken into account.
 With setMinimumWires(Int_t minwires) the algorithm can be forced to keep
 a minimum number of wires. The method will stop removing drift cells from the active
 sample if the minimum number is reached.


HMdcDeDx2(const char* name,const char* title, const char* context) : HParCond(name,title,context)

~HMdcDeDx2()
 destructor

void clear()

void printParam(TString opt)
 prints the parameters of HMdcDeDx2 to the screen.

Bool_t init(HParIo* inp,Int_t* set)
 intitializes the container from an input

Bool_t createMaxPar(Bool_t print)
 create the maximum allowed value of ToT to keep
 the numerical calulations inside the range of double.
 The procedure loops over all s,m,abin,dbin to set
 the boundary automatically. This function needs
 the parameter array for the fit functions to be initialized
 before. The boundary value are stored inside the container.

Bool_t initContainer()

Int_t write(HParIo* output)
 writes the container to an output

void putParams(HParamList* l)
 Puts all params of HMdcDeDx2 to the parameter list of
 HParamList (which ist used by the io);

Bool_t getParams(HParamList* l)

void sort(Int_t nHits)
 Puts the measurement values into decreasing order.

Float_t calcDeDx(HMdcSeg* seg[2], Float_t* meanold,Float_t* sigmaold,UChar_t* nwire, Float_t* sigmanew,UChar_t* nwiretrunc, Int_t vers,Int_t mod, Bool_t useTruncMean, Float_t truncMeanWindow)
 calculates dedx of a MDC segment by doing normalization for
 impact angle and distance from wire for the fired cells of
 the segment. The measurements are sorted in decreasing order,
 the mean and sigma is calculated. Afterwards the truncated mean
 according to the set window (in sigma units) arround the mean
 is calculated and returned as dedx. Mean,sigma,number of wires before
 truncating,truncated mean,truncated mean sigma and number of wires
 after truncating can be obtained by the returned pointers.
 Different methods can be selected:
 vers==0 : Input filling from inner HMdcSeg
 vers==1 : Input filling from outer HMdcSeg
 vers==2 : Input filling from inner+outer HMdcSeg (default)
 mod==0  : calc dedx from 1st module in segment
 mod==1  : calc dedx from 2nd module in segment
 mod==2  : calc dedx from whole segment (default)
 useTruncMean: kTRUE (default) apply truncated Mean
 truncMeanWindow (unit: SIGMA RMS of mean TOT): -99 (default) use standard window
 returns -99 if nothing was calculated

UChar_t fillingInput(HMdcSeg* seg[2])
 Fills array of measurements from HMdcSeg and returns the number of fired wires in Segment
 The single measurements are normalized by impact angle and distance from wire. Only the first
 wires < MAX_WIRES are accepted.

UChar_t select(Float_t mean,Float_t sigma,UChar_t nWire,Float_t wind)
 loops over the measurement array with nWire values and puts values to
 -99 if the measurements are outside the wanted window arround the mean value mean
 (set with setWindow() (sigma units)).
 The procedure keeps at least a number of minimum wires set with setMinimumWires().

void setFuncPar(Int_t s,Int_t m,Int_t abin,Int_t dbin,Double_t* p,Int_t size)
 set the values for the dedx functions by s, m, angle bin, dist bin. Needs pointer
 to parameter array[N_PARAM]

void setFuncPar(Double_t* p)
 set all values for the dedx functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST][N_PARAM]

void getFuncPar(Double_t* p)
 set all values for the dedx functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST][N_PARAM]

void setFuncMaxPar(Int_t s,Int_t m,Int_t abin,Int_t dbin,Double_t p)
 set the values for the dedx functions by s, m, angle bin, dist bin.

void setFuncMaxPar(Double_t* p)
 set all values for the dedx functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST]

void getFuncMaxPar(Double_t* p)
 set all values for the dedx functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST]

void setFuncWidthPar(Int_t s,Int_t m,Int_t abin,Int_t dbin,Double_t* p,Int_t size)
 set the values fpr the dedx width functions by s, m, angle bin, dist bin. Needs pointer
 to parameter array[N_SHIFT_PARAM]

void setFuncWidthPar(Double_t* p)
 set all values for the dedx width functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST][N_SHIFT_PARAM]

void getFuncWidthPar(Double_t* p)
 set all values for the dedx width functions. Needs pointer
 to parameter array[6][4][N_ANGLE][N_DIST][N_SHIFT_PARAM]

void calcSegPoints(HMdcSeg * seg,Double_t& x1, Double_t& y1, Double_t& z1, Double_t& x2, Double_t& y2, Double_t& z2)
 calculates 2 coordinate points from segment

void findBin(Int_t m,Double_t* angle,Double_t* mindist,Int_t* abin,Int_t* dbin)

Double_t toTSigma(Int_t s,Int_t m,Double_t angle,Double_t mindist,Int_t shift)
 returns asymmetric width of ToT for single drift cell
 shift == 0 (default) : returns 0
       ==-1           : low   bound width of ToT distribution (sigma)
       == 1           : upper bound width of ToT distribution (sigma)

Double_t toTTodEdX(Int_t s,Int_t m,Double_t angle,Double_t mindist,Double_t ToT)
 calculation of ToT -> dEdX for single drift cell

Double_t dEdXToToT(Int_t s,Int_t m,Double_t angle,Double_t mindist,Double_t dEdX)
 calculation of dEdX -> ToT for single drift cell

Double_t normalize(Int_t s,Int_t m,Double_t angle,Double_t mindist,Double_t ToT)
 calibrate ToT :
 dEdX=TotTodEdX(t2-t1) ----> dEdXToToT(dEdX) for reference module,angle,distbin

Double_t energyLoss(Int_t id,Double_t p,Double_t hefr)
 Calculates the dEdx (MeV 1/g cm2) of an particle with GEANT ID id
 and momentum p (MeV) for He/i-butan gas mixture with He fraction hefr
 (He (hefr) / i-butan (1-hefr)). The fomular is taken from PDG and doesn't
 include the density correction term. The values for the mean excitation
 energy are taken from Sauli.

TGraph* energyLossGraph(Int_t id,Double_t hefr,TString opt,Bool_t exchange,Int_t markerstyle,Int_t markercolor,Float_t markersize)
 creates a TGraph for particle with GEANT ID id and
 and He/i-butan gas mixture with He fraction hefr
 (He (hefr) / i-butan (1-hefr) dEdx vs p,beta,1/beta2 or betagamma
 depending on opt (p,beta,1/beta2,betagamma). exchange=kTRUE
 will exchange x-axis and y-axis.

Double_t scaledTimeAboveThreshold(HGeantKine* kine, Double_t p, Float_t t1, Float_t t2, Float_t* t2err, Int_t s,Int_t m, Double_t alpha,Double_t mindist)
 calculated the t2 of drift cell measurent
 from dedx of a given particle with momentum p.
 The assumption is that all drift cell measurements
 belonging to one track can be treated independent.
 return t2 and smeared error to pointer t2err for an
 single measurement. If the result t2-t1 would
 be < 0 a random value for t2 0-20 ns larger than
 t1(already including error) will be created.
 The smeared error will be returned to the pointer t2err.



Inline Functions


               void getFuncPar(Double_t* p)
               void getFuncMaxPar(Double_t* p)
               void getFuncWidthPar(Double_t* p)
              Int_t getN_Param()
              Int_t getN_Shift_Param()
              Int_t getN_Angle()
              Int_t getN_Dist()
            Float_t getWindow()
               void setWindow(Float_t win)
              Int_t getMinimumWires()
               void setMinimumWires(Int_t minwires)
               void setDebugMode(Bool_t dodebug)
               void setUseCalibration(Bool_t ok)
          Double_t* getArray(Int_t& size)
               void setHeFraction(Double_t fr)
           Double_t getHeFraction()
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void Streamer(TBuffer& b)
               void StreamerNVirtual(TBuffer& b)
          HMdcDeDx2 HMdcDeDx2(const HMdcDeDx2&)
         HMdcDeDx2& operator=(const HMdcDeDx2&)


Last update: Fri Jan 26 12:05:45 2007


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.