]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSNeuralTracking.C
Transition to NewIO
[u/mrichter/AliRoot.git] / ITS / AliITSNeuralTracking.C
1 // ARGUMENTS:
2 // 1. number of azymuthal sectors (it's better not to go under 8 or over 40)
3 // 2. the ROOT file to read (WITHOUT exstension)
4 // 3. event number
5 // 4. if specified a string, a fstream named like the argument is opened and
6 //    the elapsed CPU time is stored (not useful)
7
8 // the macro will save a file named, for example "galice_<nsecs>.root"
9 // containing may AliITSneuralTrack objects
10
11 void AliITSNeuralTracking
12 (Int_t nsecs = 12, const char* rfile = "galice", Int_t event = 0, const char* save = 0)
13 {
14         TStopwatch timer;
15         Double_t CONVERT = TMath::Pi() / 180.0;
16         const char* wfile = Form("%s_%d.root", rfile, nsecs);
17         cout << "Reading file " << rfile << ".root and saving in " << wfile << endl;
18         
19 // ==================================
20 // ==== CURVATURE CUT DEFINITION ====
21 // ==================================
22
23         // These values define the curvature cuts for all steps 
24         // within a sector.
25         // For a greater clarity, the cuts are indicated in units
26         // of transverse momentum (GeV/c) but these value have no
27         // exact physical meaning, but are useful to understand
28         // well what means a choice in the value of a certain 
29         // curvature constraint
30         // NOTE: becareful to make sure that the 'ncuts' variable
31         //       have the same value of the dimension of the allocated arrays
32         
33         Int_t ncuts;
34         Double_t *p, *cut;
35         
36         ncuts = 5;
37         p = new Double_t[5];
38         cut = new Double_t[5];  
39         p[0] = 2.0;
40         p[1] = 1.0;
41         p[2] = 0.7;
42         p[3] = 0.5;
43         p[4] = 0.3;
44         
45         for (Int_t i = 0; i < ncuts; i++) cut[i] = 0.003 * 0.2 / p[i];
46         
47
48 // ==========================
49 // ==== OTHER PARAMETERS ====
50 // ==========================
51             
52         Bool_t   flag   = kFALSE; // for now, don't change this line, please...
53         
54         Double_t diff   = 0.02;   // helicoidal cut 
55         Double_t dtheta = 1.0;    // delta-theta cut
56         Double_t temp   = 1.0;    // temperature parameter
57         Double_t var    = 0.0001; // stabilization threshold
58         
59         Double_t exp    = 7.0;    // straight-line excitator
60         Double_t gtoc   = 3.0;    // gain/cost contribution ratio
61         
62         Double_t min    = 0.4;    // minimum in random activation initialization
63         Double_t max    = 0.6;    // maximum in random activation initialization
64         
65         
66 // =========================
67 // ==== NEURAL TRACKING ====
68 // =========================
69         
70         AliITSneuralTracker *ANN = new AliITSneuralTracker(nsecs, ncuts, cut, CONVERT*dtheta);
71                 
72         TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(Form("%s.root", rfile));
73         if (!file) file = new TFile(Form("%s.root", rfile),"UPDATE");
74
75         //Double_t Xv = -0.001097;
76         //Double_t Yv = -0.00347647;
77         //Double_t Zv =  0.000631345;
78         //ANN->SetVertex(Xv, Yv, Zv);
79         // You should find the vertex with VertexMacro.C
80         // and then put by hand the found values with
81         // the above method.
82         
83         Int_t points = ANN->ReadFile(Form("%s.root", rfile), event);
84                 
85         ANN->SetTemperature(temp);
86         ANN->SetVariationLimit(var);
87         ANN->SetGainToCostRatio(gtoc);
88         ANN->SetExponent(exp);
89         ANN->SetInitLimits(min, max);
90         ANN->SetDiff(diff);
91         
92         cout << points << " points found " << endl << endl;
93         
94         TStopwatch timer;
95         timer.Start();
96         
97         ANN->Go(wfile, flag);
98         
99         timer.Stop();
100         cout << endl;
101         timer.Print();
102         
103         if (save) {
104                 fstream ksave(save, ios::app);
105                 ksave << nsecs << " " << timer->CpuTime() << endl;
106                 ksave.close();
107         }
108         
109 //      delete gAlice;
110 //      gAlice = 0;
111 }