FEDRA emulsion software from the OPERA Collaboration
EdbLinking.h
Go to the documentation of this file.
1#ifndef ROOT_EdbLinking
2#define ROOT_EdbLinking
3
4#include "EdbScanCond.h"
5#include "TEnv.h"
6#include "EdbAlignmentV.h"
7#include "EdbSEQ.h"
8
11{
12 public:
13
20
22 float eCHI2Pmax;
23 float eChi2Acorr;
26 float eShr0;
27 float eDShr;
28 float eBinOK;
30
31 int eNshr[2];
33
34 TObjArray eSegCouples;
35
37
39
41 int remove;
42 float dr;
43 float dt;
46
47 public:
48 EdbLinking();
49 virtual ~EdbLinking(){}
50
51 bool VerifyShrinkageCorr( int side );
52 void GetPar(TEnv &env);
53 void GetDoubletsPar(TEnv &env);
54 void GetPreselectionPar(EdbSEQ &seq, TEnv &env);
55 void SaveCouplesTree(const char *file=0);
56 void FullLinking(TObjArray &p1, TObjArray &p2);
57 void FullLinking(EdbPattern &p1, EdbPattern &p2);
58 void FillCombinationsAtMeanZ(TObjArray &p1, TObjArray &p2);
59 void CorrectShrinkage(TObjArray &p1, TObjArray &p2, float dshr);
60 void CorrectShrinkage(EdbPattern &p1, EdbPattern &p2, float dshr);
61 void CorrectShrinkage(float dshr);
62 void Link(EdbPattern &p1, EdbPattern &p2, EdbLayer &l1, EdbLayer &l2, TEnv &env, Double_t area1=-1, Double_t area2=-1);
63 void CorrectAngles(TObjArray &p1, TObjArray &p2);
65
67 void RankCouples( TObjArray &arr1,TObjArray &arr2 );
68
69 void DoubletsFilterOut(TObjArray &p1, TObjArray &p2, bool fillhist=0);
70 void DumpDoubletsTree(EdbAlignmentV &adup, const char *name);
71
72 void CloneCouplesTree( const char *ifile, const char *ofile, EdbAffine2D *aff=0, TCut *cut=0 );
73
74 void SetApplyCorr(bool corr) {eCorr[0].eApplyCorr=corr;eCorr[1].eApplyCorr=corr;eCorr[2].eApplyCorr=corr;}
75
76 void ProduceReport();
77
78 ClassDef(EdbLinking,1) // microtracks linking in one plate
79};
80#endif /* ROOT_EdbLinking */
Definition: EdbAffine.h:17
universal basic alignment class
Definition: EdbAlignmentV.h:13
EdbSegCorr eCorr[2]
corrections for side 1 and 2 (v[7]) - the result of the alignment
Definition: EdbAlignmentV.h:24
fast 2-dim histogram class (used as a basis for EdbCell2)
Definition: EdbCell2.h:19
Definition: EdbLayer.h:39
microtracks linking in one plate
Definition: EdbLinking.h:11
float eNsigmaEQlnk
equalization cut for linking
Definition: EdbLinking.h:25
void DoubletsFilterOut(TObjArray &p1, TObjArray &p2, bool fillhist=0)
Definition: EdbLinking.cxx:594
int eNcorrMin
min number of segments for the correction calculation
Definition: EdbLinking.h:32
float eNsigmaEQshr
equalization cut for shrinkage
Definition: EdbLinking.h:24
void GetDoubletsPar(TEnv &env)
Definition: EdbLinking.cxx:69
bool eDoSaveCouples
Definition: EdbLinking.h:14
float eChi2Acorr
acceptance to save into couples tree
Definition: EdbLinking.h:23
EdbLayer eL1
Definition: EdbLinking.h:36
void RankCouples(TObjArray &arr1, TObjArray &arr2)
Definition: EdbLinking.cxx:438
bool eDoFullLinking
Definition: EdbLinking.h:17
TObjArray eSegCouples
segment couples objects to fill couples format tree
Definition: EdbLinking.h:34
EdbScanCond eCond
scanning conditions for couples ranking
Definition: EdbLinking.h:29
void Link(EdbPattern &p1, EdbPattern &p2, EdbLayer &l1, EdbLayer &l2, TEnv &env, Double_t area1=-1, Double_t area2=-1)
Definition: EdbLinking.cxx:96
bool eDoCorrectShrinkage
Definition: EdbLinking.h:16
float eDTfull
acceptance for the full linking
Definition: EdbLinking.h:21
Double_t EstimatePatternArea(EdbPattern &p)
Definition: EdbLinking.cxx:90
void GetPreselectionPar(EdbSEQ &seq, TEnv &env)
Definition: EdbLinking.cxx:78
virtual ~EdbLinking()
Definition: EdbLinking.h:49
void GetPar(TEnv &env)
Definition: EdbLinking.cxx:34
void SetApplyCorr(bool corr)
Definition: EdbLinking.h:74
void CorrectAngles(TObjArray &p1, TObjArray &p2)
Definition: EdbLinking.cxx:233
float eCHI2Pmax
acceptance to save into couples tree
Definition: EdbLinking.h:22
void WriteShrinkagePlots()
Definition: EdbLinking.cxx:418
int eCPRankingAlg
couples ranking algorithm (0-default, 1-likelihood used)
Definition: EdbLinking.h:19
bool eDoDumpDoubletsTree
Definition: EdbLinking.h:18
EdbH2 eHdxyShr[2]
service histograms used for the shrinkage correction
Definition: EdbLinking.h:38
bool eDoCorrectAngles
Definition: EdbLinking.h:15
bool VerifyShrinkageCorr(int side)
Definition: EdbLinking.cxx:202
void DumpDoubletsTree(EdbAlignmentV &adup, const char *name)
Definition: EdbLinking.cxx:624
void ProduceReport()
Definition: EdbLinking.cxx:505
EdbLayer eL2
layers with the geometry and corrections
Definition: EdbLinking.h:36
void FillCombinationsAtMeanZ(TObjArray &p1, TObjArray &p2)
Definition: EdbLinking.cxx:221
void CorrectShrinkage(TObjArray &p1, TObjArray &p2, float dshr)
Definition: EdbLinking.cxx:367
struct EdbLinking::RemoveDoublets eRemoveDoublets
int eNshr[2]
number of coins found for shrinkage correction
Definition: EdbLinking.h:31
float eBinOK
mean cell occupancy
Definition: EdbLinking.h:28
void CloneCouplesTree(const char *ifile, const char *ofile, EdbAffine2D *aff=0, TCut *cut=0)
Definition: EdbLinking.cxx:640
EdbLinking()
Definition: EdbLinking.cxx:29
void SaveCouplesTree(const char *file=0)
Definition: EdbLinking.cxx:296
float eDShr
range for the shrinkage correction
Definition: EdbLinking.h:27
void FullLinking(TObjArray &p1, TObjArray &p2)
Definition: EdbLinking.cxx:265
float eDRfull
Definition: EdbLinking.h:21
float eShr0
starting value for shrinkage correction search
Definition: EdbLinking.h:26
Definition: EdbPattern.h:273
Definition: EdbSEQ.h:12
Definition: EdbScanCond.h:10
Bool_t eApplyCorr
do correction
Definition: EdbSegCorr.h:13
TCut cut
Definition: check_shower.C:6
TFile * file
Definition: write_pvr.C:3
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24
Definition: EdbLinking.h:40
int remove
Definition: EdbLinking.h:41
float dt
Definition: EdbLinking.h:43
int checkview
Definition: EdbLinking.h:44
float dr
Definition: EdbLinking.h:42
p
Definition: testBGReduction_AllMethods.C:8