ROOT logo
HYDRA - THE HADES ANALYSIS PACKAGE » (UNKNOWN) » HParticlePairMaker

class HParticlePairMaker: public TObject

_HADES_CLASS_DESCRIPTION



 HParticlePairMaker

 Build vectors of freference candidates, fothers candidates and fpairs
 for the current event. As reference candidates marked as kIsLepton
 will be used.

 nectEvent() to be called inside eventloop after HParticleSorter selection.
 1. clear all internal structures
 2. fill vector freference by candidates which are flagged kIsLepton ( or kIsUsed if it is setup)
    All other candidates (disregarding kIsLepton or kIsUsed) are added to vector fothers
    To avoid unneeded pairs the user can set a userFilter by set setUserFilter(Bool_t (*userfilter)(HParticleCand*))
    In this case only candidates passing this function (return kTRUE) will be taken into account
    for freference and fothers.
 3. created all pair combinations between freference x freference (no double counting) and freference x fothers.
       If setDoSkippedFullCandPairs(kTRUE) is  set (default) in addition all pairs of full reconstructed
       particles ( inner+outer MDC + META) of fOthers are build (ffullrecoOthers x ffullrecoOthers +
       ffullrecoOthers x fnofullrecoOthers ). All pairs with full reco candidates of the event are build
       this way for investigation, even the ones which are skipped from reference by any criteria or double hit
       rejection.
    a. pass the selectPID1 and selectPID2 functions
       will get pid1 or pid2 assigned. Otherwise pid is -2(fake+), -1(fake-) (see HPhysicsConstants).
       if pid2 fails it will get pid assigned opposite to reference cand in pairCase10.
       fakes for pid2 will be assigned by polarity or opposite to reference cand.
    b. eClosePairSelect flags will be set for each pair (see hparticledef.h).
       This flags can be used later to classify the pairs.
    c. pairs are inserted in fpairs vector. With static HParticlePair::setDoMomentumCorrection(Bool_t doit)
       the enegryloss correction can be swithced on/off
    d. secondary vertex and vertex cut vars are calulated. Which  primary ebevnt vertex
       should be used has to be set by the user (default = kVertexParticle, see eVertex
       in hparticledef.h). The user can provide an own event vertex per event by calling setVertex().


 Bool_t  HParticleTool::evalPairsFlags(UInt_t flag,HParticlePair& pair)
 can be used to classify the pair according to the predefined pair cases
 in ePairCase + eClosePairSelect (hparticledef.h). There are dedicated filter
 functions provided to select from the pairs vector


Lepton Pair Cases Lepton Pair Cases
   Fig 1. All predefined pairCases for Leptons


Hadron Pair Cases
   Fig 2. All predefined pairCases for Hadrons


