FEDRA emulsion software from the OPERA Collaboration
EdbShowAlg_NN Class Reference

#include <EdbShowAlg_NN.h>

Inheritance diagram for EdbShowAlg_NN:
Collaboration diagram for EdbShowAlg_NN:

Public Member Functions

 ClassDef (EdbShowAlg_NN, 1)
 
TMultiLayerPerceptron * Create_NN_ALG_MLP (TTree *inputtree, Int_t inputneurons)
 
void CreateANNTree ()
 
 EdbShowAlg_NN ()
 
void Execute ()
 
void Finalize ()
 
Int_t GetMeansBeforeAndAfter (Float_t &mean_dT, Float_t &mean_dR, EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
 
Int_t GetMinsBeforeAndAfter (Float_t &min_dT, Float_t &min_dR, EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
 
Int_t GetNSegBeforeAndAfter (EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
 
TString GetWeightFileString ()
 
void Init ()
 
void Initialize ()
 
void LoadANNWeights ()
 
void LoadANNWeights (TMultiLayerPerceptron *TMlpANN, TString WeightFileString)
 
void Print ()
 
void SetANNWeightString ()
 
void SetWeightFileString (TString WeightFileString)
 
virtual ~EdbShowAlg_NN ()
 
- Public Member Functions inherited from EdbShowAlg
void AddRecoShowerArray (EdbShowerP *shower)
 
 ClassDef (EdbShowAlg, 1)
 
Double_t DeltaR_NoPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithoutPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaTheta (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaComponentwise (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaSingleAngles (EdbSegP *s1, EdbSegP *s2)
 
 EdbShowAlg ()
 
 EdbShowAlg (TString AlgName, Int_t AlgValue)
 
virtual void Execute ()
 
virtual void Finalize ()
 
TString GetAlgName () const
 
Int_t GetAlgValue () const
 
Double_t GetMinimumDist (EdbSegP *seg1, EdbSegP *seg2)
 
TObjArray * GetRecoShowerArray () const
 
Int_t GetRecoShowerArrayN () const
 
EdbShowerPGetShower (Int_t i) const
 
Double_t GetSpatialDist (EdbSegP *s1, EdbSegP *s2)
 
void Help ()
 
virtual void Initialize ()
 
Bool_t IsInConeTube (EdbSegP *sTest, EdbSegP *sStart, Double_t CylinderRadius, Double_t ConeAngle)
 
void Print ()
 
void PrintAll ()
 
void PrintMore ()
 
void PrintParameters ()
 
void PrintParametersShort ()
 
void SetActualAlgParameterset (Int_t ActualAlgParametersetNr)
 
void SetEdbPVRec (EdbPVRec *Ali)
 
void SetEdbPVRecPIDNumbers (Int_t FirstPlate_eAliPID, Int_t LastPlate_eAliPID, Int_t MiddlePlate_eAliPID, Int_t NumberPlate_eAliPID)
 
void SetInBTArray (TObjArray *InBTArray)
 
void SetParameter (Int_t parNr, Float_t par)
 
void SetParameters (Float_t *par)
 
void SetRecoShowerArray (TObjArray *RecoShowerArray)
 
void SetRecoShowerArrayN (Int_t RecoShowerArrayN)
 
void SetUseAliSub (Bool_t UseAliSub)
 
void Transform_eAli (EdbSegP *InitiatorBT, Float_t ExtractSize)
 
void UpdateShowerIDs ()
 
void UpdateShowerMetaData ()
 
virtual ~EdbShowAlg ()
 

Private Attributes

Int_t eANN_Inputtype
 
Float_t eANN_var_dR_InBT_To_TestBT
 
Float_t eANN_var_dR_NextBT_To_TestBT
 
Float_t eANN_var_dR_TestBT_To_InBT
 
Float_t eANN_var_dT_InBT_To_TestBT
 
Float_t eANN_var_dT_NextBT_To_TestBT
 
Float_t eANN_var_InBT_To_TestBT
 
Float_t eANN_var_mean_dR_2after
 
Float_t eANN_var_mean_dR_2before
 
Float_t eANN_var_mean_dR_after
 
Float_t eANN_var_mean_dR_before
 
Float_t eANN_var_mean_dR_same
 
Float_t eANN_var_mean_dT_2after
 
Float_t eANN_var_mean_dT_2before
 
Float_t eANN_var_mean_dT_after
 
Float_t eANN_var_mean_dT_before
 
Float_t eANN_var_mean_dT_same
 
Float_t eANN_var_min_dR_2after
 
Float_t eANN_var_min_dR_2before
 
Float_t eANN_var_min_dR_after
 
Float_t eANN_var_min_dR_before
 
Float_t eANN_var_min_dR_same
 
Float_t eANN_var_min_dT_2after
 
Float_t eANN_var_min_dT_2before
 
Float_t eANN_var_min_dT_after
 
Float_t eANN_var_min_dT_before
 
Float_t eANN_var_min_dT_same
 
Int_t eANN_var_nseg_1after
 
Int_t eANN_var_nseg_1before
 
Int_t eANN_var_nseg_2after
 
Int_t eANN_var_nseg_2before
 
Int_t eANN_var_nseg_3after
 
Int_t eANN_var_nseg_3before
 
Int_t eANN_var_nseg_same
 
Float_t eANN_var_SpatialDist_TestBT_To_InBT
 
Float_t eANN_var_zDiff_TestBT_To_InBT
 
Float_t eANN_var_zDist_TestBT_To_InBT
 
TTree * eANNTree
 
TMultiLayerPerceptron * eTMlpANN
 
TString eWeightFileLayoutString
 
TString eWeightFileString
 

Additional Inherited Members

- Protected Member Functions inherited from EdbShowAlg
void Set0 ()
 
- Protected Attributes inherited from EdbShowAlg
Int_t eActualAlgParametersetNr
 
TString eAlgName
 
Int_t eAlgValue
 
EdbPVReceAli
 
EdbPVReceAli_Sub
 
Int_t eAli_SubNpat
 
Int_t eAliNpat
 
Int_t eFirstPlate_eAliPID
 
TObjArray * eInBTArray
 
Int_t eInBTArrayN
 
Int_t eLastPlate_eAliPID
 
Int_t eMiddlePlate_eAliPID
 
Int_t eNumberPlate_eAliPID
 
Int_t eParaN
 
TString eParaString [10]
 
Float_t eParaValue [10]
 
EdbShowerPeRecoShower
 
TObjArray * eRecoShowerArray
 
Int_t eRecoShowerArrayN
 
Int_t eUseAliSub
 

Constructor & Destructor Documentation

◆ EdbShowAlg_NN()

EdbShowAlg_NN::EdbShowAlg_NN ( )

◆ ~EdbShowAlg_NN()

EdbShowAlg_NN::~EdbShowAlg_NN ( )
virtual
41{
42 // Default Destructor
43 Log(2,"EdbShowAlg_NN::~EdbShowAlg_NN","Default Destructor");
44 if (eANNTree) {
45 delete eANNTree;
46 eANNTree=0;
47 }
48}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
TTree * eANNTree
Definition: EdbShowAlg_NN.h:52

Member Function Documentation

◆ ClassDef()

EdbShowAlg_NN::ClassDef ( EdbShowAlg_NN  ,
 
)

◆ Create_NN_ALG_MLP()

TMultiLayerPerceptron * EdbShowAlg_NN::Create_NN_ALG_MLP ( TTree *  inputtree,
Int_t  inputneurons 
)
140{
141 Log(2,"EdbShowAlg_NN::Create_NN_ALG_MLP","Creates the TTree for the data to be filled.");
142
143 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_NN::Create_NN_ALG_MLP() inputneurons= " << inputneurons << endl;
144
145 const char* layout="";
146 if (inputneurons==5) {
147 layout="eANN_var_dT_InBT_To_TestBT,eANN_var_dR_InBT_To_TestBT,eANN_var_dR_TestBT_To_InBT,eANN_var_zDiff_TestBT_To_InBT,eANN_var_SpatialDist_TestBT_To_InBT:6:5:type";
148 }
149 if (inputneurons==10) {
150 layout="@eANN_var_dT_InBT_To_TestBT,@eANN_var_dR_InBT_To_TestBT,@eANN_var_dR_TestBT_To_InBT,@eANN_var_zDiff_TestBT_To_InBT,eANN_var_SpatialDist_TestBT_To_InBT,eANN_var_nseg_2before,eANN_var_nseg_1before,eANN_var_nseg_same,eANN_var_nseg_1after,eANN_var_nseg_2after:11:10:type";
151 }
152 if (inputneurons==20) {
153 layout="@eANN_var_dT_InBT_To_TestBT,@eANN_var_dR_InBT_To_TestBT,@eANN_var_dR_TestBT_To_InBT,@eANN_var_zDiff_TestBT_To_InBT,eANN_var_SpatialDist_TestBT_To_InBT,eANN_var_nseg_2before,eANN_var_nseg_1before,eANN_var_nseg_same,eANN_var_nseg_1after,eANN_var_nseg_2after,eANN_var_mean_dT_2before,eANN_var_mean_dT_before,eANN_var_mean_dT_same,eANN_var_mean_dT_after,eANN_var_mean_dT_2after,eANN_var_mean_dR_2before,eANN_var_mean_dR_before,eANN_var_mean_dR_same,eANN_var_mean_dR_after,eANN_var_mean_dR_2after:21:20:type";
154 }
155 if (inputneurons==30) {
156 layout="@eANN_var_dT_InBT_To_TestBT,@eANN_var_dR_InBT_To_TestBT,@eANN_var_dR_TestBT_To_InBT,@eANN_var_zDiff_TestBT_To_InBT,eANN_var_SpatialDist_TestBT_To_InBT,eANN_var_nseg_2before,eANN_var_nseg_1before,eANN_var_nseg_same,eANN_var_nseg_1after,eANN_var_nseg_2after,eANN_var_mean_dT_2before,eANN_var_mean_dT_before,eANN_var_mean_dT_same,eANN_var_mean_dT_after,eANN_var_mean_dT_2after,eANN_var_mean_dR_2before,eANN_var_mean_dR_before,eANN_var_mean_dR_same,eANN_var_mean_dR_after,eANN_var_mean_dR_2after,eANN_var_min_dT_2before,eANN_var_min_dT_before,eANN_var_min_dT_same,eANN_var_min_dT_after,eANN_var_min_dT_2after,eANN_var_min_dR_2before,eANN_var_min_dR_before,eANN_var_min_dR_same,eANN_var_min_dR_after,eANN_var_min_dR_2after:31:30:type";
157 }
158
159 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_NN::Create_NN_ALG_MLP() layout: " << layout << endl;
160
161 // Create the network:
162 TMultiLayerPerceptron* TMlpANN = new TMultiLayerPerceptron(layout,simu);
163
164 if (gEDBDEBUGLEVEL>2) {
165 cout << "EdbShowAlg_NN::Create_NN_ALG_MLP() GetStructure: " << endl;
166 cout << TMlpANN->GetStructure() << endl;
167 }
168 Log(2,"EdbShowAlg_NN::Create_NN_ALG_MLP","Creates the TTree for the data to be filled...done.");
169 return TMlpANN;
170}
TMultiLayerPerceptron * TMlpANN
Definition: ShowRec.h:340
gEDBDEBUGLEVEL
Definition: energy.C:7
TTree * simu
Definition: testBGReduction_By_ANN.C:12

◆ CreateANNTree()

void EdbShowAlg_NN::CreateANNTree ( )
88{
89 Log(2,"EdbShowAlg_NN::CreateANNTree","CreateANNTree()");
90
91 if (!eANNTree) eANNTree = new TTree("EdbShowAlg_NN_eANNTree", "EdbShowAlg_NN_eANNTree");
92
93 // Variables and things important for neural Network:
94 eANNTree->Branch("dT_InBT_To_TestBT", &eANN_var_InBT_To_TestBT, "dT_InBT_To_TestBT/F");
95 eANNTree->Branch("dR_InBT_To_TestBT", &eANN_var_dR_InBT_To_TestBT, "dR_InBT_To_TestBT/F");
96 eANNTree->Branch("dR_TestBT_To_InBT", &eANN_var_dR_TestBT_To_InBT, "dR_TestBT_To_InBT/F");
97 eANNTree->Branch("zDiff_TestBT_To_InBT", &eANN_var_zDiff_TestBT_To_InBT, "zDiff_TestBT_To_InBT/F");
98 eANNTree->Branch("SpatialDist_TestBT_To_InBT", &eANN_var_SpatialDist_TestBT_To_InBT, "SpatialDist_TestBT_To_InBT/F");
99 //---------
100 eANNTree->Branch("nseg_1before", &eANN_var_nseg_1before, "nseg_1before/I");
101 eANNTree->Branch("nseg_2before", &eANN_var_nseg_2before, "nseg_2before/I");
102 eANNTree->Branch("nseg_same", &eANN_var_nseg_same, "nseg_same/I");
103 eANNTree->Branch("nseg_1after", &eANN_var_nseg_1after, "nseg_1after/I");
104 eANNTree->Branch("nseg_2after", &eANN_var_nseg_2after, "nseg_2after/I");
105 //---------
106 eANNTree->Branch("mean_dT_2before", &eANN_var_mean_dT_2before, "mean_dT_2before/F");
107 eANNTree->Branch("mean_dT_before", &eANN_var_mean_dT_before, "mean_dT_before/F");
108 eANNTree->Branch("mean_dT_same", &eANN_var_mean_dT_same, "mean_dT_same/F");
109 eANNTree->Branch("mean_dT_after", &eANN_var_mean_dT_after, "mean_dT_after/F");
110 eANNTree->Branch("mean_dT_2after", &eANN_var_mean_dT_2after, "mean_dT_2after/F");
111 //---------
112 eANNTree->Branch("mean_dR_2before", &eANN_var_mean_dR_2before, "mean_dR_2before/F");
113 eANNTree->Branch("mean_dR_before", &eANN_var_mean_dR_before, "mean_dR_before/F");
114 eANNTree->Branch("mean_dR_same", &eANN_var_mean_dR_same, "mean_dR_same/F");
115 eANNTree->Branch("mean_dR_after", &eANN_var_mean_dR_after, "mean_dR_after/F");
116 eANNTree->Branch("mean_dR_2after", &eANN_var_mean_dR_2after, "mean_dR_2after/F");
117 //---------
118 eANNTree->Branch("min_dT_2before", &eANN_var_min_dT_2before, "min_dT_2before/F");
119 eANNTree->Branch("min_dT_before", &eANN_var_min_dT_before, "min_dT_before/F");
120 eANNTree->Branch("min_dT_same", &eANN_var_min_dT_same, "min_dT_same/F");
121 eANNTree->Branch("min_dT_after", &eANN_var_min_dT_after, "min_dT_after/F");
122 eANNTree->Branch("min_dT_2after", &eANN_var_min_dT_2after, "min_dT_2after/F");
123 //---------
124 eANNTree->Branch("min_dR_2before", &eANN_var_min_dR_2before, "min_dR_2before/F");
125 eANNTree->Branch("min_dR_before", &eANN_var_min_dR_before, "min_dR_before/F");
126 eANNTree->Branch("min_dR_same", &eANN_var_min_dR_same, "min_dR_same/F");
127 eANNTree->Branch("min_dR_after", &eANN_var_min_dR_after, "min_dR_after/F");
128 eANNTree->Branch("min_dR_2after", &eANN_var_min_dR_2after, "min_dR_2after/F");
129 //---------
130 eANNTree->Branch("type", &eANN_Inputtype, "type/I");
131 //---------
132 Log(2,"EdbShowAlg_NN::CreateANNTree","CreateANNTree()...done.");
133 return;
134}
Int_t eANN_var_nseg_2before
Definition: EdbShowAlg_NN.h:87
Float_t eANN_var_dR_TestBT_To_InBT
Definition: EdbShowAlg_NN.h:60
Float_t eANN_var_min_dR_after
Definition: EdbShowAlg_NN.h:83
Float_t eANN_var_mean_dR_same
Definition: EdbShowAlg_NN.h:71
Float_t eANN_var_mean_dR_before
Definition: EdbShowAlg_NN.h:69
Float_t eANN_var_min_dR_before
Definition: EdbShowAlg_NN.h:79
Float_t eANN_var_min_dR_2after
Definition: EdbShowAlg_NN.h:85
Float_t eANN_var_min_dR_same
Definition: EdbShowAlg_NN.h:81
Float_t eANN_var_mean_dR_after
Definition: EdbShowAlg_NN.h:73
Float_t eANN_var_InBT_To_TestBT
Definition: EdbShowAlg_NN.h:57
Float_t eANN_var_mean_dT_same
Definition: EdbShowAlg_NN.h:70
Float_t eANN_var_min_dT_before
Definition: EdbShowAlg_NN.h:78
Float_t eANN_var_min_dT_after
Definition: EdbShowAlg_NN.h:82
Float_t eANN_var_mean_dT_before
Definition: EdbShowAlg_NN.h:68
Float_t eANN_var_SpatialDist_TestBT_To_InBT
Definition: EdbShowAlg_NN.h:62
Float_t eANN_var_min_dT_2before
Definition: EdbShowAlg_NN.h:76
Int_t eANN_var_nseg_1after
Definition: EdbShowAlg_NN.h:89
Float_t eANN_var_zDiff_TestBT_To_InBT
Definition: EdbShowAlg_NN.h:63
Float_t eANN_var_mean_dT_after
Definition: EdbShowAlg_NN.h:72
Float_t eANN_var_mean_dT_2before
Definition: EdbShowAlg_NN.h:66
Int_t eANN_var_nseg_same
Definition: EdbShowAlg_NN.h:92
Float_t eANN_var_dR_InBT_To_TestBT
Definition: EdbShowAlg_NN.h:59
Int_t eANN_var_nseg_2after
Definition: EdbShowAlg_NN.h:90
Float_t eANN_var_mean_dR_2after
Definition: EdbShowAlg_NN.h:75
Int_t eANN_var_nseg_1before
Definition: EdbShowAlg_NN.h:86
Float_t eANN_var_min_dT_2after
Definition: EdbShowAlg_NN.h:84
Float_t eANN_var_min_dT_same
Definition: EdbShowAlg_NN.h:80
Int_t eANN_Inputtype
Definition: EdbShowAlg_NN.h:93
Float_t eANN_var_mean_dT_2after
Definition: EdbShowAlg_NN.h:74
Float_t eANN_var_mean_dR_2before
Definition: EdbShowAlg_NN.h:67
Float_t eANN_var_min_dR_2before
Definition: EdbShowAlg_NN.h:77

◆ Execute()

void EdbShowAlg_NN::Execute ( )
virtual

Reimplemented from EdbShowAlg.

218{
219 Log(2,"EdbShowAlg_NN::Execute","Execute() DOING MAIN SHOWER RECONSTRUCTION HERE");
220
221 EdbSegP* InBT;
222 EdbSegP* Segment;
223 EdbSegP* seg;
224 EdbShowerP* RecoShower;
225
226 Bool_t StillToLoop=kTRUE;
227 Int_t ActualPID;
228 Int_t newActualPID;
229 Int_t STEP=-1;
230 Int_t NLoopedPattern=0;
232 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- STEP for patternloop direction = " << STEP << endl;
233
234
235 Double_t params[30]; // Used for ANN Evaluation
236
237 //--- Loop over InBTs:
238
239 // Since eInBTArray is filled in ascending ordering by zpositon
240 // We use the descending loop to begin with BT with lowest z first.
241 for (Int_t i=eInBTArrayN-1; i>=0; --i) {
242
243 // CounterOutPut
244 if (gEDBDEBUGLEVEL==2) if ((i%1)==0) cout << eInBTArrayN <<" InBT in total, still to do:"<<Form("%4d",i)<< "\r\r\r\r"<<flush;
245
246 //-----------------------------------
247 // 0a) Reset characteristic variables:
248 //-----------------------------------
278
279 //-----------------------------------
280 // 1) Make eAli with cut parameters:
281 //-----------------------------------
282
283 // Create new EdbShowerP Object for storage;
284 // See EdbShowerP why I have to call the Constructor as "unique" ideficable value
285 //RecoShower = new EdbShowerP(i,eAlgValue);
286 RecoShower = new EdbShowerP(i,eAlgValue,-1);
287
288 // Get InitiatorBT from eInBTArray
289 InBT=(EdbSegP*)eInBTArray->At(i);
290 if (gEDBDEBUGLEVEL>2) InBT->PrintNice();
291
292 // Add InBT to RecoShower:
293 // This has to be done, since by definition the first BT in the RecoShower is the InBT.
294 // Otherwise, also the definition of shower axis and transversal profiles is wrong!
295 RecoShower -> AddSegment(InBT);
296 cout << "Segment (InBT) " << InBT << " was added to RecoShower." << endl;
297
298 // Transform (make size smaller, extract only events having same MC) the eAli object:
299 // See in Execute_CA for description.
300 Transform_eAli(InBT,1400);
301
302 // Add InBT to RecoShower // obsolete, since we loop also over the plate containing the InBT
303 //RecoShower --> AddSegment(InBT);
304
305 //-----------------------------------
306 // 2) Loop over (whole) eAli, check BT for Cuts
307 //-----------------------------------
308 ActualPID= InBT->PID() ;
309 newActualPID= InBT->PID() ;
310
311 while (StillToLoop) {
312 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- --- Doing patterloop " << ActualPID << endl;
313
314 // just to adapt to this nomenclature of ShowRec program:
315 int patterloop_cnt=ActualPID;
316
317
318 for (Int_t btloop_cnt=0; btloop_cnt<eAli->GetPattern(ActualPID)->GetN(); ++btloop_cnt) {
319
320 Segment = (EdbSegP*)eAli->GetPattern(ActualPID)->GetSegment(btloop_cnt);
321
322 // just to adapt to this nomenclature of ShowRec program
323 seg=Segment;
324
325 if (gEDBDEBUGLEVEL>4) Segment->PrintNice();
326
327 // Reset characteristic variables:
357
358 // Now calculate NN Inputvariables: --------------------
359 // 1) zDiff_TestBT_To_InBT and SpatialDist_TestBT_To_InBT
360 eANN_var_zDiff_TestBT_To_InBT=seg->Z()-InBT->Z();
362 // 2) dT_InBT_To_TestBT
364 // 3) dR_InBT_To_TestBT (propagation order matters)
365 // In calculation it makes a difference if InBT is extrapolated to Z-pos of seg or vice versa.
368
369 // 4) nseg_1before,nseg_2before,nseg_1after,nseg_2after:
370 if (eParaValue[0]>=10) {
371 eANN_var_nseg_2before = GetNSegBeforeAndAfter(eAli,patterloop_cnt,seg,2,-1);
372 eANN_var_nseg_1before = GetNSegBeforeAndAfter(eAli,patterloop_cnt,seg,1,-1);
373 eANN_var_nseg_same = GetNSegBeforeAndAfter(eAli,patterloop_cnt,seg,0, 1);
374 eANN_var_nseg_1after = GetNSegBeforeAndAfter(eAli,patterloop_cnt,seg,1, 1);
375 eANN_var_nseg_2after = GetNSegBeforeAndAfter(eAli,patterloop_cnt,seg,2, 1);
376 }
377
378 // 5) mean_dT,dR nseg_1before,nseg_2before,nseg_1after,nseg_2after: mean of all dTheta and dR compliment segments:
379 if (eParaValue[0]>=20) {
385 }
386
387 // 6) nseg_1before,nseg_2before,nseg_1after,nseg_2after: mean of all dTheta and dR compliment segments:
388 if (eParaValue[0]>=30) {
394 }
395 // end of calculate NN Inputvariables: --------------------
396
397 // Fill ANN params Array:
428
429
430
431 // Now apply cut conditions: NN Neural Network Alg --------------------
432 Double_t value=0;
433 value=eTMlpANN->Evaluate(0, params);
434 if (gEDBDEBUGLEVEL>3) {
435 cout << "Inputvalues: ";
436 for (int hj=0; hj<5; hj++) cout << " " << params[hj];
437 cout << " , Value= " << value << " eParaValue= " << eParaValue[1] << endl;
438 }
439 if (value<eParaValue[1]) continue;
440 // end of cut conditions: NN Neural Network Alg --------------------
441
442 // If we arrive here, Basetrack Segment has passed criteria
443 // and is then added to the RecoShower:
444 // Check if its not the InBT which is already added:
445 if (Segment->X()==InBT->X()&&Segment->Y()==InBT->Y()) {
446 ; // is InBT, do nothing;
447 }
448 else {
449 RecoShower -> AddSegment(Segment);
450 }
451 cout << "Segment " << Segment << " was added to &RecoShower : " << &RecoShower << endl;
452 } // of btloop_cnt
453
454 //------------
455 newActualPID=ActualPID+STEP;
456 ++NLoopedPattern;
457
458 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- --- newActualPID= " << newActualPID << endl;
459 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- --- NLoopedPattern= " << NLoopedPattern << endl;
460 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- --- eNumberPlate_eAliPID= " << eNumberPlate_eAliPID << endl;
461 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_NN::Execute--- --- StillToLoop= " << StillToLoop << endl;
462
463 // This if holds in the case of STEP== +1
464 if (STEP==1) {
465 if (newActualPID>eLastPlate_eAliPID) StillToLoop=kFALSE;
466 if (newActualPID>eLastPlate_eAliPID) cout << "EdbShowAlg_NN::Execute--- ---Stopp Loop since: newActualPID>eLastPlate_eAliPID"<<endl;
467 }
468 // This if holds in the case of STEP== -1
469 if (STEP==-1) {
470 if (newActualPID<eLastPlate_eAliPID) StillToLoop=kFALSE;
471 if (newActualPID<eLastPlate_eAliPID) cout << "EdbShowAlg_NN::Execute--- ---Stopp Loop since: newActualPID<eLastPlate_eAliPID"<<endl;
472 }
473 // This if holds general, since eNumberPlate_eAliPID is not dependent of the structure of the gAli subject:
474 if (NLoopedPattern>eNumberPlate_eAliPID) StillToLoop=kFALSE;
475 if (NLoopedPattern>eNumberPlate_eAliPID) cout << "EdbShowAlg_NN::Execute--- ---Stopp Loop since: NLoopedPattern>eNumberPlate_eAliPID"<<endl;
476
477 ActualPID=newActualPID;
478 } // of // while (StillToLoop)
479
480 // Obligatory when Shower Reconstruction is finished!
481 RecoShower ->Update();
482 //RecoShower ->PrintBasics();
483
484
485 // Add Shower Object to Shower Reco Array.
486 // Not, if its empty:
487 // Not, if its containing only one BT:
488 if (RecoShower->N()>1) eRecoShowerArray->Add(RecoShower);
489
490 // Set back loop values:
491 StillToLoop=kTRUE;
492 NLoopedPattern=0;
493 } // of // for (Int_t i=eInBTArrayN-1; i>=0; --i) {
494
495
496 // Set new value for eRecoShowerArrayN (may now be < eInBTArrayN).
498
499 cout << "EdbShowAlg_NN::eRecoShowerArray() Entries: " << eRecoShowerArray->GetEntries() << endl;
500 Log(2,"EdbShowAlg_NN::Execute","Execute() DOING MAIN SHOWER RECONSTRUCTION HERE...done.");
501 return;
502}
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
void PrintNice() const
Definition: EdbSegP.cxx:392
Int_t PID() const
Definition: EdbSegP.h:148
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
TMultiLayerPerceptron * eTMlpANN
Definition: EdbShowAlg_NN.h:53
Int_t GetMinsBeforeAndAfter(Float_t &min_dT, Float_t &min_dR, EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
Definition: EdbShowAlg_NN.cxx:585
Int_t GetNSegBeforeAndAfter(EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
Definition: EdbShowAlg_NN.cxx:523
Float_t eANN_var_dT_InBT_To_TestBT
Definition: EdbShowAlg_NN.h:58
Int_t GetMeansBeforeAndAfter(Float_t &mean_dT, Float_t &mean_dR, EdbPVRec *local_gAli, Int_t patterloop_cnt, EdbSegP *seg, Int_t n_patterns, Int_t BeforeOrAfter)
Definition: EdbShowAlg_NN.cxx:672
void Transform_eAli(EdbSegP *InitiatorBT, Float_t ExtractSize)
Definition: EdbShowAlg.cxx:107
Float_t eParaValue[10]
Definition: EdbShowAlg.h:52
Double_t GetSpatialDist(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:463
void SetRecoShowerArrayN(Int_t RecoShowerArrayN)
Definition: EdbShowAlg.h:117
Double_t DeltaThetaSingleAngles(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:455
TObjArray * eRecoShowerArray
Definition: EdbShowAlg.h:70
Int_t eAlgValue
Definition: EdbShowAlg.h:50
Int_t eFirstPlate_eAliPID
Definition: EdbShowAlg.h:65
Int_t eLastPlate_eAliPID
Definition: EdbShowAlg.h:66
Int_t eInBTArrayN
Definition: EdbShowAlg.h:63
Double_t DeltaR_WithPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:411
TObjArray * eInBTArray
Definition: EdbShowAlg.h:62
EdbPVRec * eAli
Definition: EdbShowAlg.h:59
Int_t eNumberPlate_eAliPID
Definition: EdbShowAlg.h:68
Definition: EdbShowerP.h:28
Int_t N() const
Definition: EdbShowerP.h:412
void Update()
Definition: EdbShowerP.cxx:975
Double_t params[3]
Definition: testBGReduction_By_ANN.C:84

◆ Finalize()

void EdbShowAlg_NN::Finalize ( )
virtual

Reimplemented from EdbShowAlg.

516{
517 Log(2,"EdbShowAlg_NN::Finalize","Finalize(). Nothing done here yet.");
518 return;
519}

◆ GetMeansBeforeAndAfter()

Int_t EdbShowAlg_NN::GetMeansBeforeAndAfter ( Float_t &  mean_dT,
Float_t &  mean_dR,
EdbPVRec local_gAli,
Int_t  patterloop_cnt,
EdbSegP seg,
Int_t  n_patterns,
Int_t  BeforeOrAfter 
)
673{
674// cout << "GetMeansBeforeAndAfter( XX, " << patterloop_cnt << " , Seg, " << n_patterns << ", " << BeforeOrAfter << endl;
675 mean_dT=-1;
676 mean_dR=-1;
677
678 float Z_minus1=0;
679 float Z_normal=local_gAli->GetPattern(patterloop_cnt)->Z();
680 float Z_plus1=0;
681
682 int npat=local_gAli->Npatterns();
683 Bool_t edge_npat_upper=kFALSE;
684 Bool_t edge_npat_lower=kFALSE;
685 Int_t Factor=-1;
686
687 if (patterloop_cnt==npat-1) {
688 edge_npat_upper=kTRUE;
689 }
690 if (patterloop_cnt==0) {
691 edge_npat_lower=kTRUE;
692 }
693
694 if (!edge_npat_lower) {
695 Z_minus1=local_gAli->GetPattern(patterloop_cnt-1)->Z();
696 if (Z_minus1<Z_normal) Factor=1;
697 }
698 if (!edge_npat_upper) {
699 Z_plus1 =local_gAli->GetPattern(patterloop_cnt+1)->Z();
700 if (Z_plus1>Z_normal) Factor=1;
701 }
702 // New PID we want to have:
703 Int_t patterloop_test=patterloop_cnt+Factor*n_patterns*BeforeOrAfter;
704
705 // Does this plate exist? If not, return 0 directly:
706 if (patterloop_test>=npat || patterloop_test<0) {
707 //cout << "So NEW n_patterns would be " << patterloop_test << " BUT IT DOES NOT MATHC in our local_gAli sceheme, which means its not existing. RETURNING 0 " << endl;
708 return 0;
709 }
710
711 // Since we have checked now for bounds we can FindCompliments:
712 Int_t n_return=0;
713 TObjArray array;
714 array.Clear();
715 EdbPattern* TestPattern= (EdbPattern*)local_gAli->GetPattern(patterloop_test);
716 TestPattern -> FillCell(20,20,0.01,0.01);
717 n_return = TestPattern->FindCompliments(*seg,array,3,3);
718 //cout << " Found " << n_return << " compliments in 2,2 sigma area:" << endl;
719
720 if (n_return==0) return n_return;
721
722 //seg->PrintNice();
723 for (int i=0; i<n_return; i++) {
724 EdbSegP* s_of_array=(EdbSegP*)array.At(i);
725 if (i==0) {
726 mean_dT=0;
727 mean_dR=0;
728 }
729 //s_of_array->PrintNice();
730 mean_dT+=DeltaThetaSingleAngles(seg,s_of_array);
731 mean_dR+=DeltaR_NoPropagation(seg,s_of_array);
732 }
733 if (n_return>0) mean_dT=mean_dT/(Double_t)n_return;
734 if (n_return>0) mean_dR=mean_dR/(Double_t)n_return;
735
736 return n_return;
737}
Int_t npat
Definition: Xi2HatStartScript.C:33
Definition: EdbPattern.h:273
int FindCompliments(EdbSegP &s, TObjArray &arr, float nsig, float nsigt)
Definition: EdbPattern.cxx:1354
Int_t Npatterns() const
Definition: EdbPattern.h:366
Float_t Z() const
Definition: EdbPattern.h:84
Double_t DeltaR_NoPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:396

◆ GetMinsBeforeAndAfter()

Int_t EdbShowAlg_NN::GetMinsBeforeAndAfter ( Float_t &  min_dT,
Float_t &  min_dR,
EdbPVRec local_gAli,
Int_t  patterloop_cnt,
EdbSegP seg,
Int_t  n_patterns,
Int_t  BeforeOrAfter 
)
586{
587// cout << "GetMinsBeforeAndAfter( XX, " << patterloop_cnt << " , Seg, " << n_patterns << ", " << BeforeOrAfter << endl;
588 min_dT=-1;
589 min_dR=-1;
590
591 float Z_minus1=0;
592 float Z_normal=local_gAli->GetPattern(patterloop_cnt)->Z();
593 float Z_plus1=0;
594
595 int npat=local_gAli->Npatterns();
596 Bool_t edge_npat_upper=kFALSE;
597 Bool_t edge_npat_lower=kFALSE;
598 Int_t Factor=-1;
599
600 if (patterloop_cnt==npat-1) {
601 edge_npat_upper=kTRUE;
602 }
603 if (patterloop_cnt==0) {
604 edge_npat_lower=kTRUE;
605 }
606
607 if (!edge_npat_lower) {
608 Z_minus1=local_gAli->GetPattern(patterloop_cnt-1)->Z();
609 if (Z_minus1<Z_normal) Factor=1;
610 }
611 if (!edge_npat_upper) {
612 Z_plus1 =local_gAli->GetPattern(patterloop_cnt+1)->Z();
613 if (Z_plus1>Z_normal) Factor=1;
614 }
615 // New PID we want to have:
616 Int_t patterloop_test=patterloop_cnt+Factor*n_patterns*BeforeOrAfter;
617
618 // Does this plate exist? If not, return 0 directly:
619 if (patterloop_test>=npat || patterloop_test<0) {
620 //cout << "So NEW n_patterns would be " << patterloop_test << " BUT IT DOES NOT MATHC in our local_gAli sceheme, which means its not existing. RETURNING 0 " << endl;
621 return 0;
622 }
623
624 // Since we have checked now for bounds we can FindCompliments:
625 Int_t n_return=0;
626 TObjArray array;
627 array.Clear();
628 EdbPattern* TestPattern= (EdbPattern*)local_gAli->GetPattern(patterloop_test);
629 TestPattern -> FillCell(20,20,0.01,0.01);
630 n_return = TestPattern->FindCompliments(*seg,array,3,3);
631 //cout << " Found " << n_return << " compliments in 2,2 sigma area:" << endl;
632
633 if (n_return==0) return n_return;
634
635 if (n_return==1) {
636 EdbSegP* s_of_array=(EdbSegP*)array.At(0);
637 min_dT=DeltaThetaSingleAngles(seg,s_of_array);
638 min_dR=DeltaR_NoPropagation(seg,s_of_array);
639 }
640
641 Float_t tmp_min_dT=-1;
642 Float_t tmp_min_dR=-1;
643 Float_t tmp2_min_dT=-1;
644 Float_t tmp2_min_dR=-1;
645 Float_t angle;
646 Float_t dist;
647
648 if (n_return>1) {
649 for (int i=0; i<n_return; i++) {
650 EdbSegP* s_of_array=(EdbSegP*)array.At(i);
651 if (i==0) {
652 min_dT=999999;
653 min_dR=9999999;
654 }
655 angle=(Float_t)DeltaThetaSingleAngles(seg,s_of_array);
656 tmp_min_dT=min_dT;
657 tmp2_min_dT=TMath::Min(angle, tmp_min_dT);
658 min_dT=tmp2_min_dT;
659
660 dist=(Float_t)DeltaR_NoPropagation(seg,s_of_array);
661 tmp_min_dR=min_dR;
662 tmp2_min_dR=TMath::Min(dist, tmp_min_dR);
663 min_dR=tmp2_min_dR;
664 }
665 }
666 return n_return;
667}

◆ GetNSegBeforeAndAfter()

Int_t EdbShowAlg_NN::GetNSegBeforeAndAfter ( EdbPVRec local_gAli,
Int_t  patterloop_cnt,
EdbSegP seg,
Int_t  n_patterns,
Int_t  BeforeOrAfter 
)
524{
525// cout << "GetNSegBeforeAndAfter( XX, " << patterloop_cnt << " , Seg, " << n_patterns << ", " << BeforeOrAfter << endl;
526
527 float Z_minus1=0;
528 float Z_normal=local_gAli->GetPattern(patterloop_cnt)->Z();
529 float Z_plus1=0;
530
531 int npat=local_gAli->Npatterns();
532 Bool_t edge_npat_upper=kFALSE;
533 Bool_t edge_npat_lower=kFALSE;
534 Int_t Factor=-1;
535
536 if (patterloop_cnt==npat-1) {
537 edge_npat_upper=kTRUE;
538 }
539 if (patterloop_cnt==0) {
540 edge_npat_lower=kTRUE;
541 }
542
543 if (!edge_npat_lower) {
544 Z_minus1=local_gAli->GetPattern(patterloop_cnt-1)->Z();
545// cout << "WHAT IS GREATER? Z_normal Z_minus1 " << Z_normal << " " << Z_minus1 << endl;
546 //Factor=(Int_t)TMath::Sign(Z_normal,Z_minus1);
547 if (Z_minus1<Z_normal) Factor=1;
548 }
549 if (!edge_npat_upper) {
550 Z_plus1 =local_gAli->GetPattern(patterloop_cnt+1)->Z();
551// cout << "WHAT IS GREATER? Z_normal Z_plus1 " << Z_normal << " " << Z_plus1 << endl;
552 if (Z_plus1>Z_normal) Factor=1;
553 }
554
555// cout << Z_minus1 << endl;
556// cout << Z_normal << endl;
557// cout << Z_plus1 << endl;
558// cout << "Is edge_npat_lower = " << edge_npat_lower << endl;
559// cout << "Is edge_npat_upper = " << edge_npat_upper << endl;
560// cout << "Factor = " << Factor << endl;
561
562 // New PID we want to have:
563 Int_t patterloop_test=patterloop_cnt+Factor*n_patterns*BeforeOrAfter;
564// cout << "So NEW n_patterns would be " << patterloop_test << endl;
565
566 // Does this plate exist? If not, return 0 directly:
567 if (patterloop_test>=npat || patterloop_test<0) {
568 //cout << "So NEW n_patterns would be " << patterloop_test << " BUT IT DOES NOT MATHC in our local_gAli sceheme, which means its not existing. RETURNING 0 " << endl;
569 return 0;
570 }
571
572 // Since we have checked now for bounds we can FindCompliments:
573 TObjArray array;
574 array.Clear();
575 EdbPattern* TestPattern= (EdbPattern*)local_gAli->GetPattern(patterloop_test);
576 TestPattern -> FillCell(20,20,0.01,0.01);
577 int n_return = TestPattern->FindCompliments(*seg,array,3,3);
578 //cout << " Found " << n_return << " compliments in 2,2 sigma area:" << endl;
579
580 return n_return;
581}

◆ GetWeightFileString()

TString EdbShowAlg_NN::GetWeightFileString ( )
inline
113 {
114 return eWeightFileString;
115 }
TString eWeightFileString
Definition: EdbShowAlg_NN.h:50

◆ Init()

void EdbShowAlg_NN::Init ( void  )
53{
54 Log(2,"EdbShowAlg_NN::EdbShowAlg_NN","Init()");
55
56 // Init with values according to NN Alg:
57 eParaValue[0]=5;
58 eParaString[0]="INPUTNEURONS";
59 eParaValue[1]=0.75;
60 eParaString[1]="OUTPUT"; // NN output value
61
62 eWeightFileString="weights.txt";
64
66 //cout << " eANNTree->Show(0) " <<endl; eANNTree->Show(0);
68
69 // Standard Weights:
72
73 return;
74}
TString eWeightFileLayoutString
Definition: EdbShowAlg_NN.h:51
void SetANNWeightString()
Definition: EdbShowAlg_NN.cxx:175
void LoadANNWeights()
Definition: EdbShowAlg_NN.cxx:193
void CreateANNTree()
Definition: EdbShowAlg_NN.cxx:87
TMultiLayerPerceptron * Create_NN_ALG_MLP(TTree *inputtree, Int_t inputneurons)
Definition: EdbShowAlg_NN.cxx:139
TString eParaString[10]
Definition: EdbShowAlg.h:53

◆ Initialize()

void EdbShowAlg_NN::Initialize ( )
virtual

Reimplemented from EdbShowAlg.

80{
81 Log(2,"EdbShowAlg_NN::EdbShowAlg_NN","Initialize()");
82 return;
83}

◆ LoadANNWeights() [1/2]

void EdbShowAlg_NN::LoadANNWeights ( )
194{
195 eTMlpANN->LoadWeights(eWeightFileString);
196 return;
197}

◆ LoadANNWeights() [2/2]

void EdbShowAlg_NN::LoadANNWeights ( TMultiLayerPerceptron *  TMlpANN,
TString  WeightFileString 
)
204{
205 TMlpANN->LoadWeights(WeightFileString);
206 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_NN::LoadANNWeights " << WeightFileString << endl;
207 return;
208}

◆ Print()

void EdbShowAlg_NN::Print ( )

◆ SetANNWeightString()

void EdbShowAlg_NN::SetANNWeightString ( )
176{
177 Log(2,"EdbShowAlg_NN::SetANNWeightString","SetANNWeightString()");
178 int inputneurons=eParaValue[0];
179 if (inputneurons==5) eWeightFileString="ANN_WEIGHTS_PARASET_0_To_20.txt";
180 if (inputneurons==10) eWeightFileString="ANN_WEIGHTS_PARASET_20_To_40.txt";
181 if (inputneurons==20) eWeightFileString="ANN_WEIGHTS_PARASET_40_To_60.txt";
182 if (inputneurons==30) eWeightFileString="ANN_WEIGHTS_PARASET_60_To_80.txt";
183
184 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_NN::SetANNWeightString " << eWeightFileString << endl;
185 Log(2,"EdbShowAlg_NN::SetANNWeightString","SetANNWeightString()...done.");
186 return;
187}

◆ SetWeightFileString()

void EdbShowAlg_NN::SetWeightFileString ( TString  WeightFileString)
inline
109 {
110 eWeightFileString=WeightFileString;
111 return;
112 }

Member Data Documentation

◆ eANN_Inputtype

Int_t EdbShowAlg_NN::eANN_Inputtype
private

◆ eANN_var_dR_InBT_To_TestBT

Float_t EdbShowAlg_NN::eANN_var_dR_InBT_To_TestBT
private

◆ eANN_var_dR_NextBT_To_TestBT

Float_t EdbShowAlg_NN::eANN_var_dR_NextBT_To_TestBT
private

◆ eANN_var_dR_TestBT_To_InBT

Float_t EdbShowAlg_NN::eANN_var_dR_TestBT_To_InBT
private

◆ eANN_var_dT_InBT_To_TestBT

Float_t EdbShowAlg_NN::eANN_var_dT_InBT_To_TestBT
private

◆ eANN_var_dT_NextBT_To_TestBT

Float_t EdbShowAlg_NN::eANN_var_dT_NextBT_To_TestBT
private

◆ eANN_var_InBT_To_TestBT

Float_t EdbShowAlg_NN::eANN_var_InBT_To_TestBT
private

◆ eANN_var_mean_dR_2after

Float_t EdbShowAlg_NN::eANN_var_mean_dR_2after
private

◆ eANN_var_mean_dR_2before

Float_t EdbShowAlg_NN::eANN_var_mean_dR_2before
private

◆ eANN_var_mean_dR_after

Float_t EdbShowAlg_NN::eANN_var_mean_dR_after
private

◆ eANN_var_mean_dR_before

Float_t EdbShowAlg_NN::eANN_var_mean_dR_before
private

◆ eANN_var_mean_dR_same

Float_t EdbShowAlg_NN::eANN_var_mean_dR_same
private

◆ eANN_var_mean_dT_2after

Float_t EdbShowAlg_NN::eANN_var_mean_dT_2after
private

◆ eANN_var_mean_dT_2before

Float_t EdbShowAlg_NN::eANN_var_mean_dT_2before
private

◆ eANN_var_mean_dT_after

Float_t EdbShowAlg_NN::eANN_var_mean_dT_after
private

◆ eANN_var_mean_dT_before

Float_t EdbShowAlg_NN::eANN_var_mean_dT_before
private

◆ eANN_var_mean_dT_same

Float_t EdbShowAlg_NN::eANN_var_mean_dT_same
private

◆ eANN_var_min_dR_2after

Float_t EdbShowAlg_NN::eANN_var_min_dR_2after
private

◆ eANN_var_min_dR_2before

Float_t EdbShowAlg_NN::eANN_var_min_dR_2before
private

◆ eANN_var_min_dR_after

Float_t EdbShowAlg_NN::eANN_var_min_dR_after
private

◆ eANN_var_min_dR_before

Float_t EdbShowAlg_NN::eANN_var_min_dR_before
private

◆ eANN_var_min_dR_same

Float_t EdbShowAlg_NN::eANN_var_min_dR_same
private

◆ eANN_var_min_dT_2after

Float_t EdbShowAlg_NN::eANN_var_min_dT_2after
private

◆ eANN_var_min_dT_2before

Float_t EdbShowAlg_NN::eANN_var_min_dT_2before
private

◆ eANN_var_min_dT_after

Float_t EdbShowAlg_NN::eANN_var_min_dT_after
private

◆ eANN_var_min_dT_before

Float_t EdbShowAlg_NN::eANN_var_min_dT_before
private

◆ eANN_var_min_dT_same

Float_t EdbShowAlg_NN::eANN_var_min_dT_same
private

◆ eANN_var_nseg_1after

Int_t EdbShowAlg_NN::eANN_var_nseg_1after
private

◆ eANN_var_nseg_1before

Int_t EdbShowAlg_NN::eANN_var_nseg_1before
private

◆ eANN_var_nseg_2after

Int_t EdbShowAlg_NN::eANN_var_nseg_2after
private

◆ eANN_var_nseg_2before

Int_t EdbShowAlg_NN::eANN_var_nseg_2before
private

◆ eANN_var_nseg_3after

Int_t EdbShowAlg_NN::eANN_var_nseg_3after
private

◆ eANN_var_nseg_3before

Int_t EdbShowAlg_NN::eANN_var_nseg_3before
private

◆ eANN_var_nseg_same

Int_t EdbShowAlg_NN::eANN_var_nseg_same
private

◆ eANN_var_SpatialDist_TestBT_To_InBT

Float_t EdbShowAlg_NN::eANN_var_SpatialDist_TestBT_To_InBT
private

◆ eANN_var_zDiff_TestBT_To_InBT

Float_t EdbShowAlg_NN::eANN_var_zDiff_TestBT_To_InBT
private

◆ eANN_var_zDist_TestBT_To_InBT

Float_t EdbShowAlg_NN::eANN_var_zDist_TestBT_To_InBT
private

◆ eANNTree

TTree* EdbShowAlg_NN::eANNTree
private

◆ eTMlpANN

TMultiLayerPerceptron* EdbShowAlg_NN::eTMlpANN
private

◆ eWeightFileLayoutString

TString EdbShowAlg_NN::eWeightFileLayoutString
private

◆ eWeightFileString

TString EdbShowAlg_NN::eWeightFileString
private

The documentation for this class was generated from the following files: