HYDRA_development_version
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
hparticlepairmaker.h
Go to the documentation of this file.
1 #ifndef __HPARTICLEPAIRMAKER_H__
2 #define __HPARTICLEPAIRMAKER_H__
3 
4 
5 #include "hparticledef.h"
6 #include "hparticlepair.h"
7 #include "hparticlecand.h"
8 #include "hparticletool.h"
9 #include "hphysicsconstants.h"
10 #include "heventheader.h"
11 
12 
13 #include "TObject.h"
14 #include "TLorentzVector.h"
15 
16 
17 #include <vector>
18 #include <map>
19 #include <iostream>
20 #include <iomanip>
21 
22 using namespace std;
23 
24 class HParticlePairMaker : public TObject
25 {
26 private:
27 
28  vector<HParticleCand*> freference; //! reference candidates (kIsLepton flagged)
29  vector<HParticleCand*> fothers; //! other candidates (not KIsLepton flagged)
30  vector<HParticleCand*> ffullrecoOthers; //! full reco cands (inner/outer MDC + META) inside others
31  vector<HParticleCand*> fnofullrecoOthers; //! not full reco cands (inner/outer MDC or META missing) inside others
32  vector<HParticlePair> fpairs; //! all pair combinations freference x fothers
33 
34 
35  map<Int_t, vector<HParticleCand*> > mTofHittoCand; //! TOF hit lookup detector hit ind -> list of candidates using this hit
36  map<Int_t, vector<HParticleCand*> > mTofClsttoCand; //! TOF cluster lookup detector hit ind -> list of candidates using this hit
37  map<Int_t, vector<HParticleCand*> > mRpcClsttoCand; //! RPC cluster lookup detector hit ind -> list of candidates using this hit
38  map<Int_t, vector<HParticleCand*> > mShowertoCand; //! SHOWER hit lookup detector hit ind -> list of candidates using this hit
39  map<Int_t, vector<HParticleCand*> > mEmctoCand; //! EMC Cluster lookup detector hit ind -> list of candidates using this hit
40  map<Int_t, vector<HParticleCand*> > mInnerMdctoCand; //! inner Seg lookup detector hit ind -> list of candidates using this hit
41  map<Int_t, vector<HParticleCand*> > mOuterMdctoCand; //! outer Seg lookup detector hit ind -> list of candidates using this hit
42  map<Int_t, vector<HParticleCand*> > mRichtoCand; //! RICH hit lookup detector hit ind -> list of candidates using this hit
43  map<HVirtualCand*,vector<HParticlePair*> > mCandtoPair; //! candidate lookup candidate -> list of pairs using this candidate
44 
45 
46 
47  Bool_t (*fselectPID1)(HParticleCand*); //! selection function pid1 (default positrons)
48  Bool_t (*fselectPID2)(HParticleCand*); //! selection function pid2 (default electrons)
49  Bool_t (*fuserFilter)(HParticleCand*); //! user filter function to avoid unneeded combinatorics
50 
51  Int_t fPID1; //! pid1 (default positrons)
52  Int_t fPID2; //! pid2 (default electrons)
53  Int_t fMotherPID; //! default dilepton
54 
55  Bool_t fuse_kIsLepton; //! == kTRUE use kIsLepton as refererence selection (default)
56  Bool_t fdoSkippedFullCandPairs; //! == kTRUE build also pairs of skipped full reco cands (inner/outer MDC+META) with others
57  static Bool_t frequireRich; //! ask for rich index in selctPos/selectNeg function
58 
59  Int_t fVertexCase; //! which eventvertex to use (see eVertex in hparticledef.h)
60  HGeomVector fVertex; // vertex for current event
61 
62  vector<UInt_t> fCaseCt; //! counter array for cases
63  vector<UInt_t> fCaseVec; //! vector for pair cases
64  Int_t richCandCt; //! counter for all pair cases with both candidates matching a Rich (check)
65 
66  void clearVectors();
67  void bookHits(HParticleCand* cand1);
68  void selectPID(HParticleCand* cand1,Int_t& pid1,Bool_t warn=kTRUE);
69 public:
70 
73 
74  //--------------------------------------------------------------------------
75  // setters
76  void setPIDs(Int_t pid1,Int_t pid2,Int_t motherpid) {
77  fPID1 = pid1;
78  fPID2 = pid2;
79  fMotherPID = motherpid;
80  }
81 
82  void setPIDsSelection(Bool_t (*selPID1)(HParticleCand*), Bool_t (*selPID2)(HParticleCand*)) {
83  fselectPID1=selPID1;
84  fselectPID2=selPID2;
85  }
86  void setUserFilter(Bool_t (*userfilter)(HParticleCand*)) {
87  fuserFilter=userfilter;
88  }
89  void setDoSkippedFullCandPairs(Bool_t doit) { fdoSkippedFullCandPairs = doit; }
90  void setUseLeptons (Bool_t use) { fuse_kIsLepton = use;}
91  static void setRequireRich(Bool_t use) { frequireRich = use;}
92  static Bool_t getRequireRich() { return frequireRich ;}
93  void setVertexCase(Particle::eVertex vertexCase) { fVertexCase = vertexCase; };
94  void setVertex(HGeomVector& v) { fVertex = v; fVertexCase = kVertexUser;}
95  void setVertex(HVertex& v) { fVertex = v.getPos(); fVertexCase = kVertexUser;}
96 
97  static Bool_t selectPos(HParticleCand*);
98  static Bool_t selectNeg(HParticleCand*);
99 
100  //--------------------------------------------------------------------------
101  // event action
102  void nextEvent();
103 
104  vector<HParticleCand*>& getReferenceVector() { return freference; }
105  vector<HParticleCand*>& getOthersVector () { return fothers; }
106  vector<HParticlePair>& getPairsVector () { return fpairs; }
107 
108  //--------------------------------------------------------------------------
109  // filter functions
110  void filterPairsVector(vector<HParticlePair*>& filterpairs,UInt_t flag=0);
111  void filterPairsVector(vector<HParticlePair*>& filterpairs,vector<UInt_t>& flags);
112  Int_t filterCandidates (HVirtualCand* cand,vector<HVirtualCand*>& candidates,UInt_t flag=0,Float_t oAngle=-1);
113  Int_t filterCandidates (HVirtualCand* cand,vector<HParticlePair*>& filterpairs,UInt_t flag=0,Float_t oAngle=-1);
114 
115  //--------------------------------------------------------------------------
116  // lookup functions
117  Int_t getSameRich (HParticleCand* cand,vector<HParticleCand*>& candidates,UInt_t flag=0,Bool_t isReference = kTRUE);
118  Int_t getSameInnerMdc(HParticleCand* cand,vector<HParticleCand*>& candidates,UInt_t flag=0,Bool_t isReference = kTRUE);
119  Int_t getSameOuterMdc(HParticleCand* cand,vector<HParticleCand*>& candidates,UInt_t flag=0,Bool_t isReference = kTRUE);
120  Int_t getSameMeta (HParticleCand* cand,vector<HParticleCand*>& candidates,UInt_t flag=0,Bool_t isReference = kTRUE);
121 
122  void plotPairCaseStat();
123  ClassDef(HParticlePairMaker,0)
124 };
125 
126 
127 #endif // __HPARTICLEPAIRMAKER_H__
static Bool_t frequireRich
== kTRUE build also pairs of skipped full reco cands (inner/outer MDC+META) with others ...
vector< UInt_t > fCaseCt
static void setRequireRich(Bool_t use)
Int_t fPID1
user filter function to avoid unneeded combinatorics
void setVertex(HVertex &v)
Int_t fVertexCase
ask for rich index in selctPos/selectNeg function
void setVertex(HGeomVector &v)
vector< HParticleCand * > ffullrecoOthers
other candidates (not KIsLepton flagged)
map< Int_t, vector< HParticleCand * > > mInnerMdctoCand
EMC Cluster lookup detector hit ind -> list of candidates using this hit.
vector< HParticlePair > fpairs
not full reco cands (inner/outer MDC or META missing) inside others
map< Int_t, vector< HParticleCand * > > mShowertoCand
RPC cluster lookup detector hit ind -> list of candidates using this hit.
map< Int_t, vector< HParticleCand * > > mRpcClsttoCand
TOF cluster lookup detector hit ind -> list of candidates using this hit.
HGeomVector fVertex
which eventvertex to use (see eVertex in hparticledef.h)
map< HVirtualCand *, vector< HParticlePair * > > mCandtoPair
RICH hit lookup detector hit ind -> list of candidates using this hit.
HGeomVector & getPos(void)
Definition: heventheader.h:23
Int_t fMotherPID
pid2 (default electrons)
static Bool_t getRequireRich()
map< Int_t, vector< HParticleCand * > > mRichtoCand
outer Seg lookup detector hit ind -> list of candidates using this hit
void setDoSkippedFullCandPairs(Bool_t doit)
void setVertexCase(Particle::eVertex vertexCase)
Int_t richCandCt
vector for pair cases
vector< HParticleCand * > & getReferenceVector()
map< Int_t, vector< HParticleCand * > > mEmctoCand
SHOWER hit lookup detector hit ind -> list of candidates using this hit.
void setPIDs(Int_t pid1, Int_t pid2, Int_t motherpid)
map< Int_t, vector< HParticleCand * > > mTofClsttoCand
TOF hit lookup detector hit ind -> list of candidates using this hit.
Int_t fPID2
pid1 (default positrons)
vector< HParticleCand * > fothers
reference candidates (kIsLepton flagged)
Bool_t fdoSkippedFullCandPairs
== kTRUE use kIsLepton as refererence selection (default)
vector< HParticleCand * > freference
vector< HParticleCand * > & getOthersVector()
map< Int_t, vector< HParticleCand * > > mOuterMdctoCand
inner Seg lookup detector hit ind -> list of candidates using this hit
void setPIDsSelection(Bool_t(*selPID1)(HParticleCand *), Bool_t(*selPID2)(HParticleCand *))
map< Int_t, vector< HParticleCand * > > mTofHittoCand
all pair combinations freference x fothers
void setUseLeptons(Bool_t use)
vector< HParticlePair > & getPairsVector()
vector< HParticleCand * > fnofullrecoOthers
full reco cands (inner/outer MDC + META) inside others
Bool_t fuse_kIsLepton
default dilepton
vector< UInt_t > fCaseVec
counter array for cases
void setUserFilter(Bool_t(*userfilter)(HParticleCand *))