FEDRA emulsion software from the OPERA Collaboration
EGraphRecProc.cxx File Reference
#include "EGraphRec.h"
#include "EdbScanProc.h"
#include "EdbVertex.h"
#include <TThread.h>
#include <iostream>
Include dependency graph for EGraphRecProc.cxx:

Functions

 ClassImp (EGraphRecProc) void *ThSBProcess(void *ptr)
 
voidThSBCheckProcess (void *ptr)
 

Function Documentation

◆ ClassImp()

ClassImp ( EGraphRecProc  )
26{
27 EGraphRec *egraph = (EGraphRec*) ptr;
28 EdbScanProc *scanProc = egraph->GetScanProc(); // Scan Proc
29 EdbScanSet *scanSet = egraph->GetScanSet(); // Scan Set
30 EdbPattern *predTracks = egraph->GetPredTracks(); // predicted tracks
31 EdbPVRec *foundTracks = egraph->GetFoundTracks(); // found tracks
32 ProcBrick_t brickToProc = egraph->GetBrickToProc(); // Brick to process
33 ProcId_t procId = egraph->GetProcId();
34
35 Int_t pID_IC[4] = {brickToProc.brickId, 0, brickToProc.ver,
36 procId.interCalib};
37 Int_t pID_PS[4] = {brickToProc.brickId, 0, brickToProc.ver,
38 procId.predScan};
39 Int_t pID_IC_prev[4] = {brickToProc.brickId, 0, brickToProc.ver,
40 procId.interCalib};
41 Int_t pID_PS_prev[4] = {brickToProc.brickId, 0, brickToProc.ver,
42 procId.predScan};
43
44 Int_t firstPlate = brickToProc.firstPlate;
45 Int_t lastPlate = brickToProc.lastPlate;
46 Int_t step=(lastPlate>=firstPlate) ? brickToProc.step : -1*brickToProc.step;
47 Int_t nFailurePlates = 0;
48
49 predTracks->Reset();
50
51 for (Int_t plate = firstPlate; plate != lastPlate + step; plate += step) {
52
53 pID_PS[1] = pID_IC[1] = plate;
54 pID_PS_prev[1] = pID_IC_prev[1] = plate - step;
55
56 // get prediction for the current plate
57
58 // scanProc->ReadPred(*predTrack, pID_PS);
59
60 // scanning intercalibration zone
61
62 // Checking directory (if there is no scanning)
63
64 TString str;
65 scanProc->MakeFileName(str, pID_PS, "root");
66
67 FileStat_t buf;
68
69 if (gSystem->GetPathInfo(gSystem->DirName(str), buf)) {
70 cout << "ERROR! Directory " << gSystem->DirName(str)
71 << " does not exist.\n";
72 cout << "Please run scanning procedure for the plate number "
73 << pID_PS[1] << endl;
74 continue;
75 }
76
77 // linking process
78
79 if (egraph->IsSBToLink() && !scanProc->LinkRunAll(pID_IC)) return 0;
80
81 // alignment process between previous and current plate
82
83 if (egraph->IsSBToAlgn() && plate != firstPlate) {
84 if (!scanProc->SetAFFDZ(pID_IC_prev, pID_IC, 1300.*step)) return 0;
85 if (scanProc->AlignAll(pID_IC_prev, pID_IC, 1, 4, "-z") < 0) return 0;
86 }
87
88 // Scan back mode. Search predictions.
89
90 scanProc->CopyPar(pID_IC, pID_PS);
91
92 if (plate != firstPlate) {
93 scanProc->CopyAFFPar(pID_IC_prev, pID_IC, pID_PS_prev, pID_PS);
94 scanProc->ProjectFound(pID_PS_prev, pID_PS);
95 }
96
97 // adding scanned plate to the brick
98 // read affine and assemble "EdbBrickP" object
99
100 if (scanSet->AddID(new EdbID(pID_PS), step)) {
101 if (scanProc->AssembleScanSet(*scanSet) < 1) return 0;
102 scanSet->SetAsReferencePlate(firstPlate);
103 }
104
105 // after scann
106
107 foundTracks->ResetTracks();
108 scanProc->ReadFoundTracks(*scanSet, *foundTracks);
109
110 // come back to previous scanned plate + 1 if there are no tracks were
111 // found
112
113 EdbPattern pat;
114 if (!scanProc->ReadFound(pat, pID_PS)) {
115 if (Abs(step) != 1) {
116 plate += -step + (step > 0 ? 1 : -1);
117 if (!nFailurePlates) plate += step;
118 }
119 nFailurePlates++;
120 cout << "Track does not found in the plate " << plate
121 << ". starting scann plate " << plate + step << endl;
122 }
123 else nFailurePlates = 0;
124
125 egraph->DrawEvent();
126
127 if (nFailurePlates == 3) {
128 cout << "There are no any tracks during " << nFailurePlates
129 << "plates. Stop scanning" << endl;
130 return 0;
131 }
132 }
133
134 return 0;
135}
Definition: EGraphRec.h:24
EdbScanSet * GetScanSet() const
Definition: EGraphRec.h:57
ProcBrick_t GetBrickToProc() const
Definition: EGraphRec.h:54
EdbPVRec * GetFoundTracks() const
Definition: EGraphRec.h:59
void DrawEvent()
Definition: EGraphRec.cxx:142
Bool_t IsSBToLink() const
Definition: EGraphRec.h:52
EdbScanProc * GetScanProc() const
Definition: EGraphRec.h:56
ProcId_t GetProcId() const
Definition: EGraphRec.h:55
EdbPattern * GetPredTracks() const
Definition: EGraphRec.h:58
Bool_t IsSBToAlgn() const
Definition: EGraphRec.h:53
Definition: EdbID.h:7
Definition: EdbPVRec.h:148
void ResetTracks()
Definition: EdbPVRec.cxx:963
Definition: EdbPattern.h:273
void Reset()
Definition: EdbPattern.cxx:1509
scanned data processing
Definition: EdbScanProc.h:12
int CopyPar(EdbID id1, EdbID id2, bool overwrite=true)
Definition: EdbScanProc.h:48
int ReadFoundTracks(EdbScanSet &ss, EdbPVRec &ali, int flag=-1)
Definition: EdbScanProc.cxx:221
bool SetAFFDZ(int id1[4], int id2[4], float dz)
Definition: EdbScanProc.cxx:1054
int AlignAll(int id1[4], int id2[4], int npre=1, int nfull=3, const char *opt="-z")
Definition: EdbScanProc.cxx:2595
int CopyAFFPar(int id1c[4], int id2c[4], int id1p[4], int id2p[4], bool overwrite=true)
Definition: EdbScanProc.cxx:976
int LinkRunAll(int id[4], int npre=3, int nfull=1, int correct_ang=1)
Definition: EdbScanProc.cxx:1065
int ReadFound(EdbPattern &pred, int id[4], int flag=-1)
Definition: EdbScanProc.h:70
void MakeFileName(TString &s, int id[4], const char *suffix, bool inplate=true)
Definition: EdbScanProc.cxx:1819
int AssembleScanSet(EdbScanSet &ss)
Definition: EdbScanProc.cxx:135
bool ProjectFound(int id1[4], int id2[4])
Definition: EdbScanProc.cxx:1096
Definition: EdbScanSet.h:11
Bool_t SetAsReferencePlate(Int_t pid)
Definition: EdbScanSet.cxx:191
Bool_t AddID(EdbID *id, Int_t step)
Definition: EdbScanSet.cxx:476
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
Definition: EGraphRecProc.h:15
Int_t brickId
Definition: EGraphRecProc.h:16
Int_t step
Definition: EGraphRecProc.h:20
Int_t firstPlate
Definition: EGraphRecProc.h:17
Int_t ver
Definition: EGraphRecProc.h:19
Int_t lastPlate
Definition: EGraphRecProc.h:18
Definition: EGraphRecProc.h:25
Int_t predScan
Definition: EGraphRecProc.h:28
Int_t interCalib
Definition: EGraphRecProc.h:26

◆ ThSBCheckProcess()

void * ThSBCheckProcess ( void ptr)

Set enable "Execute Event" button after finishing process

140{
142
143 EGraphRec *egraph = (EGraphRec*) ptr;
144 TThread *thread = egraph->GetThSBProcess();
145 egraph->GetThSBProcess()->Join(); // Join ThSBProcess
146 egraph->GetButtonSBStart()->SetEnabled(kTRUE); // Enable button
147 SafeDelete(thread);
148 return 0;
149}
TThread * GetThSBProcess() const
Definition: EGraphRec.h:60
TGTextButton * GetButtonSBStart()
Definition: EGraphRec.h:61