KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
mdweight.cpp
1//Created by KVClassFactory on Thu May 7 11:02:24 2015
2//Author: John Frankland,,,
3
4#include "mdweight.h"
5#include "TMath.h"
6#include "KVNucleusEvent.h"
7
9
10namespace MicroStat {
11
12
18
20 {
21 // energy distribution of particle in gas
22 // arg[0] = energy/available energy
23 // par[0] = number of particles in gas
24 // par[1] = mass ratio = massTot/(massTot-mPart)
25 Double_t e = arg[0];
26 Double_t N = par[0];
27 Double_t massRat = par[1];
28
29 Double_t val = 0.;
30 if ((e > 0) && (e * massRat < 1) && (N > 2)) {
31 val = TMath::Sqrt(e) * TMath::Power(1. - massRat * e, (3.*N - 8.) / 2.);
32 }
33 return val;
34 }
35
36
37
41
43 {
44 // find/create energy distribution for given number of particles
45 // N and mass ratio R.
46
47 TString fn = Form("fKE_N%d_R%f", N, R);
48 TF1* f = (TF1*)fKEDist.FindObject(fn);
49 if (!f) {
50 f = new TF1("mdweight::SPKEDist", mdweight::edist, 0, 1, 2);
51 f->SetRange(0, 1.);
52 f->SetParameters((Double_t)N, R);
53 f->SetNpx(100);
54 f->SetName(fn);
55 fKEDist.Add(f);
56 }
57 return f;
58 }
59
60
61
68
69 double CosThetaDist(double* x, double* par)
70 {
71 // parameterise angular distribution as function of ratio between
72 // P(cos theta=+/-1) and P(cos theta=0)
73 // par[0] = a = P(cos theta=+/-1)
74 // par[1] = b = P(cos theta=0)
75 // x = cosTheta
76 if (abs(par[0] - par[1]) < 1.e-6) return par[1];
77
78 double a = par[0];
79 double b = par[1];
80 double cosTheta = x[0];
81 double alpha = TMath::Pi() * cosTheta;
82 return (a - b) / 2.*(1. - TMath::Cos(alpha)) + b;
83 }
84
85
86
88
90 : fCosTheta("CosTheta", CosThetaDist, -1, 1, 2)
91 {
95 eDisp = 0.0;
96 massTot = 0.0;
97 massTot0 = 0.0;
98 px = 0.0;
99 py = 0.0;
100 pz = 0.0;
101 SetAnisotropy(1, 1);
102 }
103
104
105
108
110 {
111 // Destructor
112 }
113
114
115
119
121 {
122 // Set available energy, E, and calculate statistical weight
123 // for this event
124
125 if (E <= 0) {
127 setWeight(0.);
128 return;
129 }
131 Double_t N = e->GetMult();
132 Double_t logmass_sum, mass_sum;
133 logmass_sum = mass_sum = 0.;
134 for (auto& n : EventIterator(e)) {
135 Double_t m = n.GetMass();
136 logmass_sum += TMath::Log(m);
137 mass_sum += m;
138 }
139 A = (3 * (N - 1) / 2.) * log2pi - TMath::LnGamma(3 * (N - 1) / 2.)
140 + 1.5 * (logmass_sum - TMath::Log(mass_sum));
141 B = (3 * N - 5) / 2.;
142
144 setWeight(w);
145 }
146
147
148
152
154 {
155 // Call before generating an event with StatWeight::GenerateEvent
156 // using the given partition and available energy
157
158 massTot0 = 0;
159
160 for (auto& e : EventIterator(partition)) massTot0 += e.GetMass();
161
163 }
164
165
166
170
172 {
173 // Called by StatWeight::GenerateEvent before generating another
174 // event using the same partition as the last
177 px = py = pz = 0.0;
178 }
179
180
181
187
189 {
190 // Called by StatWeight::GenerateEvent when adding a particle to the event
191 // N is the number of particles still to add including this one
192 //
193 // The algorithm was written by Daniel Cussol (LPC Caen, France).
194
195 Double_t mPart = part->GetMass();
196 Double_t rR = mPart / massTot;
197 Double_t ratio;
198 if (rR < 1.) ratio = 1. / (1. - rR);
199 else ratio = 1.;
200 Double_t ec = 0.; // kinetic energy of particle
201 Double_t ppx, ppy, ppz; // components of particle momentum
202 ppx = ppy = ppz = 0;
203 if (N >= 2) {
204 Double_t p = 0.; //momentum to give particle
205 if (N > 2) {
206 // draw random KE from 1-particle distribution for given N & ratio
207 ec = eDisp * getKEdist(N, ratio)->GetRandom();
208 p = sqrt(2.*mPart * ec);
209 } else {
210 // last 2 particles: share remaining available energy
211 p = sqrt(2.*(massTot - mPart) * mPart * eDisp / massTot);
212 ec = p * p / 2. / mPart;
213 }
214 Double_t ct = fCosTheta.GetRandom(-1, 1); //1. - 2.*gRandom->Rndm();
215 Double_t st = TMath::Sqrt(1. - ct * ct);
216 Double_t phi = gRandom->Rndm() * 2.*TMath::Pi();
217 ppz = ct * p;
218 ppx = st * TMath::Cos(phi) * p;
219 ppy = st * TMath::Sin(phi) * p;
220 ppx += px * rR;
221 ppy += py * rR;
222 ppz += pz * rR;
223 } else if (N == 1) {
224 ppx = px;
225 ppy = py;
226 ppz = pz;
227 ec = 0;
228 }
229 part->SetMomentum(ppx, ppy, ppz);
230
231 eDisp -= ec * ratio;
232 massTot -= mPart;
233 px -= ppx;
234 py -= ppy;
235 pz -= ppz;
236 }
237
238
239}/* namespace MicroStat */
240
int Int_t
#define f(i)
#define e(i)
double Double_t
#define N
winID w
winID h TVirtualViewer3D TVirtualGLPainter p
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 b
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
Class for iterating over nuclei in events accessed through base pointer/reference.
Abstract base class container for multi-particle events.
Definition KVEvent.h:67
Description of properties and kinematics of atomic nuclei.
Definition KVNucleus.h:126
void SetMomentum(const TVector3 *v)
Definition KVParticle.h:542
Double_t GetMass() const
Definition KVParticle.h:574
virtual TObject * FindObject(const char *name) const
virtual void SetOwner(Bool_t enable=kTRUE)
virtual void Add(TObject *obj)
void setAvailableEnergy(Double_t e)
Definition StatWeight.h:34
Double_t GetAvailableEnergy() const
Definition StatWeight.h:48
void setWeight(Double_t w)
Definition StatWeight.h:30
Calculate molecular dynamics ensemble weights for events .
Definition mdweight.h:20
Double_t log10twelve
Definition mdweight.h:22
void initGenerateEvent(KVEvent *partition)
Definition mdweight.cpp:153
virtual void SetWeight(KVEvent *e, Double_t E)
Definition mdweight.cpp:120
TF1 * getKEdist(Int_t, Double_t)
function used to draw random CosTheta values
Definition mdweight.cpp:42
static Double_t edist(Double_t *, Double_t *)
Definition mdweight.cpp:19
Double_t massTot0
Definition mdweight.h:23
virtual ~mdweight()
Destructor.
Definition mdweight.cpp:109
virtual void nextparticleGenerateEvent(Int_t, KVNucleus *)
Definition mdweight.cpp:188
void SetAnisotropy(double a, double b)
Definition mdweight.h:38
KVHashList fKEDist
Definition mdweight.h:24
virtual Double_t GetRandom(Double_t xmin, Double_t xmax, TRandom *rng=nullptr, Option_t *opt=nullptr)
Double_t Rndm() override
RVec< PromoteType< T > > abs(const RVec< T > &v)
Double_t x[n]
const Int_t n
double CosThetaDist(double *x, double *par)
Definition mdweight.cpp:69
Expr< UnaryOp< Sqrt< T >, Expr< A, T, D, D2, R >, T >, T, D, D2, R > sqrt(const Expr< A, T, D, D2, R > &rhs)
Double_t Exp(Double_t x)
Double_t Power(Double_t x, Double_t y)
constexpr Double_t E()
Double_t Log(Double_t x)
Double_t Sqrt(Double_t x)
Double_t Cos(Double_t)
constexpr Double_t Pi()
constexpr Double_t R()
Double_t LnGamma(Double_t z)
Double_t Sin(Double_t)
constexpr Double_t TwoPi()
TMarker m
TArc a
ClassImp(TPyArg)