CbmRoot
Loading...
Searching...
No Matches
CbmKresSelectAnn.cxx
Go to the documentation of this file.
1/* Copyright (C) 2017-2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Ievgenii Kres, Florian Uhlig [committer] */
4
18#include "CbmKresSelectAnn.h"
19
20#include "TMath.h"
21#include "TSystem.h"
22#include "TTree.h"
23
24#include <boost/assign/list_of.hpp>
25
26#include <iostream>
27#include <string>
28#include <vector>
29
30
31using namespace std;
32
33CbmKresSelectAnn::CbmKresSelectAnn() : fAnnWeights(), fNN(nullptr) {}
34
36
38{
39 TTree* simu = new TTree("MonteCarlo", "MontecarloData");
40 Double_t x[6];
41 Double_t xOut;
42
43 simu->Branch("x0", &x[0], "x0/D");
44 simu->Branch("x1", &x[1], "x1/D");
45 simu->Branch("x2", &x[2], "x2/D");
46 simu->Branch("x3", &x[3], "x3/D");
47 simu->Branch("x4", &x[4], "x4/D");
48 simu->Branch("x5", &x[5], "x5/D");
49 simu->Branch("xOut", &xOut, "xOut/D");
50
51 fNN = new TMultiLayerPerceptron("x0,x1,x2,x3,x4,x5:12:xOut", simu);
52
53 fAnnWeights = string(gSystem->Getenv("VMCWORKDIR")) + "/analysis/conversion2/KresAnalysis_ann_weights.txt";
54 cout << "-I- CbmKresSelectAnn: get ANN weight parameters from: " << fAnnWeights << endl;
55 fNN->LoadWeights(fAnnWeights.c_str());
56}
57
58double CbmKresSelectAnn::DoSelect(double InvariantMass, double OpeningAngle, double PlaneAngle_last, double ZPos,
59 TVector3 Momentum1, TVector3 Momentum2)
60{
61 double AnnValue = 0;
62
63 double p1 =
64 TMath::Sqrt(Momentum1.X() * Momentum1.X() + Momentum1.Y() * Momentum1.Y() + Momentum1.Z() * Momentum1.Z());
65 double p2 =
66 TMath::Sqrt(Momentum2.X() * Momentum2.X() + Momentum2.Y() * Momentum2.Y() + Momentum2.Z() * Momentum2.Z());
67
68 if (InvariantMass > 0.5 || OpeningAngle > 45 || ZPos > 100 || p1 > 10 || p2 > 10) { return AnnValue = -1; }
69
70
71 double params[6];
72
73 params[0] = InvariantMass / 0.1;
74 params[1] = OpeningAngle / 30;
75 params[2] = PlaneAngle_last / 30;
76 params[3] = ZPos / 100;
77 params[4] = p1 / 5;
78 params[5] = p2 / 5;
79
80 if (params[0] > 1.0) params[0] = 1.0;
81 if (params[1] > 1.0) params[1] = 1.0;
82 if (params[2] > 1.0) params[2] = 1.0;
83 if (params[3] > 1.0) params[3] = 1.0;
84 if (params[4] > 1.0) params[4] = 1.0;
85 if (params[5] > 1.0) params[5] = 1.0;
86
87 AnnValue = fNN->Evaluate(0, params);
88
89 return AnnValue;
90}
double DoSelect(double InvariantMass, double OpeningAngle, double PlaneAngle_last, double ZPos, TVector3 Momentum1, TVector3 Momentum2)
std::string fAnnWeights
TMultiLayerPerceptron * fNN
Hash for CbmL1LinkKey.