#ifndef __ClassTraj_H
#define __ClassTraj_H
#include "KVRungeKutta.h"
#include "KVWilckeReactionParameters.h"
Double_t EtotDeb, EtotFin, TKEL, TKE, VPOT;
Bool_t fused, inelastic, incoming, interacting, outgoing;
public:
virtual ~ClassTraj();
{
}
{
return sqrt(
pow(rp[0] - rt[0], 2) +
pow(rp[1] - rt[1], 2));
}
{
return sqrt(
pow(vp[0] - vt[0], 2) +
pow(vp[1] - vt[1], 2));
}
{
return (rp[0] - rt[0]) / dr(rp, rt);
}
{
return (rp[1] - rt[1]) / dr(rp, rt);
}
{
return (vp[0] - vt[0]) / dv(vp, vt);
}
{
return (vp[1] - vt[1]) / dv(vp, vt);
}
{
(VP[0] + VCM[0]) /
sqrt(
pow(VP[0] + VCM[0], 2) +
pow(VP[1] + VCM[1], 2))), VP[1] + VCM[1]);
};
{
};
{
Double_t lz = dv_x(VP, VT) * dr_y(RP, RT) - dv_y(VP, VT) * dr_x(RP, RT);
return lz;
}
void Reset()
{
if (CMPositionPad) delete CMPositionPad;
if (LabPositionPad) delete LabPositionPad;
if (!ProjARC) {
}
if (!TargARC) {
}
if (!ProjARClab) {
ProjARClab =
new TArc(0, 0, 2 * fWilcke.
GetCP());
}
if (!TargARClab) {
TargARClab =
new TArc(0, 0, 2 * fWilcke.
GetCT());
}
CMPositionPad =
new TCanvas(
"CMPOS",
"CM Positions", 800, 500);
if (infox) delete infox;
UpdateInfos(0, 0);
LabPositionPad =
new TCanvas(
"LABPOS",
"Lab Positions", 800, 500);
LabPositionPad->
DrawFrame(-50, -300, 550, 300);
}
{
if (infox) {
if (outgoing)infox->
SetText(-40, -30,
Form(
"b=%.1ffm t=%.0ffm/c OUTGOING",
b, time));
else if (incoming)infox->
SetText(-40, -30,
Form(
"b=%.1ffm t=%.0ffm/c INCOMING",
b, time));
else if (interacting)infox->
SetText(-40, -30,
Form(
"b=%.1ffm t=%.0ffm/c INTERACTING",
b, time));
else infox->
SetText(-40, -30,
Form(
"b=%.1ffm t=%.0ffm/c",
b, time));
}
}
};
#endif
#define ClassDef(name, id)
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void on
char * Form(const char *fmt,...)
Description of properties and kinematics of atomic nuclei.
static Double_t hbar
hbar*c in MeV.fm
Adaptive step-size 4th order Runge-Kutta ODE integrator from Numerical Recipes.
virtual void CalcDerivs(Double_t X, Double_t *Y, Double_t *DYDX)=0
Reaction parameters for heavy-ion collisions from systematics of Wilcke et al.
virtual void SetFillColor(Color_t fcolor)
virtual void SetFillStyle(Style_t fstyle)
virtual void SetLineStyle(Style_t lstyle)
virtual void SetMarkerSize(Size_t msize=1)
void Draw(Option_t *option="") override
void Draw(Option_t *option="") override
virtual void Draw(Option_t *option="")
TH1F * DrawFrame(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, const char *title="") override
virtual void SetText(Double_t x, Double_t y, const char *text)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
RVec< PromoteType< T > > acos(const RVec< T > &v)
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 Sign(Double_t a, Double_t b)
constexpr Double_t RadToDeg()
#include "ClassTraj.h"
using namespace std;
#include "KVTreeAnalyzer.h"
{
alpha = 0;
RP = coords;
RT = coords + 2;
VP = coords + 4;
VT = coords + 6;
CMPositionPad = LabPositionPad = 0;
ProjARC = TargARC = ProjARClab = TargARClab = 0;
infox = 0;
if (fVisualization) Reset();
}
void ClassTraj::UpdateARCPositions(
Double_t t)
{
ProjARC->SetX1(RP[0]);
ProjARC->SetY1(RP[1]);
TargARC->SetX1(RT[0]);
TargARC->SetY1(RT[1]);
ProjTRAJ->SetPoint(nTRAJPoints++, RP[0], RP[1]);
ProjARClab->SetX1(RP[0] + RCM[0] + VCM[0]*t);
ProjARClab->SetY1(RP[1] + RCM[1] + VCM[1]*t);
TargARClab->SetX1(RT[0] + RCM[0] + VCM[0]*t);
TargARClab->SetY1(RT[1] + RCM[1] + VCM[1]*t);
ProjTRAJlab->SetPoint(nTRAJPointslab++, RP[0] + RCM[0] + VCM[0]*t, RP[1] + RCM[1] + VCM[1]*t);
}
ClassTraj::~ClassTraj()
{
delete [] coords;
}
{
if (fVisualization)UpdateInfos(b, 0);
alpha = fric;
if (fVisualization) {
}
RT[0] = RT[1] = 0;
for (int i = 0; i < 2; i++) {
RCM[i] = (fWilcke.GetAP() * RP[i] + fWilcke.GetAT() * RT[i]) / (1.*fWilcke.GetAC());
RP[i] -= RCM[i];
RT[i] -= RCM[i];
}
E0 = fWilcke.ECM(ELAB) - fWilcke.GetTotalPotential()->Eval(dr(RP, RT));
VP[0] = V0;
VP[1] = VT[0] = VT[1] = 0.;
for (int i = 0; i < 2; i++) {
VCM[i] = (fWilcke.GetAP() * VP[i] + fWilcke.GetAT() * VT[i]) / (1.*fWilcke.GetAC());
VP[i] -= VCM[i];
VT[i] -= VCM[i];
}
if (fVisualization) {
nTRAJPoints = 0;
nTRAJPointslab = 0;
UpdateARCPositions(0);
CMPositionPad->cd();
ProjARC->Draw();
TargARC->Draw();
ProjTRAJ->Draw("l");
CMPositionPad->Modified();
CMPositionPad->Update();
LabPositionPad->cd();
ProjARClab->Draw();
TargARClab->Draw();
ProjTRAJlab->Draw("l");
LabPositionPad->Modified();
LabPositionPad->Update();
}
}
{
Double_t Force = -fWilcke.GetTotalPotential()->Derivative(distance);
for (int i = 0; i < 2; i++) {
drp[i] = vp[i];
drt[i] = vt[i];
}
dvp[0] = fp * uv_x;
dvp[1] = fp * uv_y;
dvt[0] = ft * uv_x;
dvt[1] = ft * uv_y;
if (fInitialDeriv) {
if (distance <= minDist)minDist = distance;
if (incoming && distance < fWilcke.GetRint() + 10.) {
}
else if (interacting && distance > fWilcke.GetRint() + 10.) {
}
if (tcon < 0 && distance < fWilcke.GetRint()) tcon = t;
if (tcon > 0 && tsep < 0) {
if (t > tcon && distance > fWilcke.GetRint()) tsep = t;
}
}
if (alpha > 0) {
if (fInitialDeriv) {
}
if (distance < fWilcke.GetRint()) {
Double_t vpar = vrel * (vv_x * uv_x + vv_y * uv_y);
Double_t fforce = FRICTION * randomX + FRICTION;
dvp[0] += fp * (uv_x);
dvp[1] += fp * (uv_y);
dvt[0] += ft * (uv_x);
dvt[1] += ft * (uv_y);
}
}
}
{
tcon = tsep = -1.;
while (Time < tmax) {
if (outgoing) {
dt *= 1.5;
if (Time + dt > tmax) {
dt = tmax - Time;
}
}
Time += dt;
if (fVisualization) {
UpdateInfos(ImpactParameter, Time);
UpdateARCPositions(Time);
CMPositionPad->Modified();
CMPositionPad->Update();
LabPositionPad->Modified();
LabPositionPad->Update();
}
}
if (fVisualization) {
ProjTRAJ->SetLineWidth(1);
ProjTRAJlab->SetLineWidth(1);
ProjTRAJ->SetLineStyle(2);
ProjTRAJlab->SetLineStyle(2);
}
}
{
fWilcke.GetAP() * (VP[0] * VP[0] + VP[1] * VP[1])
+ fWilcke.GetAT() * (VT[0] * VT[0] + VT[1] * VT[1]));
Double_t epot = fWilcke.GetTotalPotential()->Eval(dr(RP, RT));
return ekin + epot;
}
{
alpha = friction;
InitTree(e);
DoOneTraj(tmin, friction, b, e, tmax, dt);
}
theFile->Write();
theFile->Close();
}
{
TString simTitle =
Form(
"%s %f MeV/u (alpha=%f)", sysName.Data(), e, alpha);
theFile =
new TFile(filename,
"recreate");
theTree =
new TTree(
"ClassTrajSim", simTitle);
theTree->Branch("ImpactParameter", &ImpactParameter);
theTree->Branch("AngMom", &AngMom);
theTree->Branch("TKEL", &TKEL);
theTree->Branch("TKE", &TKE);
theTree->Branch("VPOT", &VPOT);
theTree->Branch("minDist", &minDist);
theTree->Branch("ThetaProjLab", &ThetaProjLab);
theTree->Branch("absThetaProjLab", &absThetaProjLab);
theTree->Branch("EProjLab", &EProjLab);
theTree->Branch("fused", &fused);
theTree->Branch("inelastic", &inelastic);
theTree->Branch("RotationTime", &RotationTime);
theTree->Branch("RotationAngle", &RotationAngle);
theTree->Branch("tcon", &tcon);
theTree->Branch("tsep", &tsep);
}
{
minDist = 500;
InitTrajectory(b, e, friction);
AngMom = fWilcke.k(e) *
b;
Run(tmin, tmax, dt);
ThetaProjLab = ThetaLabProj();
EProjLab = ELabProj();
absThetaProjLab =
abs(ThetaProjLab);
TKE = TotalEnergy();
TKEL = fWilcke.ECM(ELAB) - TKE;
RotationTime = (tsep > 0 ? tsep - tcon : tmax - tcon);
VPOT = fWilcke.GetTotalPotential()->Eval(dr(RP, RT));
if (dr(RP, RT) < fWilcke.GetRint()) fused = kTRUE;
if (TKEL > 1.e-4 && !fused) inelastic =
kTRUE;
theTree->Fill();
}
{
alpha = friction;
InitTree(e);
while (ntraj--) {
if (!(ntraj % 5000)) cout << ntraj << " more to go..." << endl;
DoOneTraj(tmin, friction, b, e, tmax, dt);
}
theFile->Write();
theTree->StartViewer();
}
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
R__EXTERN TRandom * gRandom
virtual void Print(Option_t *option="") const
const Char_t * GetSymbol(Option_t *opt="") const
static Double_t kAMU
atomic mass unit in MeV
virtual void SetLineWidth(Width_t lwidth)
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
virtual Double_t Uniform(Double_t x1, Double_t x2)
RooCmdArg Integrate(bool flag)
RVec< PromoteType< T > > abs(const RVec< T > &v)
#define R0(v, w, x, y, z, i)