KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVClust3D.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Tue Sep 27 10:46:46 2016
2 //Author: Eric BONNET
3 
4 #include "KVClust3D.h"
5 #include "KVNumberList.h"
6 #include "TMath.h"
7 
9 
10 //____________________________________________________________________________//
11 
12 
13 
16 void KVClust3D::init()
17 {
18  //initialization method
19  fNclusters = 0;
20  fVoisins = new TH2D("hvoisins", "FromKVClust3D", 1000, 0.5, 1000.5, 1000, 0.5, 1000.5);
21 
22  fNcells = new TArrayI(10);
23  fPop = new TArrayD(10);
24 
25  fThreshold = 0.0;
26  fNvoisins = 1;
27  fKeepLonelyCells = kTRUE;
28 
29  fCluster = ProduceTH3D("cluster_index");
30 }
31 
32 
33 
34 
37 
39  : TH3D()
40 {
41  // Default constructor
42 }
43 
44 
45 
46 
49 
50 KVClust3D::KVClust3D(const char* name, const char* title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup, Int_t nbinsz, Double_t zlow, Double_t zup)
51  : TH3D(name, title, nbinsx, xlow, xup, nbinsy, ylow, yup, nbinsz, zlow, zup)
52 {
53  // Constructor inherited from TH3D
54  init();
55 }
56 
57 
58 
59 
62 
63 KVClust3D::KVClust3D(const char* name, const char* title, Int_t nbinsx, const Float_t* xbins, Int_t nbinsy, const Float_t* ybins, Int_t nbinsz, const Float_t* zbins)
64  : TH3D(name, title, nbinsx, xbins, nbinsy, ybins, nbinsz, zbins)
65 {
66  // Constructor inherited from TH3D
67  init();
68 }
69 
70 
71 
72 
75 
76 KVClust3D::KVClust3D(const char* name, const char* title, Int_t nbinsx, const Double_t* xbins, Int_t nbinsy, const Double_t* ybins, Int_t nbinsz, const Double_t* zbins)
77  : TH3D(name, title, nbinsx, xbins, nbinsy, ybins, nbinsz, zbins)
78 {
79  // Constructor inherited from TH3D
80  init();
81 }
82 
83 
84 
85 
92 
94 {
95  // Destructor
96  /*
97  delete fCluster;
98  delete fVoisins;
99  delete fNcells;
100  delete fPop;
101  */
102 }
103 
104 
105 
106 
108 
110 {
111  fKeepLonelyCells = val;
112 }
113 
114 
115 
116 
118 
120 {
121  return fKeepLonelyCells;
122 }
123 
124 
125 
126 
128 
130 {
131  fNvoisins = val;
132 }
133 
134 
135 
136 
138 
140 {
141  return fNvoisins;
142 }
143 
144 
145 
146 
148 
150 {
151  return GetXaxis()->GetBinWidth(1) * GetYaxis()->GetBinWidth(1) * GetZaxis()->GetBinWidth(1);
152 }
153 
154 
155 
156 
158 
160 {
161  Double_t vol = GetVolumeCell();
162  Int_t ncell = TMath::Power(2 * fNvoisins + 1, 3.);
163  return ncell * vol;
164 }
165 
166 
167 
168 
172 
174 {
175  //Reset object associated to the KVClust3D object
176  //Don't reset the object it self
177  fNcells->Reset();
178  fPop->Reset();
179 
180  fCluster->Reset();
181  fVoisins->Reset();
182  fNclusters = 0;
183 }
184 
185 
186 
187 
190 
192 {
193  //Set threshold on the content of a cell to be considered has filled
194  fThreshold = value;
195 
196 }
197 
198 
199 
200 
203 
205 {
206  //Set density threshold on the content of a cell to be considered has filled
207  fThreshold = value * GetVolumeCell();
208 
209 }
210 
211 
212 
213 
216 
218 {
219  //return the value of the threshold for the content of a cell to be considered has filled
220  return fThreshold;
221 
222 }
223 
224 
225 
226 
228 
230 {
231 
232  printf("-----------------------------------------\n");
233  printf("volume cellule : %lf\n", GetVolumeCell());
234  printf("nombre de cellules voisines : %lf\n", TMath::Power(2.*fNvoisins + 1, 3.));
235  printf("volume voisin: %lf\n", GetVolumeVoisin());
236  printf("-----\nseuil :\n\t- contenu : %lf\n", GetThreshold());
237  printf("\t- densite : %lf\n", GetThreshold() / GetVolumeCell());
238 
239 
240 
241 }
242 
243 
244 
245 
250 
251 TH3D* KVClust3D::ProduceTH3D(const Char_t* name) const
252 {
253  //Create TH3D object, clone of the KVClust3D object
254  //user has to delete it after its use
255  //
256  TH3D* histo = new TH3D(name, "FromKVClust3D",
257  GetNbinsX(),
258  GetXaxis()->GetBinLowEdge(1),
259  GetXaxis()->GetBinLowEdge(GetNbinsX() + 1),
260  GetNbinsY(),
261  GetYaxis()->GetBinLowEdge(1),
262  GetYaxis()->GetBinLowEdge(GetNbinsY() + 1),
263  GetNbinsZ(),
264  GetZaxis()->GetBinLowEdge(1),
266  );
267 
268  return histo;
269 
270 }
271 
272 
273 
274 
284 
286 {
287  //Apply clusterization of the 3D histogram
288  //For each cell, with a content greater than the threshold defined before
289  // using the SetThreshold() method, test the content of the closest neighbors (8)
290  // to see if they are connected
291  // At the end, a clone histogram is available, using this method : GetClusterIndex().
292  // The content of its cells corresponds
293  // to the index of each found cluster
294  //
295  Int_t Ntemp = 1;
296  Int_t NCellsTot = 0;
297  Double_t PopTot = 0;
298 
299  ResetObjects();
300 
301 // if (!fCluster)
302 // fCluster = ProduceTH3D("cluster_index");
303 
304  //Remplissage de l histo 3D fCluster :
305  // Contenu de l objet KVClust3D au dessus du seuil -> contenu 1
306  //
307  for (Int_t nx = 1; nx <= GetNbinsX(); nx += 1) {
308  for (Int_t ny = 1; ny <= GetNbinsY(); ny += 1) {
309  for (Int_t nz = 1; nz <= GetNbinsZ(); nz += 1) {
310  if (GetBinContent(nx, ny, nz) >= fThreshold) {
311  fCluster->SetBinContent(nx, ny, nz, 1);
312  NCellsTot += 1;
313  PopTot += GetBinContent(nx, ny, nz);
314  }
315  else {
316  fCluster->SetBinContent(nx, ny, nz, 0);
317  }
318  }
319  }
320  }
321 
322  for (Int_t nx = 2; nx <= GetNbinsX() - 1; nx += 1) {
323  for (Int_t ny = 2; ny <= GetNbinsY() - 1; ny += 1) {
324  for (Int_t nz = 2; nz <= GetNbinsZ() - 1; nz += 1) {
325  Int_t idx = fCluster->GetBinContent(nx, ny, nz);
326  if (idx >= 1) { //Cellule avec la bonne densite
327  if (idx > 1) { //cellule deja affectee
328  Int_t ncv = -1;
329  for (Int_t nxv = nx - fNvoisins; nxv <= nx + fNvoisins; nxv += 1) { //on boucle sur les cellules voisines
330  for (Int_t nyv = ny - fNvoisins; nyv <= ny + fNvoisins; nyv += 1) {
331  for (Int_t nzv = nz - fNvoisins; nzv <= nz + fNvoisins; nzv += 1) {
332  Int_t idxv = fCluster->GetBinContent(nxv, nyv, nzv);
333  //cellule a affilier
334  if (idxv >= 1) {
335  ncv += 1;
336  if (idxv == 1) { //cellule mise dans le cluster idx
337  fNcells->AddAt(1 + fNcells->At(idx), idx);
338  fPop->AddAt(GetBinContent(nxv, nyv, nzv) + fPop->At(idx), idx);
339  fCluster->SetBinContent(nxv, nyv, nzv, idx);
340  }
341  else {
342  // clusters idx et idxv sont voisins
343  fVoisins->Fill(TMath::Min(idx, idxv), TMath::Max(idx, idxv));
344  }
345  }
346  else {
347  //cellule vide
348  //? SetBinContent(nxv,nyv,nzv,0);
349  //fCluster->SetBinContent(nxv,nyv,nzv,0);
350  }
351  }
352  }
353  }
354  if (ncv == 0)
355  printf("gros pb de coherence ... cluster %1.0lf declare avec une cellule isolee\n", fCluster->GetBinContent(nx, ny, nz));
356  }
357  else if (idx == 1) { //cellule pas encore affectee
358  Int_t ncv = -1;
359  KVNumberList nl;
360  for (Int_t nxv = nx - fNvoisins; nxv <= nx + fNvoisins; nxv += 1) { //on boucle sur les cellules voisines
361  for (Int_t nyv = ny - fNvoisins; nyv <= ny + fNvoisins; nyv += 1) {
362  for (Int_t nzv = nz - fNvoisins; nzv <= nz + fNvoisins; nzv += 1) {
363  Int_t idxv = fCluster->GetBinContent(nxv, nyv, nzv);
364  if (idxv >= 1) { //cellule voisine a affilier
365  if (idxv > 1) { // la cellule appartient a un cluster
366  //cellule voisine appartenant a un cluster deja identifie
367  nl.Add(idxv);
368  }
369  else {
370 
371  }
372  ncv += 1;
373  }
374  }
375  }
376  }
377  //au moins une cellule voisine non nulle
378  //
379  if (ncv >= 1) {
380  //les cellules voisines sont associees a un ou plusieur cluster
381  if (nl.GetNValues() >= 1) {
382  //on prend le premier cluster qui vient
383  Int_t newid = nl.First();
384  for (Int_t nxv = nx - fNvoisins; nxv <= nx + fNvoisins; nxv += 1) { //on boucle sur les cellules voisines
385  for (Int_t nyv = ny - fNvoisins; nyv <= ny + fNvoisins; nyv += 1) {
386  for (Int_t nzv = nz - fNvoisins; nzv <= nz + fNvoisins; nzv += 1) {
387  Int_t idxv = fCluster->GetBinContent(nxv, nyv, nzv);
388  //on rajoute les cellules au cluster choisi
389  if (idxv == 1) {
390  fNcells->AddAt(1 + fNcells->At(newid), newid);
391  fPop->AddAt(GetBinContent(nxv, nyv, nzv) + fPop->At(newid), newid);
392  fCluster->SetBinContent(nxv, nyv, nzv, newid);
393  }
394  }
395  }
396  }
397  if (nl.GetNValues() > 1) {
398  //on met les implementes les clusters contigus
399  nl.Begin();
400  while (!nl.End()) {
401  Int_t compt = nl.Next();
402  if (compt != newid) {
403  fVoisins->Fill(TMath::Min(compt, newid), TMath::Max(compt, newid));
404  }
405  }
406  }
407  }
408  //creation d un nouveau cluster
409  else {
410 
411  Ntemp += 1;
412  if (Ntemp > 999)
413  Warning("Clusterize", "Ntemp>999 ... should be some problems in the treatment of neigbouring clusters...");
414  //Info("Clusterize","New cluster %d",Ntemp);
415 
416  if (Ntemp >= fNcells->GetSize()) {
417  fNcells->Set(Ntemp + 1);
418  fPop->Set(Ntemp + 1);
419  }
420  for (Int_t nxv = nx - fNvoisins; nxv <= nx + fNvoisins; nxv += 1) { //on boucle sur les cellules voisines
421  for (Int_t nyv = ny - fNvoisins; nyv <= ny + fNvoisins; nyv += 1) {
422  for (Int_t nzv = nz - fNvoisins; nzv <= nz + fNvoisins; nzv += 1) {
423  Int_t idxv = fCluster->GetBinContent(nxv, nyv, nzv);
424  if (idxv == 1) {
425  fNcells->AddAt(1 + fNcells->At(Ntemp), Ntemp);
426  fPop->AddAt(GetBinContent(nxv, nyv, nzv) + fPop->At(Ntemp), Ntemp);
427  fCluster->SetBinContent(nxv, nyv, nzv, Ntemp);
428  }
429  }
430  }
431  }
432  }
433  }
434  else {
435  if (IsLonelyCellsAreKept()) {
436  Ntemp += 1;
437  if (Ntemp > 999)
438  Warning("Clusterize", "Ntemp>999 ... should be some problems in the treatment of neigbouring clusters...");
439  //Info("Clusterize","Lonely cell %d",Ntemp);
440  if (Ntemp >= fNcells->GetSize()) {
441  fNcells->Set(Ntemp + 1);
442  fPop->Set(Ntemp + 1);
443  }
444  fNcells->AddAt(1 + fNcells->At(Ntemp), Ntemp);
445  fPop->AddAt(GetBinContent(nx, ny, nz) + fPop->At(Ntemp), Ntemp);
446  fCluster->SetBinContent(nx, ny, nz, Ntemp);
447  }
448  else {
449  fCluster->SetBinContent(nx, ny, nz, 0);
450  }
451  }
452  }
453  }
454  else {
455 
456  }
457  }
458  }
459  }
460 
461  //nclus commence a 1
462  //le nombre actuel de cluster est nclus-1
463  //
464  //boucle sur les clusters voisins
465  //on rassemble les clusters avec des cellules contigues
466  //et on decrement en consequence le nombre de cluster
467  //
468  fNclusters = Ntemp - 1;
469 
470  for (Int_t ii = Ntemp; ii >= 2; ii -= 1) {
471  for (Int_t jj = ii + 1; jj <= Ntemp; jj += 1) {
472  Int_t contenu = fVoisins->GetBinContent(ii, jj);
473  if (contenu > 0) {
474  if (fNcells->At(jj) > 0) {
475 
476  fNcells->AddAt(fNcells->At(jj) + fNcells->At(ii), ii);
477  fPop->AddAt(fPop->At(jj) + fPop->At(ii), ii);
478  for (Int_t nx = 1; nx <= GetNbinsX(); nx += 1) {
479  for (Int_t ny = 1; ny <= GetNbinsY(); ny += 1) {
480  for (Int_t nz = 1; nz <= GetNbinsZ(); nz += 1) {
481  if (fCluster->GetBinContent(nx, ny, nz) == jj)
482  fCluster->SetBinContent(nx, ny, nz, ii);
483  }
484  }
485  }
486  fNcells->SetAt(0, jj);
487  fPop->SetAt(0, jj);
488 
489  fNclusters -= 1;
490  }
491  }
492  }
493  }
494 
495  //on passe les contenus en négatifs si ils appartiennent a aucun cluster
496  //
497  for (Int_t nx = 1; nx <= GetNbinsX(); nx += 1) {
498  for (Int_t ny = 1; ny <= GetNbinsY(); ny += 1) {
499  for (Int_t nz = 1; nz <= GetNbinsZ(); nz += 1) {
500  if (fCluster->GetBinContent(nx, ny, nz) < 2) {
501  fCluster->SetBinContent(nx, ny, nz, -1);
502  }
503  }
504  }
505  }
506 
507  //on reordonne les clusters, liste possible avec des clusters vides
508  //on reprend l indexage en partant de 1 jusqu a nclust
509  // pour avoir en sortie les contenus suivants dans hclust :
510  // <=0 la cellule appartient a aucun cluster
511  // contenu ii>0, la cellule appartient a iieme cluster
512  //
513  Int_t nreel = 1;
514  for (Int_t ii = 0; ii <= Ntemp; ii += 1) {
515  if (fNcells->At(ii) > 0) {
516  for (Int_t nx = 1; nx <= GetNbinsX(); nx += 1) {
517  for (Int_t ny = 1; ny <= GetNbinsY(); ny += 1) {
518  for (Int_t nz = 1; nz <= GetNbinsZ(); nz += 1) {
519  if (fCluster->GetBinContent(nx, ny, nz) == ii) {
520  fCluster->SetBinContent(nx, ny, nz, nreel);
521  }
522  }
523  }
524  }
525  fNcells->SetAt(fNcells->At(ii), nreel);
526  fPop->SetAt(fPop->At(ii), nreel);
527  nreel += 1;
528  }
529  }
530 
531  if (nreel != fNclusters + 1)
532  printf("pb de coherence %d %d\n", nreel, fNclusters + 1);
533 
534 }
535 
536 
537 
538 
542 
544 {
545  //return the number of clusters found by the Clusterize() method
546  //
547  return fNclusters;
548 }
549 
550 
551 
556 
558 {
559  //return the array of size (or number of cells) associated to each cluster
560  // found by the Clusterize() method
561  //
562  return fNcells->GetArray();
563 }
564 
565 
566 
571 
573 {
574  //return the array of population (or total number of contents) associated to each cluster
575  // found by the Clusterize() method
576  //
577  return fPop->GetArray();
578 }
579 
580 
581 
582 
586 
588 {
589  //print results of the Clusterize() method
590  //
591  Info("PrintResults", "%d clusters found", fNclusters);
592  for (Int_t ii = 0; ii < fNclusters; ii += 1) {
593  printf("%d %d %lf\n", ii, fNcells->At(ii), fPop->At(ii));
594  }
595 
596 }
597 
598 
599 
600 
607 
609 {
610  //return the clone 3D histogram where each cell has for content
611  //the index of the cluster it is associated to
612  //It's allow, afterward, for the user to calculate other observable associated to each
613  //clusters
614  //
615  return fCluster;
616 }
617 
618 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
char Char_t
bool Bool_t
double Double_t
float Float_t
const Bool_t kTRUE
TH3 object which allow clusterization in cell density.
Definition: KVClust3D.h:17
void init()
initialization method
Definition: KVClust3D.cpp:16
Int_t GetNclusters() const
Definition: KVClust3D.cpp:543
void PrintResults() const
Definition: KVClust3D.cpp:587
KVClust3D()
Default constructor.
Definition: KVClust3D.cpp:38
void ResetObjects()
Definition: KVClust3D.cpp:173
TH3D * ProduceTH3D(const Char_t *name) const
Definition: KVClust3D.cpp:251
Double_t * GetPop() const
Definition: KVClust3D.cpp:572
void SetDensityThreshold(Double_t)
Set density threshold on the content of a cell to be considered has filled.
Definition: KVClust3D.cpp:204
void Clusterize()
Definition: KVClust3D.cpp:285
virtual ~KVClust3D()
Definition: KVClust3D.cpp:93
Bool_t fKeepLonelyCells
tells if we keep or not lonely cells, par default kTRUE
Definition: KVClust3D.h:25
TArrayD * fPop
array where total contents of the clusters are stored
Definition: KVClust3D.h:29
Int_t fNvoisins
Definition: KVClust3D.h:26
Int_t GetNvoisins() const
Definition: KVClust3D.cpp:139
Double_t GetThreshold() const
return the value of the threshold for the content of a cell to be considered has filled
Definition: KVClust3D.cpp:217
void SetThreshold(Double_t)
Set threshold on the content of a cell to be considered has filled.
Definition: KVClust3D.cpp:191
TH2D * fVoisins
internal histogram to gather possible connected clusters
Definition: KVClust3D.h:23
void PrintInputs() const
Definition: KVClust3D.cpp:229
TArrayI * fNcells
exemple sur une dimension les cellules de nx-fNvoisins a nx+fNvoisins sont prises en compte,...
Definition: KVClust3D.h:28
void SetNvoisins(Int_t)
Definition: KVClust3D.cpp:129
void SetLonelyCells(Bool_t)
Definition: KVClust3D.cpp:109
TH3D * fCluster
3D histogram where index of each cluster is putted as content in the associated cells
Definition: KVClust3D.h:24
TH3D * GetClusterIndex() const
Definition: KVClust3D.cpp:608
Double_t GetVolumeVoisin() const
Definition: KVClust3D.cpp:159
Bool_t IsLonelyCellsAreKept() const
Definition: KVClust3D.cpp:119
Double_t GetVolumeCell() const
Definition: KVClust3D.cpp:149
Double_t fThreshold
threshold on the content of each cell to be considered has filled
Definition: KVClust3D.h:22
Int_t fNclusters
number of clusters found by the Clusterize method
Definition: KVClust3D.h:21
Int_t * GetSize() const
Definition: KVClust3D.cpp:557
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
Int_t First() const
Returns smallest number included in list.
Bool_t End(void) const
Definition: KVNumberList.h:197
Int_t GetNValues() const
void Begin(void) const
void Add(Int_t)
Add value 'n' to the list.
Int_t Next(void) const
void SetAt(Double_t v, Int_t i)
Double_t At(Int_t i) const
void Set(Int_t n)
void AddAt(Double_t c, Int_t i)
void Reset()
Double_t * GetArray()
void Set(Int_t n)
void SetAt(Double_t v, Int_t i)
void Reset()
Int_t * GetArray()
Int_t At(Int_t i) const
void AddAt(Int_t c, Int_t i)
Int_t GetSize() const
virtual Double_t GetBinWidth(Int_t bin) const
virtual Int_t GetNbinsY() const
virtual Int_t GetNbinsZ() const
TAxis * GetXaxis()
TAxis * GetYaxis()
virtual Int_t GetNbinsX() const
virtual Double_t GetBinLowEdge(Int_t bin) const
TAxis * GetZaxis()
virtual void Reset(Option_t *option="")
virtual Int_t Fill(const char *namex, const char *namey, Double_t w)
virtual Double_t GetBinContent(Int_t bin) const
virtual void Reset(Option_t *option="")
virtual Double_t GetBinContent(Int_t bin) const
virtual void SetBinContent(Int_t bin, Double_t content)
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
Double_t Min(Double_t a, Double_t b)
Double_t Power(Double_t x, Double_t y)
Double_t Max(Double_t a, Double_t b)
const double xbins[xbins_n]