KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
MCSampler.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Thu May 7 14:35:05 2015
2 //Author: John Frankland,,,
3 
4 #include "MCSampler.h"
5 #include "TRandom.h"
6 
7 #include <TGraph.h>
8 #include <TMultiGraph.h>
9 #include "TPad.h"
10 #include "TH1.h"
11 
13 
14 
15 namespace MicroStat {
16 
17 
20 
21  void MCSampler::init()
22  {
23  // default initialisations
24  fWeightList = 0;
25  fLastPicked = 0;
26  fTheLegend = 0;
27  fLegendProbaMin = 0;
28  fModifyMasses = kFALSE;
29  }
30 
31 
32 
35 
36  void MCSampler::initialiseWeightList()
37  {
38  // set up the TClonesArray
39 
40  fWeightList = new TClonesArray(fWeight);
41  }
42 
43 
44 
47 
48  MCSampler::MCSampler()
49  {
50  // Default constructor
51  init();
52  }
53 
54 
55 
56 
59 
60  MCSampler::MCSampler(const Char_t* name, const Char_t* title) : KVBase(name, title)
61  {
62  // Write your code here
63  init();
64  }
65 
66 
67 
70 
72  {
73  // Destructor
74 
76  }
77 
78 
79 
82 
83  void MCSampler::SetEventList(TTree* t, const TString& branchname)
84  {
85  // Define the TTree or TChain containing all possible events (partitions).
86 
87  fPartitions = t->GetEntries();
88  fBranch = t->GetBranch(branchname);
89  if (!fBranch) {
90  Error("SetEventList", "cannot find branch %s", branchname.Data());
91  }
92  fPartition = 0;
94  }
95 
96 
97 
101 
103  {
104  // Set the kind of statistical weight to be used
105  // This is the name of a class derived from MicroStat::StatWeight
106 
108  if (!fWeight) {
109  Error("SetStatWeight", "class %s not found", w.Data());
110  }
111  }
112 
113 
114 
118 
120  {
121  // if nuclear masses are modified, we have to update those in the
122  // partition read from file
123  KVNucleus* n;
124  while ((n = fPartition->GetNextParticle())) n->SetA(n->GetA());
125  }
126 
127 
128 
133 
134  void MCSampler::CalculateWeights(Double_t excitation_energy)
135  {
136  // calculate weights of all partitions for the given excitation energy
137  // (in MeV) of the initial compound nucleus
138 
139  //Info("CalculateWeights","Calculating channel weights for E*=%f",excitation_energy);
140 
142 
143  fSumWeights = 0;
144 
145  for (int i = 0; i < fPartitions; i++) {
147  GetPartition(i);
149  w->SetWeight(fPartition, excitation_energy + fPartition->GetChannelQValue());
150  w->SetIndex(i);
151  fSumWeights += w->GetWeight();
152 
153  //std::cout << i << " Q=" << fPartition->GetChannelQValue() << " weight=" << w->GetWeight() << endl;
154  }
155 
156  // sort weights in decreasing order
157  fWeightList->Sort();
158  }
159 
160 
161 
169 
171  {
172  // Return the TTree index of a random channel according to the
173  // previously calculated weights
174  // The corresponding weight for the chosen channel can be retrieved
175  // by calling GetRandomChannelWeight()
176  // If no channel is open (i.e. all weights = 0, E* < Q value of first channel),
177  // we return -1.
178 
180 
181  int i;
182  for (i = 0; i < fPartitions; i++) {
183  Double_t w = ((StatWeight*)(*fWeightList)[i])->GetWeight();
184  if (x < w) break;
185  x -= w;
186  }
187  if (i == fPartitions) {
188  fLastPicked = nullptr;
189  return -1;
190  }
192  return fLastPicked->GetIndex();
193  }
194 
195 
196 
198 
199  void MCSampler::SetBranch(TTree* theTree, const TString& bname, void* variable, const TString& vartype)
200  {
201  TString leaflist = bname + "/";
202  leaflist += vartype;
203  TBranch* b = theTree->GetBranch(bname);
204  if (!b) theTree->Branch(bname, variable, leaflist);
205  else b->SetAddress(variable);
206  }
207 
208 
209 
219 
220  void MCSampler::SetUpTreeBranches(KVEvent*& event, TTree* theTree, const TString& bname)
221  {
222  // If branch 'bname' does not exist in 'theTree', it will be created and linked to
223  // 'event' which must contain the address of a valid KVEvent object.
224  // If branch already exists, we link it to 'event'.
225  //
226  // Branches will also be created/filled with the following information:
227  // ESTAR/D the excitation energy (Exx)
228  // EDISP/D the available kinetic energy
229  // IPART/I the partition index
230  TBranch* b = theTree->GetBranch(bname);
231  if (!b) {
232  theTree->Branch(bname, "KVEvent", event, 10000000, 0)->SetAutoDelete(kFALSE);
233  }
234  else {
235  b->SetAddress(event);
236  }
237  SetBranch(theTree, "ESTAR", &ESTAR, "D");
238  SetBranch(theTree, "EDISP", &EDISP, "D");
239  SetBranch(theTree, "IPART", &IPART, "L");
240  }
241 
242 
243 
252 
253  void MCSampler::GenerateEvents(TTree* theTree, KVEvent* event, Double_t Exx, Long64_t npartitions, Long64_t nev_part)
254  {
255  // Generate (npartitions*nev_part) events for a given value of excitation energy 'Exx' [MeV].
256  // For efficiency, for each chosen partition we generate nev_part events.
257  // This means :
258  // - calculating the statistical weight of all available decay channels/partitions
259  // - picking a channel at random
260  // - generating momenta of all nuclei in chosen channel
261  // - filling the TTree with the new event
262 
263  Info("GenerateEvents", "Generating events for E*=%f", Exx);
264 
265  if (!SetExcitationEnergy(Exx)) {
266  Error("GenerateEvents", "Excitation energy is too low, no channels are open");
267  return;
268  }
269 
270  // generate events
271  while (npartitions--) {
272 
273  SetDecayChannel();
274 
275 
276  for (Long64_t iev = 0; iev < nev_part; iev++) {
277 
278  // generate nev_part events for chosen partition
279 
280  GenerateEvent(theTree, event);
281 
282  }
283 
284  }
285 
286  }
287 
288 
289 
294 
296  {
297  // Define excitation energy for random event generation
298  // This will calculate the channel weights for all partitions
299  // If sumOfWeights = 0 i.e. no channels are open, return kFALSE
300 
301  ESTAR = Exx;
302  CalculateWeights(Exx);
303  if (fSumWeights == 0) return kFALSE;
304  return kTRUE;
305  }
306 
307 
312 
314  {
315  // Pick a random decay channel with excitation energy given to
316  // previous call to SetExcitationEnergy method.
317  // Returns kFALSE if no channel is open.
318 
320  if (IPART < 0) return kFALSE;
321 
323 
326 
327  return kTRUE;
328 
329  }
330 
331 
341 
343  {
344  // Generate 1 kinematical event for the partition previously
345  // picked at random after calling:
346  // SetExcitationEnergy(E*);
347  // SetDecayChannel();
348  // theTree->Fill() is automatically called.
349  // This method can be called many times before either
350  // 1) picking a new decay channel
351  // 2) changing the energy
352 
354 
355  theTree->Fill();
356 
358 
359  }
360 
361 
362 
372 
373  void MCSampler::PlotProbabilities(double emin, double emax, double estep, Option_t* opt)
374  {
375  // Plot probability of each channel as a function of E* (default)
376  // or E*/A (opt="E*/A")
377  //
378  // If SetLegendProbaMin has previously been called with a >0 value
379  // for the minimum probability, we generate a TLegend which contains
380  // only the channels whose probability reaches at least this minimum
381  // value in the given energy range. Call ShowLegend immediately
382  // after this method in order to display in current pad.
383 
384  Bool_t makeLegend = (fLegendProbaMin > 0.0);
385 
386  TMultiGraph* mg = new TMultiGraph();
387  TGraph* g[fPartitions];
388  Long64_t i = 0;
389  int mark = 20;
390  Double_t fac = 1.;
391 
392  for (; i < fPartitions; i++) {
393  KVEvent* pt = GetPartition(i);
394  if (!i && !strcmp(opt, "E*/A")) fac = pt->GetSum("GetA");
395  g[i] = new TGraph();
396  g[i]->SetName(pt->GetPartitionName());
397  g[i]->SetMarkerStyle(mark);
398  g[i]->SetMarkerColor((i % 9) + 1);
399  g[i]->SetLineColor((i % 9) + 1);
400  g[i]->SetFillColor(0);
401  g[i]->ResetBit(BIT(20));
402  mg->Add(g[i], "pl");
403 
404  if (!(i % 9)) {
405  mark++;
406  if (mark > 30) mark = 20;
407  }
408 
409  }
410  // build TGraph
411  for (double E = emin; E <= emax; E += estep) {
412 
413  CalculateWeights(E * fac);
414  for (i = 0; i < fPartitions; i++) {
415 
416  if (GetSumWeights() > 0.0) {
417  StatWeight* wt = GetWeight(i);
418 
419  Double_t proba = wt->GetWeight() / GetSumWeights();
420  Long64_t voie = wt->GetIndex();
421 
422  g[voie]->SetPoint(g[voie]->GetN(), E, proba * 100.);
423  // if we are building a TLegend, we 'mark' each TGraph which
424  // contains probabilities higher than the required minimum
425  // by setting a bit in the TObject status bitfield
426  if (makeLegend && proba > fLegendProbaMin) g[voie]->SetBit(BIT(20));
427  }
428  }
429  }
430  if (gPad) gPad->Clear();
431  mg->Draw("apl");
432  if (!strcmp(opt, "E*/A")) {
433  mg->GetHistogram()->SetXTitle("E*/A (MeV)");
434  }
435 
436 
438 
439  else {
440  mg->GetHistogram()->SetXTitle("E* (MeV)");
441  }
442 
443 
445 
446  mg->GetHistogram()->SetYTitle("Probability");
447  if (makeLegend) {
448  mg->GetHistogram()->SetAxisRange(fLegendProbaMin * 100, 100., "Y");
449  // now add all sufficiently probable decay channels to the TLegend
450  fTheLegend = new TLegend(.7, 0.12, .88, .88);
451  fTheLegend->SetHeader(Form("Channels with P>%4.1f%%", fLegendProbaMin * 100));
452  for (i = 0; i < mg->GetListOfGraphs()->GetSize(); i++) {
453  if (g[i]->TestBit(BIT(20))) {
454  fTheLegend->AddEntry(g[i], g[i]->GetName(), "pl");
455  }
456  }
457  }
458 
459  gPad->Modified();
460  gPad->Update();
461 
462  }
463 
464 
469 
470  void MCSampler::PlotMultiplicities(double emin, double emax, double estep, Option_t* opt)
471  {
472  // Plot multiplicities of n,p,d,t,3He,4He,... as a function of E* (default)
473  // or E*/A (opt="E*/A")
474  //
475 
476  TMultiGraph* mg = new TMultiGraph();
477  const int nparticles = 7;
478  TGraph* g[nparticles];
479 
480  TString particles[] = {"1n", "1H", "2H", "3H", "3He", "4He", "Z>2"};
481 
482  int mark = 20;
483  Double_t fac = 1.;
484 
485  if (!strcmp(opt, "E*/A")) {
486  KVEvent* pt = GetPartition(0);
487  fac = pt->GetSum("GetA");
488  }
489 
490 
491 
493 
494  for (int i = 0; i < nparticles; i++) {
495  g[i] = new TGraph();
496  g[i]->SetName(particles[i]);
497  g[i]->SetMarkerStyle(mark);
498  g[i]->SetMarkerColor((i % 9) + 1);
499  g[i]->SetLineColor((i % 9) + 1);
500  g[i]->SetFillColor(0);
501  g[i]->ResetBit(BIT(20));
502  mg->Add(g[i], "pl");
503 
504  if (!(i % 9)) {
505  mark++;
506  if (mark > 30) mark = 20;
507  }
508 
509  }
510 
511  // build TGraph
512 
514 
515  for (double E = emin; E <= emax; E += estep) {
516 
517  CalculateWeights(E * fac);
518  Double_t multiplicities[nparticles];
519  memset(multiplicities, 0, sizeof(double)*nparticles);
520 
521  for (Long64_t i = 0; i < fPartitions; i++) {
522 
523  if (GetSumWeights() > 0.0) {
524  StatWeight* wt = GetWeight(i);
525 
526  Double_t proba = wt->GetWeight() / GetSumWeights();
527  if (proba > 1.e-06) {
528  KVEvent* pt = GetPartition(wt->GetIndex());
529  KVNucleus* n;
530  while ((n = pt->GetNextParticle())) {
531  for (int j = 0; j < nparticles - 1; j++) {
532  if (!strcmp(n->GetSymbol(), particles[j].Data())) multiplicities[j] += proba;
533  }
534  if (n->GetZ() > 2) multiplicities[nparticles - 1] += proba;
535  }
536  }
537 
538  }
539  }
540  for (int i = 0; i < nparticles; i++) {
541  g[i]->SetPoint(g[i]->GetN(), E, multiplicities[i]);
542  }
543  }
544 
545 
547 
548  if (gPad) gPad->Clear();
549  mg->Draw("apl");
550  if (!strcmp(opt, "E*/A")) {
551  mg->GetHistogram()->SetXTitle("E*/A (MeV)");
552  }
553  else {
554  mg->GetHistogram()->SetXTitle("E* (MeV)");
555  }
556 
557  mg->GetHistogram()->SetYTitle("Multiplicity");
558  gPad->Modified();
559  gPad->Update();
560 
561  }
562 
563 }
564 
565 
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define SafeDelete(p)
#define b(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
#define BIT(n)
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
#define gPad
Base class for KaliVeda framework.
Definition: KVBase.h:135
Base class container for multi-particle events.
Definition: KVEvent.h:176
Double_t GetChannelQValue() const
Definition: KVEvent.cpp:1044
KVNucleus * GetNextParticle(Option_t *opt="") const
Definition: KVEvent.cpp:564
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Monte-Carlo sampling of events with statistical weights .
Definition: MCSampler.h:21
virtual ~MCSampler()
Destructor.
Definition: MCSampler.cpp:71
Bool_t fModifyMasses
the partition index
Definition: MCSampler.h:42
KVEvent * fPartition
branch containing events
Definition: MCSampler.h:33
void PlotProbabilities(double emin=0., double emax=100., double estep=1., Option_t *opt="")
Definition: MCSampler.cpp:373
Double_t fLegendProbaMin
weight of channel picked by call to PickRandomChannel()
Definition: MCSampler.h:25
Bool_t SetDecayChannel()
Definition: MCSampler.cpp:313
void PlotMultiplicities(double emin=0., double emax=100., double estep=1., Option_t *opt="")
Definition: MCSampler.cpp:470
void SetBranch(TTree *theTree, const TString &name, void *variable, const TString &vartype)
automatically generated legend for PlotProbabilities
Definition: MCSampler.cpp:199
Double_t ESTAR
variables for TTree branches
Definition: MCSampler.h:38
TLegend * fTheLegend
minimum probability for which channels are included in automatically generated TLegend when PlotProba...
Definition: MCSampler.h:26
void initialiseWeightList()
if nuclear masses are modified
Definition: MCSampler.cpp:36
StatWeight * GetWeight(Int_t i) const
Definition: MCSampler.h:72
void SetUpTreeBranches(KVEvent *&event, TTree *theTree, const TString &bname)
Definition: MCSampler.cpp:220
Double_t EDISP
the excitation energy (Exx)
Definition: MCSampler.h:39
Long64_t IPART
the available kinetic energy
Definition: MCSampler.h:40
void SetEventList(TTree *t, const TString &branchname)
Define the TTree or TChain containing all possible events (partitions).
Definition: MCSampler.cpp:83
StatWeight * fLastPicked
Definition: MCSampler.h:24
TClonesArray * fWeightList
statistical weight class
Definition: MCSampler.h:35
void GenerateEvents(TTree *, KVEvent *event, Double_t, Long64_t npartitions, Long64_t nev_part=10)
Definition: MCSampler.cpp:253
Bool_t SetExcitationEnergy(Double_t Exx)
Definition: MCSampler.cpp:295
Double_t fSumWeights
list of weights for all events
Definition: MCSampler.h:36
Long64_t PickRandomChannel()
Definition: MCSampler.cpp:170
void CalculateWeights(Double_t excitation_energy)
Definition: MCSampler.cpp:134
KVEvent * GetPartition(Long64_t i)
Definition: MCSampler.h:54
void init()
default initialisations
Definition: MCSampler.cpp:21
void SetStatWeight(const TString &)
Definition: MCSampler.cpp:102
Long64_t fPartitions
Definition: MCSampler.h:31
void GenerateEvent(TTree *theTree, KVEvent *event)
Definition: MCSampler.cpp:342
TBranch * fBranch
number of partitions in TTree/TChain
Definition: MCSampler.h:32
TClass * fWeight
event read from fPartitionList tree
Definition: MCSampler.h:34
Double_t GetSumWeights() const
Definition: MCSampler.h:76
Abstract base class for calculating statistical weights for events .
Definition: StatWeight.h:17
void SetIndex(Long64_t i)
Definition: StatWeight.h:49
virtual void SetWeight(KVEvent *e, Double_t E)=0
void GenerateEvent(KVEvent *partition, KVEvent *event)
Definition: StatWeight.cpp:78
virtual void initGenerateEvent(KVEvent *partition)=0
Double_t GetAvailableEnergy() const
Definition: StatWeight.h:44
virtual void resetGenerateEvent()=0
Double_t GetWeight() const
Definition: StatWeight.h:40
Long64_t GetIndex() const
Definition: StatWeight.h:53
virtual void SetAddress(void *add)
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
TObject * ConstructedAt(Int_t idx)
virtual void Sort(Int_t upto=kMaxInt)
virtual void SetHeader(const char *header="", Option_t *option="")
TLegendEntry * AddEntry(const char *name, const char *label="", Option_t *option="lpf")
virtual const char * GetName() const
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Double_t Uniform(Double_t x1, Double_t x2)
const char * Data() const
virtual Int_t Fill()
virtual TBranch * GetBranch(const char *name)
virtual Long64_t GetEntries() const
virtual Int_t Branch(const char *folder, Int_t bufsize=32000, Int_t splitlevel=99)
long long Long64_t
TPaveText * pt
Double_t x[n]
const Int_t n
const long double mg
Definition: KVUnits.h:74
const long double g
masses
Definition: KVUnits.h:72
constexpr Double_t E()
#define mark(osub)