KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVSpiderLine.cpp
Go to the documentation of this file.
1 #include "KVSpiderLine.h"
2 #include <TMath.h>
3 
4 using namespace std;
5 
7 
8 
9 
12 {
13  SetName(Form("Z=%d", z_));
14  _z = z_;
15  _a = -1;
16  _filled = false;
17 
18  _line = new TGraph();
19  _line->SetName(GetName());
20  _iline = new TGraph();
21  _iline->SetName(Form("I%s", GetName()));
22 
23  _ff = 0;
24  _pow = 0;
25  _fitStatus = 1;
26  _pdy = pdy;
27  ResetCounter();
28  _nAcceptedPoints = 100000;
29 }
30 
31 
32 
34 
36 {
37  SetName(Form("Z=%d,A=%d", z_, a_));
38  _z = z_;
39  _a = a_;
40  _filled = false;
41 
42  _line = new TGraph();
43  _line->SetName(GetName());
44  _iline = new TGraph();
45  _iline->SetName(Form("I%s", GetName()));
46 
47  _ff = 0;
48  _pow = 0;
49  _fitStatus = 1;
50  _pdy = 0;
51  ResetCounter();
52  _nAcceptedPoints = 100000;
53 }
54 
55 
56 
58 
60 {
61  _line = 0;
62  _iline = 0;
63  ResetCounter();
64  _nAcceptedPoints = 100000;
65 }
66 
67 
68 
69 
71 
72 void KVSpiderLine::SetZ(int z_)
73 {
74  SetName(Form("Z=%d", z_));
75  if (_a > 0) SetName(Form("%s,A=%d", GetName(), _a));
76  _z = z_;
77  _line->SetName(GetName());
78  _iline->SetName(Form("I%s", GetName()));
79 }
80 
81 
82 
84 
85 void KVSpiderLine::SetA(int a_)
86 {
87  SetName(Form("Z=%d,A=%d", _z, a_));
88  _a = a_;
89  _line->SetName(GetName());
90  _iline->SetName(Form("I%s", GetName()));
91 }
92 
93 
94 
95 
97 
98 bool KVSpiderLine::AddPoint(double x_, double y_, bool test_, int n_)
99 {
100  if (!CheckStatus()) return false;
101  int n;
102  bool valid = true;
103 
104  if (n_ == -1) n = _line->GetN();
105  else n = n_;
106 
107  if (test_) valid = TestPoint(x_, y_);
108  if (valid) {
109  _line->SetPoint(n, x_, y_);
110  _iline->SetPoint(_iline->GetN(), x_, y_);
111  }
112 
113  _pointsCounter++;
114  return valid;
115 }
116 
117 
118 
119 
121 
122 bool KVSpiderLine::AddInterpolatePoint(double x_, double y_, bool test_, int n_)
123 {
124  if (!CheckStatus()) return false;
125  int n;
126  bool valid = true;
127 
128  if (n_ == -1) n = _iline->GetN();
129  else n = n_;
130 
131  if (test_) valid = TestPoint(x_, y_);
132  if (valid) _iline->SetPoint(n, x_, y_);
133 
134  return valid;
135 }
136 
137 
138 
139 
141 
143 {
144  if (!f) return;
145  if (!CheckStatus()) return;
146 
147  _line->Apply(f);
148  _iline->Apply(f);
149 
150  _ff = 0;
151 
152  return;
153 }
154 
155 
156 
157 
159 
160 bool KVSpiderLine::ReplaceLastPoint(double x_, double y_)
161 {
162  if (!CheckStatus()) return false;
163 
164  _line->SetPoint(GetN() - 1, x_, y_);
165  _iline->SetPoint(GetInterpolateN() - 1, x_, y_);
166 
167  return true;
168 }
169 
170 
171 
172 
174 
175 void KVSpiderLine::Sort(bool ascending_)
176 {
177  if (!CheckStatus()) return;
178 
179  _line->Sort(&TGraph::CompareX, ascending_);
180  _iline->Sort(&TGraph::CompareX, ascending_);
181 
182  return;
183 }
184 
185 
186 
187 
189 
191 {
192  return _filled;
193 }
194 
195 
196 
197 
199 
200 void KVSpiderLine::SetStatus(bool filled_)
201 {
202  _filled = filled_;
203 }
204 
205 
206 
207 
209 
210 TF1* KVSpiderLine::GetFunction(double min_, double max_)
211 {
212  if (!CheckStatus()) return 0;
213 
214  double min;
215  double max;
216 
217  double xtest = GetX();
218  if (GetX(0) > GetX()) xtest = GetX(0);
219  double ytest = GetY();
220  if (GetY(0) > GetY()) ytest = GetY(0);
221 
222  double p0 = TMath::Power(xtest, 0.4) * ytest;
223 
224  if (min_ == -1.) {
225  if (GetX(0) < GetX()) min = GetX(0) - 10;
226  else min = GetX() - 10;
227  if (min <= 0.) min += 10.;
228  }
229  else min = min_;
230  if (max_ == -1.) {
231  if (GetX(0) < GetX()) max = GetX() + 10.;
232  else max = GetX(0) + 10.;
233  }
234  else max = max_;
235 
236  if (!_ff) {
237  _ff = new TF1(GetName(), Form("[0]*TMath::Power(x,%lf)/(TMath::Power((x+[1]),[2]))", _pow), min, max);
238  _ff->SetParameters(p0, 100., 0.4);
239  }
240  else {
241  double fmin, fmax;
242  _ff->GetRange(fmin, fmax);
243  if ((min <= fmin) || (max >= fmax)) {
244  _ff->SetRange(min, max);
245  }
246  }
247 
248  _fitStatus = _line->Fit(_ff, "QN");
249  return _ff;
250 }
251 
252 
253 
255 
256 bool KVSpiderLine::TestPoint(double x_, double y_, double dy_, bool fit)
257 {
258  if (!CheckStatus()) return false;
259  if (_pointsCounter >= _nAcceptedPoints) return false;
260 
261  TF1* locf = 0;
262  if ((GetN() >= 10) && (fit)) {
263  locf = GetFunction();
264  if (_fitStatus != 0) fit = false;
265  }
266 
267  bool valid = true;
268 
269  double dy;
270  if (dy_ > 0.) dy = dy_;
271  else {
272  dy = ((GetInterpolateY() - _pdy) / _z) * 1.;
273  }
274 
275  if ((GetN() >= 10) && (fit)) {
276  double yext = locf->Eval(x_);
277  if ((TMath::Abs(yext - y_) >= dy)) {
278  valid = false;
279  }
280  }
281  else if ((TMath::Abs(GetInterpolateY() - y_) >= dy)) valid = false;
282  return valid;
283 }
284 
285 
286 
288 
289 double KVSpiderLine::GetDistance(double x_, double y_)
290 {
291  double xx = x_;
292  double yy = y_;
293  double ox = GetInterpolateX();
294  double oy = GetInterpolateY();
295 
296  double dist = TMath::Sqrt((xx - ox) * (xx - ox) + (yy - oy) * (yy - oy));
297  return dist;
298 }
299 
300 
301 
303 
305 {
306  if (!CheckStatus()) return false;
307  return _line->GetN();
308 }
309 
310 
311 
312 
314 
316 {
317  if (!CheckStatus()) return false;
318  return _iline->GetN();
319 }
320 
321 
322 
323 
324 
326 
327 double KVSpiderLine::GetX()const
328 {
329  if (!CheckStatus()) return false;
330  return GetX(_line->GetN() - 1);
331 }
332 
333 
334 
336 
337 double KVSpiderLine::GetX(int n_)const
338 {
339  if (!CheckStatus()) return false;
340  int n;
341  if ((n_ <= _line->GetN() - 1) && (n_ >= 0)) n = n_;
342  else {
343  cout << "WARNING: KVSpiderLine::GetX(): Invalid index n = '" << n_ << "'." << endl;
344  return -1.;
345  }
346  return _line->GetX()[n];
347 }
348 
349 
350 
351 
353 
354 double KVSpiderLine::GetY()const
355 {
356  if (!CheckStatus()) return false;
357  return GetY(_line->GetN() - 1);
358 }
359 
360 
361 
363 
364 double KVSpiderLine::GetY(int n_)const
365 {
366  if (!CheckStatus()) return false;
367  int n;
368  if ((n_ <= _line->GetN() - 1) && (n_ >= 0)) n = n_;
369  else {
370  cout << "WARNING: KVSpiderLine::GetY(): Invalid index n = '" << n_ << "'." << endl;
371  return -1.;
372  }
373  return _line->GetY()[n];
374 }
375 
376 
377 
378 
379 
381 
383 {
384  if (!CheckStatus()) return false;
385  return GetInterpolateX(_iline->GetN() - 1);
386 }
387 
388 
389 
391 
392 double KVSpiderLine::GetInterpolateX(int n_)const
393 {
394  if (!CheckStatus()) return false;
395  int n;
396  if ((n_ <= _iline->GetN() - 1) && (n_ >= 0)) n = n_;
397  else {
398  cout << "WARNING: KVSpiderLine::GetInterpolateX(): Invalid index n = '" << n_ << "'." << endl;
399  return -1.;
400  }
401  return _iline->GetX()[n];
402 }
403 
404 
405 
406 
408 
410 {
411  if (!CheckStatus()) return false;
412  return GetInterpolateY(_iline->GetN() - 1);
413 }
414 
415 
416 
418 
419 double KVSpiderLine::GetInterpolateY(int n_)const
420 {
421  if (!CheckStatus()) return false;
422  int n;
423  if ((n_ <= _iline->GetN() - 1) && (n_ >= 0)) n = n_;
424  else {
425  cout << "WARNING: KVSpiderLine::GetX(): Invalid index n = '" << n_ << "'." << endl;
426  return -1.;
427  }
428  return _iline->GetY()[n];
429 }
430 
431 
432 
433 
435 
437 {
438  TString opt(opt_);
439  TString com("");
440  if (!CheckStatus()) return;
441 
442  if (opt.Contains("A")) com += "A";
443  if (opt.Contains("same")) com += "same";
444 
445  if (opt.Contains("I")) {
446  _iline->SetMarkerColor(kGreen);
447  _iline->SetMarkerSize(0.9);
448  _iline->SetMarkerStyle(21);
449  _iline->SetLineWidth(1);
450  _iline->SetLineColor(kGreen);
451  _iline->Draw(Form("PL%s", com.Data()));
452  if (com.Contains("A")) com.ReplaceAll("A", "");
453  }
454 
455  if (opt.Contains("L")) {
456  _line->SetMarkerColor(kRed);
457  _line->SetMarkerSize(1.3);
458  _line->SetMarkerStyle(20);
459  _line->SetLineWidth(2);
460  _line->SetLineColor(kRed);
461  if (_line->GetN() <= 1) _line->Draw(Form("P%s", com.Data()));
462  else _line->Draw(Form("PL%s", com.Data()));
463  }
464 
465  if (opt.Contains("N")) {
466  TLatex* zname = 0;
467  if (GetX(0) > GetX()) zname = new TLatex(GetX(0), GetY(0), GetName());
468  else zname = new TLatex(GetX(), GetY(), GetName());
469  zname->SetTextSize(0.02);
470  zname->Draw("same");
471  }
472 }
473 
474 
475 
476 
478 
480 {
481  if ((!_line) || (!_iline)) return false;
482  return true;
483 }
484 
485 
486 
487 
488 
489 
490 
491 
492 
493 
494 
495 
496 
497 
498 
499 
500 
501 
502 
503 
504 
505 
506 
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define f(i)
double Double_t
const char Option_t
kRed
kGreen
char * Form(const char *fmt,...)
Part of Spider Identification.
Definition: KVSpiderLine.h:18
void Sort(bool ascending_=true)
double GetY() const
bool ReplaceLastPoint(double x_, double y_)
void Draw(Option_t *opt_="")
double GetInterpolateX() const
virtual bool TestPoint(double x_, double y_, double dy_=-1., bool fit=true)
int GetInterpolateN() const
void SetZ(int z_)
void SetA(int a_)
int GetN() const
double GetDistance(double x_, double y_)
bool CheckStatus() const
double GetX() const
bool AddInterpolatePoint(double x_, double y_, bool test_=false, int n_=-1)
bool AddPoint(double x_, double y_, bool test_=false, int n_=-1)
virtual TF1 * GetFunction(double min_=-1., double max_=-1.)
void SetStatus(bool filled_=true)
void Apply(TF1 *f)
double GetInterpolateY() const
virtual void SetTextSize(Float_t tsize=1)
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
static Bool_t CompareX(const TGraph *gr, Int_t left, Int_t right)
virtual void Draw(Option_t *option="")
const char * Data() const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
const Int_t n
gr SetName("gr")
def fit(model, train_loader, val_loader, num_epochs, batch_size, optimizer, criterion, save_best, scheduler)
double dist(AxisAngle const &r1, AxisAngle const &r2)
std::function< T(T)> GetFunction(const std::string &name)
Double_t Power(Double_t x, Double_t y)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)