Mixed Lepton Hadron Pair Cases
   Fig 2. All predefined pairCases for mixed Lepton Hadrons



 USAGE:
 #include  "hparticlepairmaker.h"
 #include  "hparticledef.h"
 #include  "hparticletool.h"
 #include  <vector>
 using namespace std;



 {
 // setup HLoop etc ....


  HEnergyLossCorrPar enLossCorr;
  enLossCorr.setDefaultPar("apr12"); //  "jan04","nov02","aug04","apr12"

 //---------------------------------------------
 // HPARTICLEPAIRMAKER SETUP :
 HParticlePairMaker pairmaker;
 // pairmaker.setPIDs(2,3,51) // default
 // pairmaker.setPIDsSelection(mySelectionPID1,mySelectionPID2) // set your own selection function for pid1+pid2 (signature :  Bool_t myfunction(HParticleCand*) )
 // pairmaker.setUseLeptons(kTRUE);// kTRUE : use kIsLepton, kFALSE: kIsUsed for selection of reference  (default kTRUE)
 // pairmake.setRequireRich(kFALSE); // do not ask for rich Hit in selection functions ( default = kTRUE)
 // pairmaker.setVertexCase(Particle::kVertexParticle); // select which event vertex to use (default = kVertexParticle)
 // HParticlePair::setDoMomentumCorrection(kTRUE); // default : kTRUE
 //---------------------------------------------

 vector<HParticlePair*> pairs;

 for (Int_t i=0; i < entries; i++) {
    Int_t nbytes =  loop.nextEvent(i);             // get next event. categories will be cleared before
    if(nbytes <= 0) { cout<<nbytes<<endl; break; } // last event reached


    //--------------------------------------------------------------------------
    // HPARTICLETRACKSORTER SETUP
    sorter.cleanUp();
    sorter.resetFlags(kTRUE,kTRUE,kTRUE,kTRUE);     // reset all flags for flags (0-28) ,reject,used,lepton
    Int_t nCandLep     = sorter.fill(HParticleTrackSorter::selectLeptons);   // fill only good leptons
    Int_t nCandLepBest = sorter.selectBest(HParticleTrackSorter::kIsBestRKRKMETA,HParticleTrackSorter::kIsLepton);
    Int_t nCandHad     = sorter.fill(HParticleTrackSorter::selectHadrons);   // fill only good hadrons (already marked good leptons will be skipped)
    Int_t nCandHadBest = sorter.selectBest(HParticleTrackSorter::kIsBestRKRKMETA,HParticleTrackSorter::kIsHadron);
    //--------------------------------------------------------------------------


    //###########################################################################
    // HPARTICLEPAIRMAKER ACTION:
    pairmaker.nextEvent();  // call it after your HParticleSorter action!



    //--------------------------------------------------------------------------
    // example  1:
    // fill list of pairs for pairCase1 and require
    // fitted inner+out+RK for second leg in addition
    pairmaker.filterPairsVector(pairs,kPairCase1|kFittedInnerMDC2|kFittedInnerMDC2|kRK2);

    for(UInt_t j = 0; j < pairs.size(); j++ ){ // print the pairs
       pairs[j]->print();
    }
    //--------------------------------------------------------------------------


    //--------------------------------------------------------------------------
    // example 2:
    // get list of all pairs and loop over the list
    // plot different pairCases

    vector<HParticlePair>& pairsVec =  pairmaker.getPairsVector();
    for(UInt_t l = 0; l < pairsVec.size(); l++ ){
        HParticlePair& pair = pairsVec[l];

        if(HParticleTool::evalPairsFlags(kPairCase1,pair)){
    	cout<<"case1"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase2,pair)){
    	cout<<"case2"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase3,pair)){
    	cout<<"case3"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase4,pair)){
    	cout<<"case4"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase5,pair)){
    	cout<<"case5"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase6,pair)){
    	cout<<"case6"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase7,pair)){
    	cout<<"case7"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase8,pair)){
    	cout<<"case8"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase9,pair)){
    	cout<<"case9"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase10,pair)){
    	cout<<"case10"<<endl;
        } else if(HParticleTool::evalPairsFlags(kPairCase11,pair)){
   	cout<<"case11"<<endl;
        }
    }
    //--------------------------------------------------------------------------
    //###########################################################################
   }  // end event loop
 }

Function Members (Methods)

