KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVRandomizor.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Fri Jun 23 12:02:13 2017
2 //Author: Eric BONNET
3 
4 #include "KVRandomizor.h"
5 #include "TMath.h"
6 
7 #include "TH1.h"
8 #include "TH2.h"
9 #include "TH3.h"
10 #include "TRandom3.h"
11 
13 
14 
15 //-------------------------
17  fNd(ndim), fNdMax(6), fMin(ndim), fMax(ndim)
18 //-------------------------
19 
22 
23 {
24  // Default constructor
25  if (fNd > fNdMax) {
26  Warning("KVRandomizor", "Too high dimensions (max : %d)", fNdMax);
27  //return; MAYBE THROW AN EXCEPTION?
28  }
29 }
30 
31 
32 //-------------------------
34 //-------------------------
35 
36 
39 {
40  // Destructor
41 }
42 
43 
44 //-------------------------
46 //-------------------------
47 
48 
50 {
51  fVmax = vmax;
52  fVmin = vmin;
53 }
54 
55 
56 //-------------------------
58 //-------------------------
59 
60 
62 {
63  for (Int_t ii = 0; ii < fNd; ii += 1) {
64  fMin[ii] = min[ii];
65  fMax[ii] = max[ii];
66  }
67 }
68 
69 
70 //-------------------------
71 std::vector<Double_t> KVRandomizor::GetPosition()
72 //-------------------------
73 
74 
76 {
77 
78  std::vector<Double_t> pos(fNd);
79  for (Int_t ii = 0; ii < fNd; ii += 1) {
80  pos[ii] = GetPosition(ii);
81  }
82  return pos;
83 
84 }
85 
86 
87 //-------------------------
89 //-------------------------
90 
91 
93 {
94 
95  return gRandom->Uniform(fMin[idx], fMax[idx]);
96 
97 }
98 
99 
100 //-------------------------
102 //-------------------------
103 
104 
106 {
107 
108  Info("ComputeValue", "To be defined in child class");
109  return 1;
110 
111 }
112 
113 
114 //-------------------------
116 //-------------------------
117 
118 
120 {
121  Double_t prob = gRandom->Uniform(fVmin, fVmax);
122 
123  return (prob <= val);
124 }
125 
126 
127 //-------------------------
129 //-------------------------
130 
131 
133 {
134  fNtest = 0;
135  TH1* h1 = 0;
136  if (fNd == 1) {
137  h1 = FillHisto1D(ntimes);
138  }
139  else if (fNd == 2) {
140  h1 = FillHisto2D(ntimes);
141  }
142  else if (fNd == 3) {
143  h1 = FillHisto3D(ntimes);
144  }
145  else {
146  Info("FillHisto", "method implemented for dimansion <= 3\nDo nothing ...");
147  return h1;
148  }
149  Info("FillHisto", "stat=%d in %d tentatives => %lf", ntimes, fNtest, Double_t(ntimes) / fNtest);
150 
151  return h1;
152 }
153 
154 
155 
156 //-------------------------
158 //-------------------------
159 
160 
162 {
163  TH1* h1 = new TH1F(
164  Form("histo%dD", fNd),
165  Form("histo%dD", fNd),
166  (fMax[0] - fMin[0]) * 100, fMin[0], fMax[0]
167  );
168  for (Int_t ii = 0; ii < ntimes; ii += 1) {
169  fNtest += 1;
170  std::vector<Double_t> pos = GetPosition();
171  if (TestValue(ComputeValue(&pos[0]))) {
172  h1->Fill(pos[0]);
173  }
174  else {
175  ii -= 1;
176  }
177  }
178 
179  return h1;
180 }
181 
182 
183 //-------------------------
185 //-------------------------
186 
187 
189 {
190  TH2* h1 = new TH2F(
191  Form("histo%dD", fNd),
192  Form("histo%dD", fNd),
193  (fMax[0] - fMin[0]) * 100, fMin[0], fMax[0],
194  (fMax[1] - fMin[1]) * 100, fMin[1], fMax[1]
195  );
196  for (Int_t ii = 0; ii < ntimes; ii += 1) {
197  fNtest += 1;
198  std::vector<Double_t> pos = GetPosition();
199  if (TestValue(ComputeValue(&pos[0]))) {
200  h1->Fill(pos[0], pos[1]);
201  }
202  else {
203  ii -= 1;
204  }
205  }
206  return h1;
207 }
208 
209 
210 //-------------------------
212 //-------------------------
213 
214 
216 {
217  TH3* h1 = new TH3F(
218  Form("histo%dD", fNd),
219  Form("histo%dD", fNd),
220  (fMax[0] - fMin[0]) * 20, fMin[0], fMax[0],
221  (fMax[1] - fMin[1]) * 20, fMin[1], fMax[1],
222  (fMax[2] - fMin[2]) * 20, fMin[2], fMax[2]
223  );
224  for (Int_t ii = 0; ii < ntimes; ii += 1) {
225  fNtest += 1;
226  std::vector<Double_t> pos = GetPosition();
227  if (TestValue(ComputeValue(&pos[0]))) {
228  h1->Fill(pos[0], pos[1], pos[2]);
229  }
230  else {
231  ii -= 1;
232  }
233  }
234 
235  return h1;
236 }
237 
238 
239 //-------------------------
241 //-------------------------
242 
243 
245 {
246  Double_t min[3] = { -5, -5, -5};
247  Double_t max[3] = {5, 5, 5};
248 
249  this->SetRange(min, max);
250  this->SetExtrema(1, 0);
251 
252  return this->FillHisto(ntimes);
253 }
254 
255 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
bool Bool_t
double Double_t
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
Test of generic class to perform sample on functions.
Definition: KVRandomizor.h:18
virtual Double_t ComputeValue(Double_t *pos)
std::vector< Double_t > fMax
Definition: KVRandomizor.h:28
TH1 * Test(Int_t ntimes=1000)
Double_t fVmin
Definition: KVRandomizor.h:26
virtual ~KVRandomizor()
Destructor.
TH1 * FillHisto2D(Int_t)
void SetRange(Double_t *min, Double_t *max)
Bool_t TestValue(Double_t)
TH1 * FillHisto(Int_t ntimes=1000)
std::vector< Double_t > fMin
Definition: KVRandomizor.h:27
Double_t fVmax
Definition: KVRandomizor.h:25
TH1 * FillHisto3D(Int_t)
TH1 * FillHisto1D(Int_t)
void SetExtrema(Double_t, Double_t vmin=0)
KVRandomizor(const char* name, const char* title);.
std::vector< Double_t > GetPosition()
virtual Int_t Fill(const char *name, Double_t w)
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Double_t Uniform(Double_t x1, Double_t x2)
TH1F * h1
void Warning(const char *location, const char *va_(fmt),...)