KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVBreakUp.cpp
1//Created by KVClassFactory on Thu Mar 25 11:25:27 2010
2//Author: bonnet
3
4#include "KVBreakUp.h"
5#include "TMath.h"
6#include "TMethodCall.h"
7#include "TRandom3.h"
8#include "KVNucleusEvent.h"
9
10using namespace std;
11
13
14
15
16
19
20void KVBreakUp::init(void)
21{
22 //Initialisation des variables
23 //appele par le constructeur
24
25 Ztotal = 0;
26 bound = 0;
27 SetConditions(0, 0, 0);
28
29 //generateurs de nombre aleatoire
30 alea = 0;
31 RedefineTRandom("TRandom3");
32
33 BreakUpMethod = "BreakUsingChain";
34 SetBreakUpMethod(BreakUpMethod);
35
36 lhisto = new KVHashList();
37 lhisto->SetOwner(kTRUE);
38 //Rempli dans les methodes type Break...()
39 hzz = new TH1F("KVBreakUp_hzz", "distri z", 200, -0.5, 199.5);
40 lhisto->Add(hzz);
41 nraffine = 0;
42
43 lobjects = new KVHashList();
44 lobjects->SetOwner(kTRUE);
45
46 DefineHistos();
47
48 niter_tot = 0;
49 tstart = tstop = tellapsed = 0;
50
51 StorePartitions();
52 current_event = KVEvent::Factory("SimEvent");
53}
54
55
56
62
64{
65 //Constructeur, l'argument taille_max correspond
66 //a la taille maximale du tableau ou sont enregistrees
67 //les produits de la cassure (par defaut 5000)
68 //Ce qui correspond a la cassure en 5000 morceaux max
69 size_max = taille_max;
70 size = new Int_t[size_max];
71 init();
72}
73
74
75
78
80{
81 //Destructeur
82 delete alea;
83 delete lhisto;
84 delete lobjects;
85
86 //delete current_event;
87
88 if (bound) delete [] bound;
89 delete [] size;
90
91 TCanvas* c1 = 0;
92 if ((c1 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("BreakUp_Control"))) delete c1;
93
94}
95
96
101
103{
104 //Remet a zero les compteurs
105 //Efface les partitions enregistrees
106 //Le contenu des histos
109 ResetHistos();
110
111}
112
113
114
119
121{
122 //Definition des histogrammes
123 //A redefinir si besoin dans les classes filles
124 //A Remplir dans la methode TreatePartition()
125
126 hzt = new TH1F("KVBreakUp_hzt", "h ztotal", 200, -0.5, 199.5);
127 lhisto->Add(hzt);
128 hmt = new TH1F("KVBreakUp_hmt", "h mtotal", 200, -0.5, 199.5);
129 lhisto->Add(hmt);
130}
131
132
133
136
138{
139 //Protected method
140 if (zt != Ztotal) {
141 Ztotal = zt;
142 //zt correspond au nombre de possibles liens a casser
143 if (bound) delete bound;
144 bound = new Int_t[Ztotal];
145 for (Int_t nn = 0; nn < Ztotal; nn += 1) bound[nn] = 1;
147 }
148}
149
150
151
154
156{
157 //Protected method
158 if (mt > size_max) {
159 Error("SetMtot", "%d -> La multiplicite max (%d) est depassee", mt, size_max);
160 exit(EXIT_FAILURE);
161 }
162 Mtotal = mt;
164}
165
166
167
170
172{
173 //Protected method
174
175 Zmin = zlim;
177
178}
179
180
181
184
186{
187 //Protected method
188 BreakUpMethod = bup_method;
189
190}
191
192
193
195
197{
198
199 if (current_event) delete current_event;
200 current_event = evt;
201
202}
203
204
205
210
212{
213 //Permet de definir une classe pour le tirage aleatoire pour la
214 //cassure en entiers
215 //La class doit derivee de TRandom
216
217 if (alea) delete alea;
218 TClass* cl = new TClass(TRandom_Method.Data());
219 alea = (TRandom*)cl->New();
220 delete cl;
221}
222
223
224
228
230{
231 //Definition du nombre a casser (zt), en plusieurs nombres entiers (mt)
232 //avec une valeur minimale (zmin)
233 SetZtot(zt);
234 SetMtot(mt);
235 SetZmin(zmin);
236}
237
238
239
243
245{
246 //Permet de definir une methode
247 //de cassure
248
249 if (bup_method == "") {
250 Info("DefineBreakUpMethod", "Available methods are");
251 cout << "BreakUsingChain" << endl;
252 cout << "BreakUsingPile" << endl;
253 cout << "BreakUsingIndividual" << endl;
254 cout << "BreakUsingLine" << endl;
255 cout << "Make your choice" << endl;
256 } else {
257 SetBreakUpMethod(bup_method);
258 }
259
260}
261
262
263
266
268{
269 //Retourne le nombre entier a casser
270 return Ztotal;
271
272}
273
274
275
278
280{
281 //Retourne le nombre d'entiers apres cassure (la multiplicite)
282 return Mtotal;
283
284}
285
286
287
290
292{
293 //Retourne la taille minimale des entiers apres cassure
294 return Zmin;
295
296}
297
298
299
302
304{
305 //Retourne methode de cassure
306 return BreakUpMethod;
307
308}
309
310
311
314
316{
317 //si choix=kTRUE, on enregistre les partitions
318 SetBit(kStorePartitions, choix);
319
320}
321
322
323
328
330{
331 //Methode de cassure
332
333 //Conditions de depart
334 //Mtotal clusters de taille minimale Zmin
335 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
336 size[mm] = Zmin - 1;
337 }
338 //On initilise la taille a Zmin-1
339 //Le tirage aleatoire se fait
340 //sur la taille restante Ztotal_corr
341 Int_t Ztotal_corr = Ztotal - (Zmin - 1) * Mtotal;
342
343 nl = "";
344 Int_t bb = 0;
345 Int_t mtire = 0;
346 while (mtire < Mtotal) {
347 //Tirage uniform d'entier entre 0 et Ztotal_corr - 1
348 bb = TMath::Nint(alea->Uniform(0, (Ztotal_corr) - 0.5));
349 //test si ce lien a deja ete casse
350 if (bound[bb] == 1) {
351 nl.Add(bb);
352 mtire += 1;
353 bound[bb] = 0;
354 }
355 }
356
357 //La cassure a reussie
358 //on recupere les tailles additionelles
359 IntArray val = nl.GetArray();
360 Int_t mtot_corr = val.size();
361
362 Int_t zc = 0;
363 //boucle sur les liens casses
364 Int_t taille = 0;
365 for (Int_t ii = 1; ii < mtot_corr; ii += 1) {
366 taille = val[ii] - val[ii - 1];
367 size[ii] += taille; //mise a jour des tailles des clusters
368 hzz->Fill(size[ii]);
369 bound[val[ii]] = 1; //Reset du lien pour le prochain tirage
370 zc += size[ii]; //incrementation de la taille totale de controle
371 }
372 taille = Ztotal_corr - val[mtot_corr - 1] + val[0];
373
374 size[0] += taille;
375 hzz->Fill(size[0]);
376 bound[val[0]] = 1;
377 zc += size[0];
378
379 //cout << zc << " " << Ztotal << endl;
380 if (zc == Ztotal && mtot_corr == Mtotal) return 1;
381 else return 0;
382
383}
384
385
386
391
393{
394 //Methode de cassure
395
396 //Conditions de depart
397 //Mtotal clusters de taille minimale Zmin
398 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
399 size[mm] = Zmin - 1;
400 }
401 //On initilise la taille a Zmin-1
402 //Le tirage aleatoire se fait
403 //sur la taille restante Ztotal_corr
404 Int_t Ztotal_corr = Ztotal - (Zmin - 1) * Mtotal;
405
406 nl = "";
407 Int_t bb = 0;
408 Int_t mtire = 0;
409 while (mtire < Mtotal - 1) {
410 //Tirage uniform d'entier entre 0 et Ztotal_corr - 1
411 bb = TMath::Nint(alea->Uniform(1, (Ztotal_corr - 1) - 0.5));
412 //test si ce lien a deja ete casse
413 if (bound[bb] == 1) {
414 //cout << "bb="<<bb << endl;
415 nl.Add(bb);
416 mtire += 1;
417 bound[bb] = 0;
418 }
419 }
420
421 //cout << "mtire="<<mtire << endl;
422 //
423 nl.Add(0);
424 nl.Add(Ztotal_corr);
425 //cout << "nl="<<nl.AsString() << endl;;
426 //La cassure a reussie
427 //on recupere les tailles additionelles
428 IntArray val = nl.GetArray();
429 Int_t mtot_corr = val.size();
430 //cout << "mtot_corr="<<mtot_corr << endl;
431 Int_t zc = 0;
432 //boucle sur les liens casses
433 Int_t taille = 0;
434 for (Int_t ii = 1; ii < mtot_corr; ii += 1) {
435 taille = val[ii] - val[ii - 1];
436 //cout << "ii="<< ii <<" taille="<<taille<<endl;
437 size[ii - 1] += taille; //mise a jour des tailles des clusters
438 hzz->Fill(size[ii - 1]);
439 bound[val[ii]] = 1; //Reset du lien pour le prochain tirage
440 zc += size[ii - 1]; //incrementation de la taille totale de controle
441 }
442
443 //cout << "zc="<<zc << " " << Ztotal << endl;
444 if (zc == Ztotal && mtot_corr - 1 == Mtotal) return 1;
445 else return 0;
446
447}
448
449
450
455
457{
458 //Methode de cassure
459
460 //Conditions de depart
461 // Mtotal clusters de taille minimale Zmin
462 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
463 size[mm] = Zmin;
464 }
465 //On initilise la taille a Zmin
466 //Le tirage aleatoire se fait
467 //sur la taille restante Ztotal_corr
468 Int_t Ztotal_corr = Ztotal - (Zmin) * Mtotal;
469
470 Double_t val[Mtotal], sum_val = 0;
471 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
472 //Tirage des Mtotal valeurs entre [0; 1[
473 val[mm] = alea->Uniform(0, 1);
474 sum_val += val[mm];
475 }
476 Int_t zc = 0;
477 Int_t surplus[Mtotal];
478 //normalisation a une charge totale de Ztotal_corr
479 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
480 surplus[mm] = TMath::Nint(val[mm] * Ztotal_corr / sum_val);
481 zc += surplus[mm];
482 }
483 //Test effets spurieux pasage reel -> entier
484 //Redistribution eventuelle de charge manquante ou en exces
485 if (zc != Ztotal_corr) {
486 Int_t diff = Ztotal_corr - zc;
487 nraffine += 1;
488 Int_t bb = 0;
489 if (diff > 0) {
490 for (Int_t mm = 0; mm < TMath::Abs(diff); mm += 1) {
491 bb = TMath::Nint(alea->Uniform(0, (Mtotal) - 0.5));
492 surplus[bb] += 1;
493 }
494 } else {
495 for (Int_t mm = 0; mm < TMath::Abs(diff); mm += 1) {
496 bb = TMath::Nint(alea->Uniform(0, (Mtotal) - 0.5));
497 if (surplus[bb] > 0)
498 surplus[bb] -= 1;
499 else
500 mm -= 1;
501 }
502 }
503 }
504 zc = 0;
505 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
506 size[mm] += surplus[mm];
507 zc += size[mm];
508 hzz->Fill(size[mm]);
509 }
510
511 if (zc == Ztotal) return 1;
512 else return 0;
513
514}
515
516
517
520
522{
523 //Methode de cassure
524 Int_t bb = 0;
525
526 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
527 size[mm] = Zmin;
528 }
529 //on distribue aleatoirement les charges
530 //disponibles une par une sur les differents fragments
531 for (Int_t mm = 0; mm < nbre_nuc; mm += 1) {
532 bb = TMath::Nint(alea->Uniform(0, (Mtotal) - 0.5));
533 size[bb] += 1;
534 }
535
536 Int_t mc = 0, zc = 0;
537 for (Int_t mm = 0; mm < Mtotal; mm += 1) {
538 Int_t taille = size[mm];
539 hzz->Fill(taille);
540 mc += 1;
541 zc += taille;
542 }
543 if (mc == Mtotal && zc == Ztotal) return 1;
544 else return 0;
545
546
547}
548
549
550
558
560{
561 //Remplissage des histogrammes predefinis
562 //A redefinir si besoin dans les classes filles
563 //
564 //Enregistrement de la partition (si demande, via StorePartitions(kTRUE) par defaut)
565 //Si une partition identique est deja presente, on incremente sa population
566 //sinon on enregistre celle-ci, voir KVPartitionList
567 hmt->Fill(Mtotal);
568 hzt->Fill(Ztotal);
569
570
571 partition = new KVIntegerList();
574
576 for (Int_t ii = 0; ii < Mtotal; ii += 1)
577 current_event->AddNucleus()->SetZ(tab[ii]);
578
579 delete [] tab;
580
582 if (Fill(partition))
583 delete partition;
584 } else {
585 delete partition;
586 }
587}
588
589
590
591
593
595{
596
597 return current_event;
598
599}
600
601
602
610
612{
613 //On realise times fois la cassure
614 //suivant les conditions definis vis la methode SetConditions et DefineBreakUpMethod
615 //Plusieurs series de cassures peuvent etre ainsi realise
616 //
617 //Si l utilisateur veut changer de conditions et ne pas melanger les partitions
618 //il faut appeler la methode Clear()
619
620 Start();
621
622 TMethodCall meth;
623 meth.InitWithPrototype(this->IsA(), BreakUpMethod.Data(), "");
624 Long_t ret;
625 if (meth.IsValid() && meth.ReturnType() == TMethodCall::kLong) {
626 for (Int_t nn = 0; nn < times; nn += 1) {
627 if (nn % 10000 == 0) Info("BreakNtimes", "%d partitions generees sur %d", nn, times);
628 meth.Execute(this, "", ret);
629 if (ret == 1)
631 else {
632 Info("BreakNtimes", "%s retourne %ld", BreakUpMethod.Data(), ret);
633 nn -= 1;
634 }
635 }
636 }
637 Info("BreakNtimes", "Tirage termine");
638 niter_tot += times;
639
640 Stop();
641 Info("BreakNtimes", "Temps ecoule en secondes : %d", GetDeltaTime());
642
643}
644
645
646
657
658void KVBreakUp::BreakNtimesOnGaussian(Int_t times, Double_t Ztot_moy, Double_t Ztot_rms, Double_t Mtot_moy, Double_t Mtot_rms, Int_t zmin)
659{
660 //On realise times fois la cassure sur une double gaussienne
661 //la methode SetCondition n'a ici aucune influence
662 //les valeurs ztot et mtot sont tirees aleatoirement sur une gaussienne
663 //a chaque iteration
664 //
665 //Plusieurs series de cassures peuvent etre ainsi realise
666 //
667 //Si l utilisateur veut changer de conditions et ne pas melanger les partitions
668 //il faut appeler la methode Clear()
669
670 Start();
671
672 TMethodCall meth;
673 meth.InitWithPrototype(this->IsA(), BreakUpMethod.Data(), "");
674 Long_t ret;
675 TRandom3* gaus = new TRandom3();
676 if (meth.IsValid() && meth.ReturnType() == TMethodCall::kLong) {
677 for (Int_t nn = 0; nn < times; nn += 1) {
678 if (nn % 1000 == 0) printf("%d partitions generees sur %d\n", nn, times);
679
680 Int_t zt = TMath::Nint(gaus->Gaus(Ztot_moy, Ztot_rms));
681 Int_t mt = TMath::Nint(gaus->Gaus(Mtot_moy, Mtot_rms));
682
683 if (mt > 0 && zt > 0 && zt >= mt * zmin) {
684 SetConditions(zt, mt, zmin);
685 meth.Execute(this, "", ret);
686 if (ret == 1)
688 else {
689 //cout << BreakUpMethod << " retourne " << ret << endl;
690 }
691 } else {
692 nn -= 1;
693 }
694 }
695 }
696 delete gaus;
697 niter_tot += times;
698
699 Stop();
700 Info("BreakNtimesOnGaussian", "Temps ecoule en secondes : %d", GetDeltaTime());
701
702}
703
704
705
713
714KVEvent* KVBreakUp::BreakOnGaussian(Double_t Ztot_moy, Double_t Ztot_rms, Double_t Mtot_moy, Double_t Mtot_rms, Int_t zmin)
715{
716 //On realise la cassure sur une double gaussienne
717 //
718 //les valeurs ztot et mtot sont tirees aleatoirement sur une gaussienne
719 //
720 //Plusieurs series de cassures peuvent etre ainsi realise
721 //
722
723 TMethodCall meth;
724 meth.InitWithPrototype(this->IsA(), BreakUpMethod.Data(), "");
725 Long_t ret;
726
727 if (meth.IsValid() && meth.ReturnType() == TMethodCall::kLong) {
728 Int_t zt = -1;
729 Int_t mt = -1;
730 while (!(mt > 0 && zt > 0 && zt >= mt * zmin)) {
731 zt = TMath::Nint(alea->Gaus(Ztot_moy, Ztot_rms));
732 mt = TMath::Nint(alea->Gaus(Mtot_moy, Mtot_rms));
733 }
734
735 SetConditions(zt, mt, zmin);
736 meth.Execute(this, "", ret);
737 if (ret == 1) {
739 return current_event;
740 } else {
741 return 0;
742 }
743 }
744 return 0;
745}
746
747
748
762
763void KVBreakUp::BreakFromHisto(TH2F* hh_zt_VS_mt, Int_t zmin)
764{
765 //On realise times fois la cassure suivant un histogramme a deux dimensions
766 //avec en definition :
767 // mtot -> axe X
768 // ztot -> axe Y
769 //la methode SetCondition n'a ici aucune influence
770 //les valeurs ztot et mtot sont tirees aleatoirement sur cette histogramme
771 //a chaque iteration
772 //
773 //Plusieurs series de cassures peuvent etre ainsi realise
774 //
775 //Si l utilisateur veut changer de conditions et ne pas melanger les partitions
776 //il faut appeler la methode Clear()
777
778 TH2F* h2 = hh_zt_VS_mt;
779 if (!h2) return;
780
781 Start();
782
783 TMethodCall meth;
784 meth.InitWithPrototype(this->IsA(), BreakUpMethod.Data(), "");
785 Long_t ret;
786
787 Int_t zt, mt;
788 Int_t stat_tot = Int_t(h2->Integral());
789 Int_t stat_par = 0;
790 for (Int_t nx = 1; nx <= h2->GetNbinsX(); nx += 1)
791 for (Int_t ny = 1; ny <= h2->GetNbinsY(); ny += 1) {
792 Int_t stat = TMath::Nint(h2->GetBinContent(nx, ny));
793 if (stat > 0) {
794 mt = TMath::Nint(h2->GetXaxis()->GetBinCenter(nx));
795 zt = TMath::Nint(h2->GetYaxis()->GetBinCenter(ny));
796 if (mt > 0 && zt > 0 && zt >= mt * zmin) {
797
798 SetConditions(zt, mt, zmin);
799 for (Int_t nn = 0; nn < stat; nn += 1) {
800 meth.Execute(this, "", ret);
801 if (ret == 1)
803 stat_par += 1;
804 if (stat_par % 1000 == 0) printf("%d partitions generees sur %d\n", stat_par, stat_tot);
805 }
806 } else {
807 cout << zt << " " << mt << endl;
808 }
809 }
810 }
811
812 niter_tot += stat_par;
813
814 Stop();
815 Info("BreakNtimesFromHisto", "Temps ecoule en secondes : %d", GetDeltaTime());
816
817}
818
819
820
824
826{
827 //Retourne le nombre total d iterations
828 //depuis le dernier clear
829 return niter_tot;
830
831}
832
833
834
839
841{
842 //Retourne la liste des histogrammes
843 //definis dans DefineHistos()
844 //si l utilisateur a utilise lhisto->Add(TObject * )
845 return lhisto;
846
847}
848
849
850
853
855{
856 //Retourne la liste d'objects autres ...
857 return lobjects;
858
859}
860
861
862
865
867{
868 //Comme c'est ىcrit
869 niter_tot = 0;
870
871}
872
873
874
877
879{
880 //Met a zero le contenu des histogrammes
881 lhisto->Execute("Reset", "");
882
883}
884
885
886
890
892{
893 //Trace les histos definis
894 //A redefinir si besoin dans les classes filles
895
896 TCanvas* c1 = 0;
897 if (!(c1 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("BreakUp_Control"))) {
898
899 c1 = new TCanvas("BreakUp_Control", "Control", 0, 0, 900, 900);
900 c1->Divide(2, 2, 0.001, 0.001, 10);
901 c1->cd(1);
902 if (hzz->GetEntries() > 0) gPad->SetLogy(1);
903 hzz->Draw();
904 c1->cd(2);
905 if (hzt->GetEntries() > 0) gPad->SetLogy(1);
906 hzt->Draw();
907 c1->cd(3);
908 if (hmt->GetEntries() > 0) gPad->SetLogy(1);
909 hmt->Draw();
910
911 }
912}
913
914
915
926
927void KVBreakUp::SaveHistos(KVString filename, KVString suff, Option_t* option)
928{
929 //Permet la sauvegarde des histogrammes
930 //dans un fichier (option est l'option classique de TFile::TFile())
931 //
932 //Si filename=="" (defaut) -> nom du fichier = "KVBreakUp_Ouput.root";
933 //Si suff=="" (defaut) -> definition d'un suffixe en utisant les arguments de SetConditions()
934 //Ce suffixe est ensuite utilise dans la redefinition des noms des histogrammes
935 //Exemple:
936 // nom de l'histogramme : histo -> histo_suff
937 //
938 if (filename == "") filename = "KVBreakUp_Ouput.root";
939 TFile* file = new TFile(filename.Data(), option);
940
941 if (suff == "")
942 suff.Form("Zt%d_Mt%d_Zm%d_%s", GetZtot(), GetMtot(), GetZmin(), GetBreakUpMethod().Data());
943
944 KVString snom;
945 for (Int_t nn = 0; nn < lhisto->GetEntries(); nn += 1) {
946 snom.Form("%s_%s", lhisto->At(nn)->GetName(), suff.Data());
947 lhisto->At(nn)->Write(snom.Data());
948 }
949 file->Close();
950
951}
952
953
954
960
962{
963 // Comme c'est ىcrit
964 // Why not "Print(Option_t*)" ?
965 // - because TCollection has several 'virtual' Print methods which
966 // are overloaded (i.e. have different arguments): BAD!
967
968 Info("Print", "Configuration for the break up");
969 printf(" Ztot=%d - Mtot=%d - Zmin=%d\n", GetZtot(), GetMtot(), GetZmin());
970 printf(" Charge to be distributed %d - Biggest possible charge %d\n", nbre_nuc, Zmin + nbre_nuc);
971 printf(" Methode de cassage aleatoire %s\n", GetBreakUpMethod().Data());
972 alea->Print();
973 printf(" Partition are stored via KVPartitionList : %d\n", Int_t(TestBit(kStorePartitions)));
974 printf("------------------------------------------------------");
975}
976
977
978
982
984{
985 //protected method
986 //Signal start
987 TDatime time;
988 tstart = time.GetHour() * 3600 + time.GetMinute() * 60 + time.GetSecond();
989
990}
991
992
993
997
999{
1000 //protected method
1001 //Signal stop
1002 TDatime time;
1003 tstop = time.GetHour() * 3600 + time.GetMinute() * 60 + time.GetSecond();
1005
1006}
1007
1008
1009
1013
1015{
1016 //Retoune le temps ecoules (en seconde)
1017 //entre un appel Start() et un appel Stop()
1018 return tellapsed;
1019
1020}
1021
1022
int Int_t
long Long_t
bool Bool_t
double Double_t
constexpr Bool_t kTRUE
const char Option_t
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
#define gROOT
#define gPad
Permet de casser aleatoirement un nombre entier (ztot) en un nombre (mtot) d'entiers plus petits d'un...
Definition KVBreakUp.h:66
void BreakFromHisto(TH2F *hh_zt_VS_mt, Int_t zmin=1)
Int_t Ztotal
Definition KVBreakUp.h:70
virtual void SaveHistos(KVString filename="", KVString suff="", Option_t *option="recreate")
Int_t * size
[size_max]->
Definition KVBreakUp.h:81
Int_t GetTotalIterations(void)
KVString BreakUpMethod
Definition KVBreakUp.h:74
KVEvent * current_event
Definition KVBreakUp.h:92
Int_t GetDeltaTime()
virtual void DrawPanel()
KVEvent * BreakOnGaussian(Double_t Ztot_moy, Double_t Ztot_rms, Double_t Mtot_moy, Double_t Mtot_rms, Int_t zmin)
Int_t Mtotal
Definition KVBreakUp.h:71
Int_t tstop
Definition KVBreakUp.h:83
TH1F * hzt
Definition KVBreakUp.h:86
Int_t tellapsed
Definition KVBreakUp.h:83
virtual void ResetHistos()
Met a zero le contenu des histogrammes.
void DefineBreakUpMethod(KVString bup_method="")
virtual void TreatePartition()
Int_t GetZtot(void) const
Retourne le nombre entier a casser.
TH1F * hzz
Definition KVBreakUp.h:85
void Clear(Option_t *="")
Int_t BreakUsingChain()
Int_t BreakUsingLine()
KVHashList * GetObjects()
Retourne la liste d'objects autres ...
virtual ~KVBreakUp()
Destructeur.
Definition KVBreakUp.cpp:79
TH1F * hmt
Definition KVBreakUp.h:87
KVIntegerList * partition
Definition KVBreakUp.h:91
void init()
Definition KVBreakUp.cpp:20
void SetZmin(Int_t zlim)
Protected method.
KVString GetBreakUpMethod(void) const
Retourne methode de cassure.
Int_t BreakUsingPile()
Methode de cassure.
void Stop()
Int_t nbre_nuc
Definition KVBreakUp.h:79
Int_t BreakUsingIndividual()
KVHashList * lhisto
Definition KVBreakUp.h:90
void Start()
TRandom * alea
Definition KVBreakUp.h:76
void StorePartitions(Bool_t choix=kTRUE)
si choix=kTRUE, on enregistre les partitions
void LinkEvent(KVEvent *)
void BreakNtimes(Int_t times=1000)
KVNumberList nl
Definition KVBreakUp.h:78
KVEvent * GetCurrentEvent()
void PrintConfig() const
void BreakNtimesOnGaussian(Int_t times, Double_t Ztot_moy, Double_t Ztot_rms, Double_t Mtot_moy, Double_t Mtot_rms, Int_t zmin=1)
KVBreakUp(Int_t taille_max=1000)
Definition KVBreakUp.cpp:63
void SetZtot(Int_t zt)
Protected method.
Int_t Zmin
Definition KVBreakUp.h:72
Int_t GetZmin(void) const
Retourne la taille minimale des entiers apres cassure.
KVHashList * GetHistos()
void SetConditions(Int_t zt, Int_t mt, Int_t zmin=1)
Int_t tstart
Definition KVBreakUp.h:83
void RedefineTRandom(KVString TRandom_Method)
@ kStorePartitions
Definition KVBreakUp.h:95
Int_t nraffine
Definition KVBreakUp.h:109
Int_t size_max
Definition KVBreakUp.h:80
Int_t niter_tot
Definition KVBreakUp.h:82
KVHashList * lobjects
Definition KVBreakUp.h:89
KVString TRandom_Method
Definition KVBreakUp.h:74
virtual void ResetTotalIterations()
Comme c'est ىcrit.
Int_t * bound
[Ztotal] tableau permettant de gérer les cassures de liens
Definition KVBreakUp.h:73
Int_t GetMtot(void) const
Retourne le nombre d'entiers apres cassure (la multiplicite)
virtual void DefineHistos()
void SetMtot(Int_t mt)
Protected method.
void SetBreakUpMethod(KVString bup_method)
Protected method.
Abstract base class container for multi-particle events.
Definition KVEvent.h:67
void Clear(Option_t *opt="")
Definition KVEvent.h:238
KVNucleus * AddNucleus()
Definition KVEvent.cpp:109
static KVEvent * Factory(const char *plugin)
Definition KVEvent.h:228
Extended version of ROOT THashList.
Definition KVHashList.h:29
Handle a list of positive integers (partition)
Int_t * CreateTableOfValues()
void Fill(Int_t *tab, Int_t mult)
void Fill(TArrayI* tab);
void SetZ(Int_t z, Char_t mt=-1)
void Add(Int_t)
Add value 'n' to the list.
IntArray GetArray() const
virtual Bool_t Fill(KVIntegerList *par)
void Clear(Option_t *option="")
Mise a zero de la liste.
virtual void Execute(const char *method, const char *params, Int_t *error=0)
virtual TObject * At(Int_t idx) const
virtual void Add(TObject *obj)
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73
virtual Double_t GetBinCenter(Int_t bin) const
virtual Int_t GetEntries() const
Int_t GetHour() const
Int_t GetSecond() const
Int_t GetMinute() const
virtual Int_t GetNbinsY() const
TAxis * GetXaxis()
virtual Int_t GetNbinsX() const
TAxis * GetYaxis()
void Draw(Option_t *option="") override
virtual Int_t Fill(const char *name, Double_t w)
virtual Double_t GetEntries() const
virtual Double_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Option_t *option="") const
virtual Double_t GetBinContent(Int_t bin) const
EReturnType ReturnType()
static const EReturnType kLong
void InitWithPrototype(const char *function, const char *proto, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Bool_t IsValid() const
void Execute()
void Print(Option_t *option="") const override
virtual const char * GetName() const
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
virtual Double_t Uniform(Double_t x1, Double_t x2)
TClass * IsA() const override
const char * Data() const
void Form(const char *fmt,...)
return c1
Int_t Nint(T x)
Double_t Abs(Double_t d)
ClassImp(TPyArg)