public:
HParticlePairMaker()
HParticlePairMaker(const HParticlePairMaker&)
virtual~HParticlePairMaker()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTObject::Clear(Option_t* = "")
virtual TObject*TObject::Clone(const char* newname = "") const
virtual Int_tTObject::Compare(const TObject* obj) const
virtual voidTObject::Copy(TObject& object) const
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
Int_tfilterCandidates(HVirtualCand* cand, vector<HVirtualCand*>& candidates, UInt_t flag = 0, Float_t oAngle = -1)
Int_tfilterCandidates(HVirtualCand* cand, vector<HParticlePair*>& filterpairs, UInt_t flag = 0, Float_t oAngle = -1)
voidfilterPairsVector(vector<HParticlePair*>& filterpairs, UInt_t flag = 0)
voidfilterPairsVector(vector<HParticlePair*>& filterpairs, vector<UInt_t>& flags)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
virtual const char*TObject::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
vector<HParticleCand*>&getOthersVector()
vector<HParticlePair>&getPairsVector()
vector<HParticleCand*>&getReferenceVector()
static Bool_tgetRequireRich()
Int_tgetSameInnerMdc(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
Int_tgetSameMeta(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
Int_tgetSameOuterMdc(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
Int_tgetSameRich(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
virtual const char*TObject::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTObject::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTObject::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTObject::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
voidnextEvent()
virtual Bool_tTObject::Notify()
voidTObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
HParticlePairMaker&operator=(const HParticlePairMaker&)
virtual voidTObject::Paint(Option_t* option = "")
voidplotPairCaseStat()
virtual voidTObject::Pop()
virtual voidTObject::Print(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
static Bool_tselectNeg(HParticleCand*)
static Bool_tselectPos(HParticleCand*)
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
voidsetDoSkippedFullCandPairs(Bool_t doit)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
static voidTObject::SetObjectStat(Bool_t stat)
voidsetPIDs(Int_t pid1, Int_t pid2, Int_t motherpid)
voidsetPIDsSelection(Bool_t (*)(HParticleCand*) selPID1, Bool_t (*)(HParticleCand*) selPID2)
static voidsetRequireRich(Bool_t use)
virtual voidTObject::SetUniqueID(UInt_t uid)
voidsetUseLeptons(Bool_t use)
voidsetUserFilter(Bool_t (*)(HParticleCand*) userfilter)
voidsetVertex(HGeomVector& v)
voidsetVertex(HVertex& v)
voidsetVertexCase(Particle::eVertex vertexCase)
virtual voidShowMembers(TMemberInspector&)
virtual voidStreamer(TBuffer&)
voidStreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidTObject::MakeZombie()
private:
voidbookHits(HParticleCand* cand1)
voidclearVectors()
voidselectPID(HParticleCand* cand1, Int_t& pid1, Bool_t warn = kTRUE)

Data Members

private:
vector<UInt_t>fCaseCt! counter array for cases
vector<UInt_t>fCaseVec! vector for pair cases
Int_tfMotherPID! default dilepton
Int_tfPID1! pid1 (default positrons)
Int_tfPID2! pid2 (default electrons)
HGeomVectorfVertexvertex for current event
Int_tfVertexCase! which eventvertex to use (see eVertex in hparticledef.h)
Bool_tfdoSkippedFullCandPairs! == kTRUE build also pairs of skipped full reco cands (inner/outer MDC+META) with others
vector<HParticleCand*>ffullrecoOthers! full reco cands (inner/outer MDC + META) inside others
vector<HParticleCand*>fnofullrecoOthers! not full reco cands (inner/outer MDC or META missing) inside others
vector<HParticleCand*>fothers! other candidates (not KIsLepton flagged)
vector<HParticlePair>fpairs! all pair combinations freference x fothers
vector<HParticleCand*>freference! reference candidates (kIsLepton flagged)
static Bool_tfrequireRich! ask for rich index in selctPos/selectNeg function
Bool_t (*)(HParticleCand*)fselectPID1! selection function pid1 (default positrons)
Bool_t (*)(HParticleCand*)fselectPID2! selection function pid2 (default electrons)
Bool_tfuse_kIsLepton! == kTRUE use kIsLepton as refererence selection (default)
Bool_t (*)(HParticleCand*)fuserFilter! user filter function to avoid unneeded combinatorics
map<HVirtualCand*,vector<HParticlePair*> >mCandtoPair! candidate lookup candidate -> list of pairs using this candidate
map<Int_t,vector<HParticleCand*> >mEmctoCand! EMC Cluster lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mInnerMdctoCand! inner Seg lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mOuterMdctoCand! outer Seg lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mRichtoCand! RICH hit lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mRpcClsttoCand! RPC cluster lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mShowertoCand! SHOWER hit lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mTofClsttoCand! TOF cluster lookup detector hit ind -> list of candidates using this hit
map<Int_t,vector<HParticleCand*> >mTofHittoCand! TOF hit lookup detector hit ind -> list of candidates using this hit
Int_trichCandCt! counter for all pair cases with both candidates matching a Rich (check)

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

HParticlePairMaker()
~HParticlePairMaker()
Bool_t selectPos(HParticleCand* )
 standard function for positrons
 require richmatch (broad window,  frequireRich =kTRUE)
 Select positive polatity
Bool_t selectNeg(HParticleCand* )
 standard function for electrons
 require richmatch (broad window,  frequireRich =kTRUE)
 Select negative polatity
void bookHits(HParticleCand* cand1)
 add candidate hits to maps
void selectPID(HParticleCand* cand1, Int_t& pid1, Bool_t warn = kTRUE)
void nextEvent()
 to be called inside eventloop after HParticleSorter selection.
 1. clear all internal structures
 2. fill vector freference by candidates which are flagged kIsLepton
    All other candidates are added to vector fothers
 3. created all pair combinations between freference and fothers which
    a. pass the selectPID1 and selectPID2 functions
       will get pid1 or pid2 assigned. Otherwise pid is -2(fake+), -1(fake-) (see HPhysicsConstants).
       if pid2 fails it will get pid assigend opposite to reference cand in pairCase10.
       fakes for pid2 will be assigned by polarity or opposite to reference cand.
    b. eClosePairSelect flags will be set for each pair.
       This flags can be used later to classify the pairs.
    c. pairs are inserted in fpairs vector.
void clearVectors()
 clear all internal vectors etc
void filterPairsVector(vector<HParticlePair*>& filterpairs, UInt_t flag = 0)
 fill all pairs which fullfill the flag to filterspairs.
 filterpairs will be cleared automatically before filling.
void filterPairsVector(vector<HParticlePair*>& filterpairs, vector<UInt_t>& flags)
 fill all pairs which fullfill the flags to filterspairs.
 filterpairs will be cleared automatically before filling.
Int_t filterCandidates(HVirtualCand* cand, vector<HVirtualCand*>& candidates, UInt_t flag = 0, Float_t oAngle = -1)
 vector candidates with all candidates which have been combined into pairs.
 if flag == 0 no filtering will be done. Filtering takes place on pairs.
 Flag can be used to filter for candidates which share 1 or more detector
 hits with the candidate (see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 In addition a filter of opening Angle of the track can be applied (oAngle < 0 does no selection)
 The function returns the number of found candidates or -1 if the candidate
 has not been found at all (should not happen).
Int_t filterCandidates(HVirtualCand* cand, vector<HParticlePair*>& filterpairs, UInt_t flag = 0, Float_t oAngle = -1)
 Fills vector filterpairs with all pairs which share cand.
 if flag == 0 no filtering will be done. Filtering takes place on pairs.
 Flag can be used to filter for candidates which share 1 or more detector
 hits with the candidate (see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 In addition a filter of opening Angle of the track can be applied (oAngle < 0 does no selection)
 The function returns the number of found pairs or -1 if the candidate
 has not been found at all (should not happen).
Int_t getSameRich(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
 fills list of candidates which share the same RICH. If the RICH is not found
 at all -1 is returned, other wise the number of candidates found for the hit.
 The input candidate will be not included in candidates. Candidates can be filtered
 by flag ( 0 == no filter, see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 The input candidate will be treated as refererence if isReference == kTRUE.
Int_t getSameInnerMdc(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
 fills list of candidates which share the same inner MDC seg. If the seg is not found
 at all -1 is returned, other wise the number of candidates found for the hit.
 The input candidate will be not included in candidates. Candidates can be filtered
 by flag ( 0 == no filter, see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 The input candidate will be treated as refererence if isReference == kTRUE.
Int_t getSameOuterMdc(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
 fills list of candidates which share the same outer MDC seg. If the seg is not found
 at all -1 is returned, other wise the number of candidates found for the hit.
 The input candidate will be not included in candidates. Candidates can be filtered
 by flag ( 0 == no filter, see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 The input candidate will be treated as refererence if isReference == kTRUE.
Int_t getSameMeta(HParticleCand* cand, vector<HParticleCand*>& candidates, UInt_t flag = 0, Bool_t isReference = kTRUE)
 fills list of candidates which share the same META hit. The number of candidates
 found for the hit is returned or -1 if the hit is not found at all..
 The input candidate will be not included in candidates. Candidates can be filtered
 by flag ( 0 == no filter, see flags eClosePairSelect+ ePairCase (hparticledef.h)).
 The input candidate will be treated as refererence if isReference == kTRUE.
void plotPairCaseStat()
 draws the fraction of pair cases for
 the default pairCases1...12.
HParticlePairMaker()
void setPIDs(Int_t pid1, Int_t pid2, Int_t motherpid)
 setters
void setPIDsSelection(Bool_t (*)(HParticleCand*) selPID1, Bool_t (*)(HParticleCand*) selPID2)
void setUserFilter(Bool_t (*)(HParticleCand*) userfilter)
void setDoSkippedFullCandPairs(Bool_t doit)
void setUseLeptons(Bool_t use)
{ fuse_kIsLepton = use;}
void setRequireRich(Bool_t use)
{ frequireRich = use;}
Bool_t getRequireRich()
{ return frequireRich ;}
void setVertexCase(Particle::eVertex vertexCase)
{ fVertexCase = vertexCase; }
void setVertex(HGeomVector& v)
{ fVertex = v; fVertexCase = kVertexUser;}
void setVertex(HVertex& v)
{ fVertex = v.getPos(); fVertexCase = kVertexUser;}
vector<HParticleCand*>& getReferenceVector()
{ return freference; }
vector<HParticleCand*>& getOthersVector()
{ return fothers; }
vector<HParticlePair>& getPairsVector()
{ return fpairs; }