KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KV2Body.cpp
1/***************************************************************************
2 KV2Body.cpp - description
3 -------------------
4 begin : 28/11/2003
5 copyright : (C) 2003 by J.D. Frankland
6 email : frankland@ganil.fr
7
8$Id: KV2Body.cpp,v 1.4 2009/02/02 13:52:29 ebonnet Exp $
9 ***************************************************************************/
10
11/***************************************************************************
12 * *
13 * This program is free software; you can redistribute it and/or modify *
14 * it under the terms of the GNU General Public License as published by *
15 * the Free Software Foundation; either version 2 of the License, or *
16 * (at your option) any later version. *
17 * *
18 ***************************************************************************/
19
20#include "KV2Body.h"
21#include "TROOT.h"
22
23using namespace std;
24
26
27
28
29
31
32void KV2Body::init()
33{
34 //Default initialisations
35 WLT = WCT = 0.;
36 for (int i = 0; i < 5; i++) {
37 WC[i] = 0.;
38 VC[i] = 0.;
39 EC[i] = 0.;
40 K[i] = 0.;
41 TETAMAX[i] = 0.;
42 TETAMIN[i] = 0.;
43 fThetaLabVsThetaCM[i] = nullptr;
44 fELabVsThetaCM[i] = nullptr;
45 fELabVsThetaLab[i] = nullptr;
46 }
47 fKoxReactionXSec = nullptr;
48 fEqbmChargeState = nullptr;
49 fEqbmChargeStateShSol = nullptr;
50 fEqbmChargeStateShGas = nullptr;
51 fSetOutgoing = kFALSE;
52
53 SetIntegralPrecision(1e-10);
54}
55
56
57
60
61KV2Body::KV2Body(): fNuclei(5), fEDiss(0)
62{
63 //default ctor
64 init();
65}
66
67
68
86
87KV2Body::KV2Body(const Char_t* systemname) : fNuclei(5), fEDiss(0)
88{
89 //Set up calculation defining entrance channel following this prescription:
90 //
91 // [Projectile_Symbol]+[Target_Symbol]@[Incident_Energy]MeV/A
92 // [Projectile_Symbol]+[Target_Symbol]@[Incident_Energy]MeV/u
93 // [Projectile_Symbol]+[Target_Symbol]
94 //
95 //Any spaces will be ignored.
96 //
97 //Valid examples:
98 //
99 // 129Xe+119Sn@50.0MeV/A
100 // 58Ni + 64Ni @ 32 MeV/u
101 // U+U@5MeV/A
102 // Ta+Zn
103 //
104 //If the format is not respected, this object will be a zombie (IsZombie() returns kTRUE)
105
106 init();
107 KVString sener = "";
108 Double_t ee = 0;
109 KVString syst(systemname);
110 syst.ReplaceAll(" ", "");
111 KVString scouple(syst);
112
113 if (syst.Contains("@")) {
114 syst.Begin("@");
115 scouple = syst.Next();
116 sener = syst.Next();
117 }
118 if (sener != "") {
119 if (sener.Contains("MeV/A")) {
120 sener.ReplaceAll("MeV/A", "");
121 ee = sener.Atof();
122 }
123 else if (sener.Contains("MeV/u")) {
124 sener.ReplaceAll("MeV/u", "");
125 ee = sener.Atof();
126 }
127 else {
128 MakeZombie();
129 }
130 }
131 scouple.Begin("+");
132 KVNucleus nnuc1(scouple.Next(), ee);
133 if (nnuc1.IsZombie()) MakeZombie();
134 else {
135 fNuclei[1] = nnuc1;
136 KVNucleus nnuc2(scouple.Next());
137 if (nnuc2.IsZombie()) MakeZombie();
138 else {
139 fNuclei[2] = nnuc2;
140 }
141 }
142 if (!IsZombie()) {
143 SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), ee));
144 }
145}
146
147
148
151
152KV2Body::KV2Body(KVNucleus*, KVNucleus* cib, KVNucleus* proj_out, Double_t): fNuclei(5), fEDiss(0)
153{
154 // Deprecated. Do not use.
155
156 Obsolete("KV2Body(KVNucleus*, KVNucleus*, KVNucleus*, Double_t)", "1.11", "1.12");
157
158 if (!cib) {
159 Info("KV2Body", "Use constructor KV2Body(const KVNucleus& compound) to simulate binary decay of compound");
160 }
161 else if (!proj_out) {
162 Info("KV2Body", "Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ) to define entrance channel");
163 }
164 else {
165 Info("KV2Body", "Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ, const KVNucleus& outgoing) to define entrance & exit channels");
166 }
167 init();
168}
169
170
171
195
196KV2Body::KV2Body(const KVNucleus& compound, double Exx)
197 : fNuclei(5), fEDiss(-Exx)
198{
199 // Set up kinematics of binary decay of compound nucleus.
200 //
201 // The excitation energy of the CN (in MeV) is either given:
202 //
203 // - in the KVNucleus object itself (i.e. compound.GetExcitEnergy() gives the E* of CN)
204 // - or in the variable Exx: in this case any E* in the KVNucleus will be ignored
205 // - or by calling SetExcitEnergy or SetDissipatedEnergy after this constructor with the NEGATIVE E*
206 //
207 //Usage:
208 //~~~~~~~~~~{.cpp}
209 // KVNucleus CN("204Pb");
210 // CN.SetExcitEnergy(1.5*CN.GetA());
211 // KV2Body CNdecay(CN);
212 //
213 // /* or: KV2Body CNdecay("204Pb", 1.5*204);
214 // * or: KV2Body CNdecay("204Pb");
215 // * CNdecay.SetExcitEnergy(-1.5*204); */
216 //
217 // // calculate decay
218 // KVNucleus alpha("4He", 67.351/4.);
219 // CNdecay.SetOutgoing(alpha);
220 //~~~~~~~~~~
221
222 init();
223 fNuclei[1] = compound;
224 if (fEDiss == 0) fEDiss = compound.GetExcitEnergy();
225}
226
227
228
241
242KV2Body::KV2Body(const KVNucleus& proj, const KVNucleus& targ, double Ediss):
243 fNuclei(5), fEDiss(Ediss)
244{
245 //Set up calculation of basic binary reaction for given projectile and target.
246 //
247 //By default the dissipated energy Ediss is zero (elastic reaction).
248 //
249 //Usage:
250 //~~~~~~~~~~{.cpp}
251 // KV2Body reaction(KVNucleus("129Xe",50), "natSn");
252 //
253 // // or with C++11 (ROOT6):
254 // KV2Body reaction({"129Xe",50}, "natSn");
255 //~~~~~~~~~~
256
257 init();
258 fNuclei[1] = proj;
259 fNuclei[2] = targ;
260 SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
261}
262
263
264
281
282KV2Body::KV2Body(const KVNucleus& proj, const KVNucleus& targ, const KVNucleus& outgoing, double Ediss):
283 fNuclei(5), fEDiss(Ediss)
284{
285 //Set up calculation of basic binary reaction for given projectile and target with definition
286 //of exit channel (outgoing projectile- or target-like fragment).
287 //
288 //By default the dissipated energy is zero (elastic reaction).
289 //
290 //Any excitation energy of the outgoing fragment will be taken into account, e.g.
291 //for quasi-elastic scattering leaving the target in an excited state.
292 //
293 //Usage:
294 //~~~~~~~~~~{.cpp}
295 // KV2Body reaction(KVNucleus("129Xe",50), "natSn", KVNucleus("24Mg",35));
296 //
297 // // or with C++11 (ROOT6):
298 // KV2Body reaction({"129Xe",50}, "natSn", {"24Mg",35});
299 //~~~~~~~~~~
300
301 init();
302 fNuclei[1] = proj;
303 fNuclei[2] = targ;
304 if (outgoing.GetExcitEnergy() > 0) fEDiss += outgoing.GetExcitEnergy();
305 SetOutgoing(outgoing);
306 SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
307}
308
309
310
313
315{
316 //Set target for reaction.
317 fNuclei[2] = targ;
319}
320
321
322
325
327{
328 //Set target for reaction
329
330 fNuclei[2].SetZandA(z, a);
332}
333
334
335
338
340{
341 //Set projectile for reaction.
342 fNuclei[1] = proj;
344}
345
346
347
350
352{
353 //Set projectile for reaction
354 fNuclei[1].SetZandA(z, a);
356}
357
358
359
361
363{
366 for (int i = 0; i < 5; i++) {
367 if (fThetaLabVsThetaCM[i]) delete fThetaLabVsThetaCM[i];
368 if (fELabVsThetaCM[i]) delete fELabVsThetaCM[i];
369 if (fELabVsThetaLab[i]) delete fELabVsThetaLab[i];
370 }
371}
372
373
374
375
379
381{
382 //Static function, calculates relativistic velocity (in cm/ns) from rest mass and total
383 //energy (i.e. KE + mass) 'E' for any particle
384
385 Double_t p = TMath::Sqrt(E * E - mass * mass);
386 return KVParticle::C() * p / E;
387}
388
389
390
391
396
397void KV2Body::SetOutgoing(const KVNucleus& proj_out)
398{
399 // Set outgoing projectile-like nucleus properties.
400 //
401 // The properties of the outgoing target-like nucleus will be deduced from mass, charge and momentum/energy conservation.
402
403 SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
405 fNuclei[3] = proj_out;
407 // all nuclei should now have E* set to zero in order to use ground state masses in kinematics
408 for (unsigned i = 1; i <= 4; ++i) if (fNuclei[i].IsDefined()) fNuclei[i].SetExcitEnergy(0);
409}
410
411
412
413
421
423{
424 // Private method, used to deduce 4th nucleus (target-like) from projectile, target
425 // and outgoing projectile using conservation of mass, momentum and energy.
426 //
427 // if the outgoing nucleus set by the user is equal to the compound nucleus
428 // formed by projectile and target, but the excitation energy of the CN
429 // was not set by the user, we calculate it here.
430
432 if (GetNucleus(2))
433 sum = *GetNucleus(1) + *GetNucleus(2);
434 else sum = *GetNucleus(1);
435 KVNucleus tmp4 = sum - *GetNucleus(3);
436 if (!tmp4.IsDefined()) {
437 // nucleus 3 is the CN proj+targ
438 // has E* been defined ?
439 if (fEDiss == 0) fEDiss = sum.GetExcitEnergy();
440 }
441 fNuclei[4] = tmp4;
442}
443
444
445
446
455
457{
458 //Return pointer to nucleus i (1 <= i <= 4)
459 //
460 // Entrance channel nuclei ..... i=1 : projectile i=2 : target
461 //
462 // Exit channel nuclei ..... i=3 : projectile-like i=4 : target-like
463 //
464 // Will return `nullptr` if any nucleus is undefined
465
466 if (i > 0 && i < 5) {
467 if (fNuclei[i].IsDefined()) return (KVNucleus*) &fNuclei[i];
468 return nullptr;
469 }
470 Warning("GetNucleus(Int_t i)", "Index i out of bounds, i=%d", i);
471 return nullptr;
472}
473
474
475
476
479
481{
482 //Calculate Q-value for reaction, including dissipated (excitation) energy
483
484 if (!fSetOutgoing) {
485 Warning("GetQReaction", "Parameters for outgoing nuclei not set");
486 return 0.0;
487 }
489
490 return QR;
491}
492
493
494
495
498
500{
501 //Calculate Q-value for reaction, assuming all nuclei in ground state
502
503 if (!fSetOutgoing) {
504 Warning("GetQGroundStates",
505 "Parameters for outgoing nuclei not set");
506 return 0.0;
507 }
508 Double_t QGG =
510 (GetNucleus(2) ? GetNucleus(2)->GetMassExcess() : 0.)
511 - GetNucleus(3)->GetMassExcess() -
512 (GetNucleus(4) ? GetNucleus(4)->GetMassExcess() : 0.);
513 return QGG;
514}
515
516
517
518
521
523{
524 //Return available kinetic energy in centre of mass
525
526 return WCT - (GetNucleus(1)->GetMass() +
527 (GetNucleus(2) ? GetNucleus(2)->GetMass() : 0.));
528}
529
530
531
535
537{
538 //Returns maximum scattering angle in lab for nuclei i=3 (quasiproj)
539 //and i=4 (quasitarget)
540 if (i < 3 || i > 4) {
541 Warning("GetMaxAngleLab(Int_t i)",
542 "Returns maximum scattering angle in lab for nuclei i=3 (quasiproj) and i=4 (quasitarget)");
543 return 0.;
544 }
545 return TETAMAX[i];
546}
547
548
549
553
555{
556 //Returns minimum scattering angle in lab for nuclei i=3 (quasiproj)
557 //and i=4 (quasitarget)
558 if (i < 3 || i > 4) {
559 Warning("GetMinAngleLab(Int_t i)",
560 "Returns minimum scattering angle in lab for nuclei i=3 (quasiproj) and i=4 (quasitarget)");
561 return 0.;
562 }
563 return TETAMIN[i];
564}
565
566
567
568
571
573{
574 //Return vector velocity of centre of mass of reaction (units: cm/ns)
575
576 return VCM;
577}
578
579
580
581
586
588{
589 //Return velocity of nucleus "i" in the centre of mass of the reaction
590 // Entrance channel nuclei ..... i=1 : projectile i=2 : target
591 // Exit channel nuclei ..... i=3 : projectile-like i=4 : target-like
592 return VC[i];
593}
594
595
596
597
602
604{
605 //Calculate laboratory grazing angle.
606 // i = 1 (default) : projectile
607 // i = 2 : target
608
609 if (i < 1 || i > 2) {
610 Warning("GetLabGrazingAngle(Int_t i)",
611 "i should be 1 (proj) or 2 (targ)");
612 return 0.0;
613 }
614 if (!GetNucleus(2)) {
615 Warning("GetLabGrazingAngle(Int_t i)",
616 "No target defined for reaction");
617 return 0.0;
618 }
619 Double_t RP = 1.28 * TMath::Power(GetNucleus(1)->GetA(), 1. / 3.)
620 - 0.76 + 0.8 * TMath::Power(GetNucleus(1)->GetA(), -1. / 3.);
621 Double_t RT = 1.28 * TMath::Power(GetNucleus(2)->GetA(), 1. / 3.)
622 - 0.76 + 0.8 * TMath::Power(GetNucleus(2)->GetA(), -1. / 3.);
623 Double_t CP = RP * (1. - TMath::Power(1. / RP, 2.));
624 Double_t CT = RT * (1. - TMath::Power(1. / RT, 2.));
625 Double_t RINT = CP + CT + 4.49 - (CT + CP) / 6.35;
626 Double_t BAR =
627 1.44 * GetNucleus(1)->GetZ() * GetNucleus(2)->GetZ() / RINT;
628 if (GetCMEnergy() < BAR) {
629 //below Coulomb barrier
630 switch (i) {
631 case 1:
632 return 180.;
633 break;
634 case 2:
635 return 90.;
636 break;
637 }
638 }
639 Double_t X = 2. * RINT * GetCMEnergy();
640 Double_t Y =
641 X / (GetNucleus(1)->GetZ() * GetNucleus(2)->GetZ() * 1.44) - 1.;
642 Double_t U = 1. / Y;
643 Double_t GR = 2. * TMath::ASin(U);
644 Double_t GPART = TMath::Pi() - GR;
645 Double_t ARAZ =
646 TMath::Sin(GR) / (TMath::Cos(GR) +
647 GetNucleus(1)->GetA() / (1.0 *
648 GetNucleus(2)->GetA()));
649 Double_t GRAZ = TMath::ATan(ARAZ);
650 if (GRAZ <= 0.)
651 GRAZ += TMath::Pi();
652
653 switch (i) {
654 case 1: //projectile
655 return (GRAZ * TMath::RadToDeg());
656 break;
657 case 2: //target
658 return (GPART * TMath::RadToDeg() / 2.);
659 break;
660 }
661 return -99.;
662}
663
664
665
666
676
678{
679 // Called to make kinematics calculation for all nuclei, use once properties of
680 // entrance-channel (and, if necessary, exit-channel) nuclei have been defined.
681 //
682 // If a compound decay is to be calculated (only 1 nucleus in entrance channel),
683 // outgoing (decay) product must be defined before calling this method.
684 //
685 // For a 2-body entrance channel, if no exit-channel nuclei are defined,
686 // we assume elastic scattering (i.e. identical outgoing nuclei)
687
688 KVNucleus* Nuc1 = GetNucleus(1);
689 KVNucleus* Nuc2 = GetNucleus(2);
690
691 // compound decay ?
692 if (!Nuc2 && !fSetOutgoing) {
693 Error("CalculateKinematics", "Set outgoing (decay) product first");
694 return;
695 }
696
697 // call SetOutgoing if not already done
698 if (!fSetOutgoing) SetOutgoing(*Nuc1);
699
700 // set everything to zero
701 WLT = WCT = BCM = 0.;
702 for (int i = 0; i < 5; i++) {
703 WC[i] = 0.;
704 VC[i] = 0.;
705 EC[i] = 0.;
706 K[i] = 0.;
707 TETAMAX[i] = 0.;
708 TETAMIN[i] = 0.;
709 }
710 VCM.SetXYZ(0., 0., 0.);
711
712 //total energy (T + m) of entrance channel
713 WLT = Nuc1->E();
714 //velocity of CM
715 VCM = Nuc1->GetMomentum();
716 if (Nuc2) {
717 WLT += Nuc2->E();
718 VCM += Nuc2->GetMomentum();
719 }
720
721 //beta of CM
722 BCM = VCM.Mag() / WLT;
723
724 VCM *= (1. / WLT) * KVParticle::C();
725 Nuc1->SetFrame("CM", KVFrameTransform(VCM, kFALSE));
726 if (Nuc2)
727 Nuc2->SetFrame("CM", KVFrameTransform(VCM, kFALSE));
728
729 //total energy in CM
730 WCT = (GetCMGamma() > 0.0 ? WLT / GetCMGamma() : 0.);
731 if (WCT == 0.0)
732 return;
733
734 //total energy proj in CM
735 WC[1] = Nuc1->GetFrame("CM")->E();
736 //kinetic energy proj in CM
737 EC[1] = Nuc1->GetFrame("CM")->GetKE();
738 VC[1] = Nuc1->GetFrame("CM")->GetVelocity().Mag();
739 K[1] = VCM.Mag() / VC[1];
740 if (Nuc2) {
741 WC[2] = Nuc2->GetFrame("CM")->E();
742 EC[2] = Nuc2->GetFrame("CM")->GetKE();
743 VC[2] = Nuc2->GetFrame("CM")->GetVelocity().Mag();
744 K[2] = VCM.Mag() / VC[2];
745 }
746
747 Double_t AM3 = GetNucleus(3)->GetMass();
748 Double_t AM4 = (GetNucleus(4) ? GetNucleus(4)->GetMass() : 0.0);
749 //total cm energy of nucleus 3 (quasiproj)
750 WC[3] = (WCT - GetEDiss()) / 2. + (AM3 - AM4) * (AM3 + AM4)
751 / (2. * (WCT - GetEDiss()));
752 VC[3] = GetVelocity(AM3, WC[3]);
753 //cm kinetic energy
754 EC[3] = WC[3] - AM3;
755 if (VC[3] > 0.) K[3] = VCM.Mag() / VC[3];
756
757 Double_t T3MAX = 0.;
758 if (TMath::AreEqualAbs(K[3], 1., 1.e-16))
759 T3MAX = TMath::PiOver2();
760 else if (K[3] < 1.)
761 T3MAX = TMath::Pi();
762 else
763 T3MAX = TMath::ATan((1. / TMath::Sqrt(K[3] * K[3] - 1.)) / GetCMGamma());
764 TETAMAX[3] = T3MAX * TMath::RadToDeg();
765
766 if (!GetNucleus(4))
767 return; //only valid for binary channels
768 //total cm energy of nucleus 4 (quasitarg)
769 WC[4] = WCT - GetEDiss() - WC[3];
770 VC[4] = GetVelocity(AM4, WC[4]);
771 //cm kinetic energy
772 EC[4] = WC[4] - AM4;
773 if (VC[4] > 0.) K[4] = VCM.Mag() / VC[4];
774
775 Double_t T4MAX = 0.;
776 if (TMath::AreEqualAbs(K[4], 1., 1.e-16))
777 T4MAX = TMath::PiOver2();
778 else if (K[4] < 1.)
779 T4MAX = TMath::Pi();
780 else
781 T4MAX = TMath::ATan((1. / TMath::Sqrt(K[4] * K[4] - 1.)) / GetCMGamma());
782 if (TMath::Abs(GetQReaction()) < 1.E-08 && K[4] < 1.)
783 T4MAX = TMath::PiOver2();
784 TETAMAX[4] = T4MAX * TMath::RadToDeg();
785}
786
787
788
789
794
796{
797 //Returns kinetic energy of nucleus "i" in the centre of mass of the reaction
798 // Entrance channel nuclei ..... i=1 : projectile i=2 : target
799 // Exit channel nuclei ..... i=3 : projectile-like i=4 : target-like
800 return EC[i];
801}
802
803
804
805
812
813void KV2Body::Print(Option_t* opt) const
814{
815 // Print out characteristics of reaction.
816 //
817 // If a two-body exit channel has been defined, you can use the following options:
818 // opt = "ruth" : list Rutherford scattering cross-sections as a function of angle in laboratory frame
819 // opt = "lab" : list energies and angles in laboratory frame
820
821 cout << " ***** REACTION " << GetNucleus(1)->GetSymbol();
822 if (GetNucleus(2))
823 cout << " + " << GetNucleus(2)->GetSymbol();
824 if (GetNucleus(3)) {
825 cout << " ---> " << GetNucleus(3)->GetSymbol();
826 }
827 if (GetNucleus(4))
828 cout << " + " << GetNucleus(4)->GetSymbol();
829 cout << " ******" << endl;
830
831 cout << " E.LAB = " << GetNucleus(1)->GetEnergy() << " MEV";
832 if (GetNucleus(3)) {
833 cout << " QGG = " << GetQGroundStates() << " MEV";
834 }
835 cout << endl;
836 cout << " E.EXC = " << GetExcitEnergy() << " MEV";
837 if (GetNucleus(3)) {
838 cout << " ==> Q-REACTION = " << GetQReaction() << " MEV";
839 }
840 cout << endl;
841 cout << " AVAILABLE ENERGY IN C.M. : ECM = " << GetCMEnergy() <<
842 " MEV (" << GetCMEnergy() / (GetNucleus(1)->GetA() +
843 (GetNucleus(2) ? GetNucleus(2)->
844 GetA() : 0.)) << " MEV/A)" << endl;
845 cout << " PROJECTILE VELOCITY IN LAB " << GetNucleus(1)->GetV().Mag()
846 << " CM/NS ( " << GetNucleus(1)->Beta() << " * C )" << endl;
847 cout << " VELOCITY OF C.M. " << GetCMVelocity().
848 Mag() << " CM/NS" << endl;
849
850 for (int i = 1; i <= 4; i++) {
851 if (GetNucleus(i))
852 cout << " ENERGY - VELOCITY OF NUCLEUS " << i << " IN CM : " <<
853 GetCMEnergy(i) << " MEV " << GetCMVelocity(i) << " CM/NS (K=" << K[i] << ")" <<
854 endl;
855 }
856 if (GetNucleus(3)) {
857 cout << " MAXIMUM SCATTERING ANGLE IN LABORATORY" << endl;
858 cout << " THETA #3# " << GetMaxAngleLab(3) << " DEG." <<
859 endl;
860 if (GetNucleus(4))
861 cout << " THETA #4# " << GetMaxAngleLab(4) << " DEG." <<
862 endl;
863 }
864 if (GetNucleus(2)) {
865 cout << " GRAZING ANGLE IN LABORATORY : PROJECTILE " <<
866 GetLabGrazingAngle(1) << " DEG." << endl;
867 cout << " GRAZING ANGLE IN LABORATORY : TARGET " <<
868 GetLabGrazingAngle(2) << " DEG." << endl;
869 }
870
871 if (GetNucleus(3)) {
872 Int_t nsteps = 15;
873
874 if (!strcmp(opt, "lab") || !strcmp(opt, "LAB")
875 || !strcmp(opt, "Lab")) {
876 //laboratory angular distribution
877 cout << endl <<
878 " LABORATORY ANGULAR DISTRIBUTION" << endl
879 << endl;
880 cout <<
881 " TETA 3 EN.DIF.3 V3 TETA 4 EN.DIF.4 TETA 3"
882 << endl;
883 cout <<
884 " (LAB) (LAB) (LAB) (LAB) (LAB) (C.M)"
885 << endl;
886 Double_t dtheta = GetMaxAngleLab(3) / nsteps;
887 for (int step = 0; step <= nsteps; step++) {
888 Double_t theta = dtheta * step;
889
890 Double_t elabP1, elabP2;
891 Double_t tcmP1, tcmP2;
892 Double_t vP1, vP2;
893 Double_t tlT1, tlT2;
894 Double_t elabT1, elabT2;
895 GetELab(3, theta, 3, elabP1, elabP2);
896 GetThetaCM(theta, 3, tcmP1, tcmP2);
897 GetVLab(3, theta, 3, vP1, vP2);
898 GetThetaLab(4, theta, 3, tlT1, tlT2);
899 GetELab(4, theta, 3, elabT1, elabT2);
900 printf
901 (" %6.2f %7.2f/%7.2f %5.2f/%5.2f %6.2f/%6.2f %7.2f/%7.2f %6.2f/%6.2f\n",
902 theta, elabP1, elabP2, vP1, vP2,
903 tlT1, tlT2, elabT1, elabT2,
904 tcmP1, tcmP2);
905 }
906 }
907 if (GetNucleus(2) && (!strcmp(opt, "ruth") || !strcmp(opt, "RUTH")
908 || !strcmp(opt, "Ruth"))) {
909 //laboratory angular distribution with Rutherford x-section
910 cout << endl <<
911 " RUTHERFORD LABORATORY ANGULAR DISTRIBUTION" <<
912 endl << endl;
913 cout <<
914 " TETA 3 TETA 3 EN.DIF.3 SIG.RUT. SIG.RUT."
915 << endl;
916 cout <<
917 " (LAB) (CM) (LAB) (LAB) (b/sr) (CM) (b/sr)"
918 << endl;
919 //go up to grazing angle
920 Double_t dtheta = GetLabGrazingAngle(1) / nsteps;
921 for (int step = 0; step < nsteps; step++) {
922 Double_t theta = dtheta * (step + 1);
923 Double_t elabP1, elabP2;
924 Double_t tcm1, tcm2;
925 GetELab(3, theta, 3, elabP1, elabP2);
926 GetThetaCM(theta, 3, tcm1, tcm2);
927 printf
928 (" %6.2f %6.2f %7.2f %10.4g %10.4g\n",
929 theta, tcm1, elabP1,
930 GetXSecRuthLab(theta), GetXSecRuthCM(theta));
931 }
932 }
933 }
934}
935
936
937
944
945Int_t KV2Body::GetELab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& e1, Double_t& e2) const
946{
947 // Calculate laboratory kinetic energy of nucleus OfNucleus as a function of the lab angle of
948 // nucleus AngleNucleus.
949 //
950 // In general there may be two solutions for a given angle, therefore we return the number
951 // of solutions (0 if ThetaLab > max lab angle for nucleus in question).
952
953 e1 = e2 = -1.;
954 if (ThetaLab > TETAMAX[AngleNucleus]) return 0;
955 // calculate CM angle(s)
956 Double_t TCM1, TCM2;
957 Int_t nsol = GetThetaCM(ThetaLab, AngleNucleus, TCM1, TCM2);
958 TCM1 = GetThetaCM(OfNucleus, TCM1, AngleNucleus);
959 if (nsol > 1) TCM2 = GetThetaCM(OfNucleus, TCM2, AngleNucleus);
960 // calculate lab energie(s) of corresponding to CM angle(s)
961 e1 = GetELab(TCM1, OfNucleus);
962 if (nsol == 2) e2 = GetELab(TCM2, OfNucleus);
963 return nsol;
964}
965
966
967
971
972Int_t KV2Body::GetThetaCM(Double_t ThetaLab, Int_t OfNucleus, Double_t& t1, Double_t& t2) const
973{
974 // Calculate CM angle of nucleus OfNucleus as a function of its lab angle.
975 // Returns number of solutions (there may be <=2 solutions).
976
977 t1 = t2 = -1.;
978 KV2Body* kin = const_cast<KV2Body*>(this);
979 TF1* TLF = kin->GetThetaLabVsThetaCMFunc(OfNucleus);
980 Int_t nsol = FindRoots(TLF, 0., 180., ThetaLab, t1, t2);
981 return nsol;
982}
983
984
985
988
990{
991 // Function calculating lab energy of nucleus par[0] for any lab angle x[0]
992
993 Double_t e1, e2;
994 Int_t nsol = GetELab((Int_t)par[0], x[0], (Int_t)par[0], e1, e2);
995 if (!nsol) return 0;
996 //if(nsol>1) Warning("ELabVsThetaLab", "Two energies are possible for %f deg. : %f and %f",
997 // x[0],e1,e2);
998 return e1;
999}
1000
1001
1002
1005
1007{
1008 // Function calculating lab energy of nucleus par[0] for any CM angle x[0]
1009
1010 Int_t OfNucleus = (Int_t)par[0];
1011 Double_t TCM = x[0] * TMath::DegToRad();
1012 Double_t WL =
1013 WC[OfNucleus] * GetCMGamma() * (1. +
1014 BCM * (VC[OfNucleus] / KVParticle::C()) *
1015 TMath::Cos(TCM));
1016 return (WL - GetNucleus(OfNucleus)->GetMass());
1017}
1018
1019
1020
1024
1026{
1027 // Calculate Lab angle of nucleus as function of CM angle x[0]
1028 // par[0] = index of nucleus = 1, 2, 3, 4
1029
1030 Double_t ThetaCM = x[0] * TMath::DegToRad();
1031 Int_t OfNucleus = (Int_t)par[0];
1032 Double_t ThetaL = TMath::ATan2(TMath::Sin(ThetaCM) + 1e-10, (K[OfNucleus] + TMath::Cos(ThetaCM)) * GetCMGamma()) * TMath::RadToDeg();
1033
1034 if (ThetaL < 0.) ThetaL += 180.;
1035 return ThetaL;
1036}
1037
1038
1039
1043
1045{
1046 // Return TF1 giving lab energy of nucleus as function of CM angle
1047 // OfNucleus = 1 or 2 (entrance channel) or 3 or 4 (exit channel)
1048
1049 if (!fELabVsThetaCM[OfNucleus]) {
1050 fELabVsThetaCM[OfNucleus] = new TF1(Form("KV2Body:ELabVsThetaCM:%d", OfNucleus),
1051 const_cast<KV2Body*>(this), &KV2Body::ELabVsThetaCM, 0, 180, 1, "KV2Body", "ELabVsThetaCM");
1052 fELabVsThetaCM[OfNucleus]->SetNpx(1000);
1053 fELabVsThetaCM[OfNucleus]->SetParameter(0, OfNucleus);
1054 }
1055 return fELabVsThetaCM[OfNucleus];
1056}
1057
1058
1059
1063
1065{
1066 // Return TF1 giving lab energy of nucleus as function of its lab angle
1067 // OfNucleus = 1 or 2 (entrance channel) or 3 or 4 (exit channel)
1068
1069 if (!fELabVsThetaLab[OfNucleus]) {
1070 fELabVsThetaLab[OfNucleus] = new TF1(Form("KV2Body:ELabVsThetaLab:%d", OfNucleus),
1071 const_cast<KV2Body*>(this), &KV2Body::ELabVsThetaLab, 0, 180, 1, "KV2Body", "ELabVsThetaLab");
1072 fELabVsThetaLab[OfNucleus]->SetNpx(1000);
1073 fELabVsThetaLab[OfNucleus]->SetParameter(0, OfNucleus);
1074 }
1075 return fELabVsThetaLab[OfNucleus];
1076}
1077
1078
1079
1080
1087
1088Int_t KV2Body::GetVLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& v1, Double_t& v2) const
1089{
1090 // Calculate laboratory velocity of nucleus OfNucleus as a function of the lab angle of
1091 // nucleus AngleNucleus.
1092 //
1093 // In general there may be two solutions for a given angle, therefore we return the number
1094 // of solutions (0 if ThetaLab > max lab angle for nucleus in question).
1095
1096 Double_t e1, e2;
1097 v1 = v2 = 0.;
1098 Int_t nsol = GetELab(OfNucleus, ThetaLab, AngleNucleus, e1, e2);
1099 if (!nsol) return nsol;
1100 Double_t etot1 = e1 + GetNucleus(OfNucleus)->GetMass();
1101 Double_t etot2 = e2 + GetNucleus(OfNucleus)->GetMass();
1102 v1 = GetVelocity(GetNucleus(OfNucleus)->GetMass(), etot1);
1103 if (nsol > 1) v2 = GetVelocity(GetNucleus(OfNucleus)->GetMass(), etot2);
1104 return nsol;
1105}
1106
1107
1108
1109
1116
1117Int_t KV2Body::GetThetaLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t& t1, Double_t& t2) const
1118{
1119 // Calculate laboratory angle of nucleus OfNucleus as a function of the laboratory angle of
1120 // nucleus AngleNucleus.
1121 //
1122 // In general there may be two solutions for a given angle, therefore we return the number
1123 // of solutions (0 if ThetaLab > max lab angle for nucleus in question).
1124
1125 t1 = t2 = -1.;
1126 if (ThetaLab > TETAMAX[AngleNucleus]) return 0;
1127 if (!(TMath::Abs(OfNucleus - AngleNucleus) % 2)) {
1128 // same nucleus!
1129 t1 = ThetaLab;
1130 return 1;
1131 }
1132 // calculate CM angle(s)
1133 Double_t TCM1, TCM2;
1134 Int_t nsol = GetThetaCM(ThetaLab, AngleNucleus, TCM1, TCM2);
1135 TCM1 = GetThetaCM(OfNucleus, TCM1, AngleNucleus);
1136 if (nsol > 1) TCM2 = GetThetaCM(OfNucleus, TCM2, AngleNucleus);
1137 // calculate lab angle(s) of corresponding CM angle(s)
1138 t1 = GetThetaLab(TCM1, OfNucleus);
1139 if (nsol == 2) t2 = GetThetaLab(TCM2, OfNucleus);
1140 return nsol;
1141}
1142
1143
1144
1148
1150{
1151 // Return TF1 giving lab angle of nucleus as function of CM angle
1152 // OfNucleus = 1 or 2 (entrance channel) or 3 or 4 (exit channel)
1153
1154 if (!fThetaLabVsThetaCM[OfNucleus]) {
1155 fThetaLabVsThetaCM[OfNucleus] = new TF1(Form("KV2Body:ThetaLabVsThetaCM:%d", OfNucleus),
1156 const_cast<KV2Body*>(this), &KV2Body::ThetaLabVsThetaCM, 0, 180, 1, "KV2Body", "ThetaLabVsThetaCM");
1157 fThetaLabVsThetaCM[OfNucleus]->SetNpx(1000);
1158 fThetaLabVsThetaCM[OfNucleus]->SetParameter(0, OfNucleus);
1159 }
1160 return fThetaLabVsThetaCM[OfNucleus];
1161}
1162
1163
1164
1171
1173{
1174 //Calculate Rutherford cross-section (b/sr) in the CM as a
1175 //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1176 //
1177 // WARNING: in inverse kinematics, projectile lab angles generally have
1178 // two corresponding CM angles. We only use the most forward (smallest) CM angle.
1179
1180 Double_t par = OfNucleus;
1181 return const_cast<KV2Body*>(this)->XSecRuthCM(&ThetaLab, &par);
1182}
1183
1184
1185
1192
1194{
1195 //Calculate Rutherford cross-section (b/sr) in the CM as a
1196 //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1197 //
1198 // WARNING: in inverse kinematics, projectile lab angles generally have
1199 // two corresponding CM angles. We only use the most forward (smallest) CM angle.
1200
1201 if (!GetNucleus(2)) {
1202 Warning("GetXSecRuthCM", "No target defined for reaction");
1203 return 0.;
1204 }
1205 Double_t PB =
1206 1.44 * GetNucleus(1)->GetZ() * GetNucleus(2)->GetZ() /
1207 GetCMEnergy();
1208 // get projectile CM angle from lab angle of nucleus par[0]
1209 Double_t TCM = GetMinThetaCMFromThetaLab(1, x[0], par[0]);
1210 if (TCM < 0) return 0;
1211 Double_t D = 1. / (16. * (TMath::Power(TMath::Sin(TCM * TMath::DegToRad() / 2.), 4.)));
1212
1213 return ((TMath::Power(PB, 2.)) * D / 100.);
1214}
1215
1216
1217
1223
1225{
1226 // Return the smallest (i.e. most forward) CM angle of nucleus OfNucleus
1227 // corresponding to laboratory angle theta of nucleus OtherNucleus
1228 //
1229 // If theta>max lab angle for nucleus, returns -1.0
1230
1231 if (theta > TETAMAX[OtherNucleus]) return -1;
1232 Double_t t1, t2;
1233 Int_t nsol = GetThetaCM(theta, OtherNucleus, t1, t2);
1234 if (!nsol) return -1.;
1235 Double_t Pt1 = GetThetaCM(OfNucleus, t1, OtherNucleus);
1236 Double_t Pt2 = GetThetaCM(OfNucleus, t2, OtherNucleus);
1237 Double_t Pt;
1238 if (nsol > 1)
1239 Pt = TMath::Min(Pt1, Pt2);
1240 else
1241 Pt = Pt1;
1242 return Pt;
1243}
1244
1245
1246
1250
1252{
1253 // Calculate CM Rutherford cross-section (b/sr) in the CM as a function of
1254 // scattering angle in the CM frame for nucleus par[0]
1255
1256 if (!GetNucleus(2)) {
1257 Warning("XSecRuthCMVsThetaCM", "No target defined for reaction");
1258 return 0.;
1259 }
1260 Double_t PB =
1261 1.44 * GetNucleus(1)->GetZ() * GetNucleus(2)->GetZ() /
1262 GetCMEnergy();
1263 Double_t TCM = x[0] * TMath::DegToRad();
1264 Int_t OfNucleus = (Int_t)par[0];
1265 if (OfNucleus == 2 || OfNucleus == 4) TCM = TMath::Pi() - TCM; // get projectile CM angle from target CM angle
1266 Double_t D = 1. / (16. * (TMath::Power(TMath::Sin(TCM / 2.), 4.)));
1267 return ((TMath::Power(PB, 2.)) * D / 100.);
1268}
1269
1270
1271
1278
1280{
1281 //Calculate Rutherford cross-section (b/sr) in the Lab as a
1282 //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1283 //
1284 // WARNING: in inverse kinematics, projectile lab angles generally have
1285 // two corresponding CM angles. We only use the most forward (smallest) CM angle.
1286
1287 Double_t par = OfNucleus;
1288 return const_cast<KV2Body*>(this)->XSecRuthLab(&ThetaLab, &par);
1289}
1290
1291
1292
1299
1301{
1302 //Calculate Rutherford cross-section (b/sr) in the Lab as a
1303 //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1304 //
1305 // WARNING: in inverse kinematics, projectile lab angles generally have
1306 // two corresponding CM angles. We only use the most forward (smallest) CM angle.
1307
1308 Double_t DSIDTB = XSecRuthCM(x, par);
1309 if (DSIDTB > 0) {
1310 // get projectile CM angle from lab angle of nucleus par[0]
1311 Double_t T3CM = GetMinThetaCMFromThetaLab(par[0], x[0], par[0]);
1312 Double_t T3L = GetThetaLab(T3CM, par[0]) * TMath::DegToRad();
1313 T3CM *= TMath::DegToRad();
1314 Double_t RLC = (TMath::Power(TMath::Sin(T3CM), 3.)) /
1315 ((TMath::Power(TMath::Sin(T3L), 3.)) * GetCMGamma() *
1316 (1. + K[(int)par[0]] * TMath::Cos(T3CM)));
1317 if (DSIDTB * RLC < 0) {
1318 Warning("GetXSecRuthLab", "negative value for choosen parameters : %lf %d\n", x[0], (Int_t)par[0]);
1319 return 0;
1320 }
1321 DSIDTB *= RLC;
1322 }
1323 return DSIDTB;
1324}
1325
1326
1327
1332
1334{
1335 //Rutherford cross-section (b/sr) function in the Lab as a
1336 //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle x[0]
1337 //including 'sin theta' factor needed for integrating over solid angles.
1338
1339 return XSecRuthLab(x, par) * TMath::Sin(x[0] * TMath::DegToRad());
1340}
1341
1342
1343//Double_t KV2Body::GetIntegratedXSecRuthLab(KVTelescope* tel, Int_t OfNucleus)
1344//{
1345// //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1346// //polar and azimuthal angular range of the given KVTelescope.
1347// // if (OfNucleus==3) => X-section for scattered projectile
1348// // if (OfNucleus==4) => X-section for scattered target
1349// //
1350// //The returned value is in barns
1351//
1352// return GetIntegratedXSecRuthLab(tel->GetThetaMin(), tel->GetThetaMax(), tel->GetPhiMin(), tel->GetPhiMax(), OfNucleus);
1353//}
1354//
1356//Double_t KV2Body::GetIntegratedXSecRuthLab(KVDetector* det, Int_t OfNucleus)
1357//{
1358// //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1359// //polar and azimuthal angular range of the given KVDetector. These will be taken
1360// //from the parent KVTelescope of the detector.
1361// // if (OfNucleus==3) => X-section for scattered projectile
1362// // if (OfNucleus==4) => X-section for scattered target
1363// //
1364// //The returned value is in barns
1365//
1366// KVTelescope* tel = (KVTelescope*)det->GetParentStructure("TELESCOPE");
1367// if (!det) {
1368// Error("GetIntegratedXSecRuthLab(KVDetector*,Int_t)",
1369// "Detector has no parent telescope: it has not been positioned in a multidetector geometry");
1370// return 0;
1371// }
1372// return GetIntegratedXSecRuthLab(tel, OfNucleus);
1373//}
1374
1375
1386
1388{
1389 //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1390 //polar and azimuthal angular range expressed in degree
1391 // if (OfNucleus==3) This angular range is considered to be the scattered projectile one
1392 // if (OfNucleus==4) This angular range is considered to be the scattered target one
1393 //
1394 //If phi1 ou phi2 ==-1 the azimuthal width is set to 2pi
1395 //Else if phi1=phi2 the azimuthal width is set to 1 ie the integral is only on theta
1396 //
1397 //The returned value is in barns
1398
1399 if (th2 < th1) return 0;
1400 Double_t dphi = 0;
1401 //azimuthal width expressed in rad
1402 if (phi1 == -1 || phi2 == -1) dphi = 2 * TMath::Pi();
1403 else if (phi1 == phi2) dphi = 1;
1404 else {
1405 dphi = phi2 - phi1;
1406 dphi *= TMath::DegToRad();
1407 }
1408 Double_t theta_min = 1.;
1409 Double_t theta_max = 179.;
1410 if (th1 < theta_min) theta_min = th1;
1411 if (th2 > theta_max) theta_max = th2;
1412
1413#if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
1414 return GetXSecRuthLabIntegralFunc(OfNucleus, theta_min, theta_max)->Integral(th1, th2, fIntPrec) * TMath::DegToRad() * dphi;
1415#else
1416 const Double_t* para = 0;
1417 return GetXSecRuthLabIntegralFunc(OfNucleus, theta_min, theta_max)->Integral(th1, th2, para, fIntPrec) * TMath::DegToRad() * dphi;
1418#endif
1419}
1420
1421
1422
1427
1428TF1* KV2Body::GetXSecRuthCMFunc(Int_t OfNucleus, Double_t theta_cm_min, Double_t theta_cm_max) const
1429{
1430 // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the CM as a
1431 // function of projectile (OfNucleus=3) or target (OfNucleus=4) CM scattering angle
1432 // By default, theta_min = 1 degree & theta_max = 179 degrees
1433
1434 TString name = "RuthXSecCM: ";
1435 name += GetNucleus(1)->GetSymbol();
1436 name += " + ";
1437 name += GetNucleus(2)->GetSymbol();
1438 name += " ";
1439 Double_t elab = GetNucleus(1)->GetEnergy() / GetNucleus(1)->GetA();
1440 name += Form("%6.1f AMeV ", elab);
1441 if (OfNucleus == 3) name += "(projectile)";
1442 else name += "(target)";
1443 TF1* fXSecRuthCM = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1444 if (!fXSecRuthCM) {
1445 fXSecRuthCM = new TF1(name.Data(),
1446 const_cast<KV2Body*>(this), &KV2Body::XSecRuthCMVsThetaCM,
1447 theta_cm_min, theta_cm_max, 1, "KV2Body", "XSecRuthCMVsThetaCM");
1448 fXSecRuthCM->SetParameter(0, OfNucleus);
1449 fXSecRuthCM->SetNpx(1000);
1450 }
1451 fXSecRuthCM->SetRange(theta_cm_min, theta_cm_max); //in case TF1 already exists, but new range is required
1452 return fXSecRuthCM;
1453}
1454
1455
1456
1457
1462
1464{
1465 // calculate Bass interaction barrier B_int
1466 //
1467 // r0 = 1.07 fm
1468
1469 const Double_t r0 = 1.07;
1470 const Double_t e2 = 1.44;
1471 Double_t A1third = pow(GetNucleus(1)->GetA(), 1. / 3.);
1472 Double_t A2third = pow(GetNucleus(2)->GetA(), 1. / 3.);
1473 Double_t R12 = r0 * (A1third + A2third);
1474
1475 Double_t Bint = GetNucleus(1)->GetZ() * GetNucleus(2)->GetZ() * e2 / (R12 + 2.7)
1476 - 2.9 * A1third * A2third / (A1third + A2third);
1477
1478 return Bint;
1479}
1480
1481
1482
1483
1489
1491{
1492 // calculate Kox reaction X-section (in barns) for a given lab energy of projectile (in MeV/nucleon)
1493 //
1494 // c parameter fitted with Landau function vs. Log10(E/A) (see Fig. 12 of PRC Kox 87)
1495 // by imposing c=0.1 at E/A=10 MeV (goes to 0 as E/A->0)
1496
1497 const Double_t r0 = 1.05;
1498 const Double_t rc = 1.3;
1499 const Double_t e2 = 1.44;
1500 const Double_t a = 1.9;
1501
1502 KVNucleus* proj = GetNucleus(1);
1503 KVNucleus* targ = GetNucleus(2);
1504 proj->SetEnergy(eproj[0]*proj->GetA());
1506 Double_t ECM = GetCMEnergy();
1507
1508 Double_t c = 11.3315 * TMath::Landau(TMath::Log10(eproj[0]), 2.47498, 0.554937);
1509 Double_t A1third = pow(proj->GetA(), 1. / 3.);
1510 Double_t A2third = pow(targ->GetA(), 1. / 3.);
1511 Double_t BC = proj->GetZ() * targ->GetZ() * e2 / (rc * (A1third + A2third));
1512 Double_t EFac = 1 - BC / ECM;
1513
1514 Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1515 pow((A1third + A2third + a * A1third * A2third / (A1third + A2third) - c), 2) *
1516 EFac;
1517
1518 return Xsec / 100.;
1519}
1520
1521
1522
1523
1527
1529{
1530 // calculate Reaction Cross Section with the "Sphere Dure"
1531 // approximation
1532
1533 Double_t A1third = pow(GetNucleus(1)->GetA(), 1. / 3.);
1534 Double_t A2third = pow(GetNucleus(2)->GetA(), 1. / 3.);
1535
1536 Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1537 pow(A1third + A2third, 2);
1538
1539 return Xsec / 100.;
1540}
1541
1542
1543
1544
1548
1550{
1551
1552 // Deduce the maximum impact parameter (in fm) from a given reaction cross section (in barn)
1553 // in the approximation Xsec = \int(0,bmax) 2Pi b db
1554
1555 return 10 * TMath::Sqrt(ReacXsec / TMath::Pi());
1556
1557}
1558
1559
1560
1561
1562
1563
1567
1569{
1570
1571 // Integrate the cross section between two impact parameter (in fm)
1572 // and give the result in barn
1573
1574 return TMath::Pi() * (TMath::Power(b2, 2.) - TMath::Power(b1, 2.)) / 100;
1575
1576}
1577
1578
1579
1580
1586
1588{
1589 // Return pointer to TF1 with Kox reaction X-section in barns as a
1590 // function of projectile lab energy (in Mev/nucleon) for this reaction.
1591 // By default the range of the function is [20,100] MeV/nucleon.
1592 // Change with TF1::SetRange.
1593
1594 if (!fKoxReactionXSec) {
1596 name += " + ";
1597 name += GetNucleus(2)->GetSymbol();
1598 fKoxReactionXSec = new TF1(name.Data(),
1599 const_cast<KV2Body*>(this), &KV2Body::KoxReactionXSec, 20, 100, 0, "KV2Body", "KoxReactionXSec");
1600 fKoxReactionXSec->SetNpx(1000);
1601 }
1602 return fKoxReactionXSec;
1603}
1604
1605
1606
1607
1634
1636{
1637 // Calculate the mean charge state of the projectile after passage through the target,
1638 // assuming that the equilibrium charge state distribution is achieved*.
1639 // t[0] = energy of projectile after the target (in MeV/nucleon)
1640 //
1641 // We use the empirical parameterization of Leon et al., At. Dat. and Nucl. Dat. Tab. 69, 217 (1998)
1642 // developed for heavy ions in the GANIL energy range (it represents a fit to data measured using
1643 // GANIL beams).
1644 //
1645 // *N.B. Concerning equilibrium charge state distributions, it is not easy to know whether, for a given
1646 // combination of projectile, projectile energy, target, and target thickness, the equilibrium
1647 // distribution is reached or not. Here are some comments from the paper cited above which
1648 // may give some guidelines:
1649 //
1650 // "The energies available at the GANIL accelerator range from 24 to 95 MeV/u. Within this energy range,
1651 // the equilibrium charge state is reached only for fairly thick targets (~1 mg/cm2 for C foils)."
1652 //
1653 // "Mean Charge State as a Function of the Target Thickness
1654 // A typical example of the variation of the mean charge as a function of the foil thickness is shown ... It is seen
1655 // that the mean charge initially increases due to the ionization process. Then, the equilibrium state is reached at
1656 // a certain thickness, the so-called equilibrium thickness, due to the equilibration of electron loss and
1657 // capture processes. Finally, for foils thicker than the equilibrium thickness, the mean charge decreases due to the
1658 // slowing down of the ions in matter leading to higher capture cross sections."
1659 //
1660 // It should be noted that, according to the data published in this and other papers, the equilibrium thickness
1661 // decreases with increasing atomic number of the target, and increases with increasing energy of the projectile.
1662
1663 KVNucleus* proj = GetNucleus(1);
1664 Double_t Zp = proj->GetZ();
1665 proj->SetEnergy(t[0]*proj->GetA());
1666 Double_t beta = proj->Beta();
1667 Double_t vp = beta * KVParticle::C();
1668 Double_t Zt = GetNucleus(2)->GetZ();
1669
1670 Double_t q = Zp * (1. - TMath::Exp(-83.275 * beta / pow(Zp, 0.477)));
1671
1672 // correction for target Z
1673 Double_t g = 0.929 + 0.269 * TMath::Exp(-0.16 * Zt) + (0.022 - 0.249 * TMath::Exp(-0.322 * Zt)) * vp / pow(Zp, 0.477);
1674 q *= g;
1675
1676 if (Zp > 54) {
1677 // f(Zp) - correction for projectiles with Z>54
1678 Double_t f = 1. - TMath::Exp(-12.905 + 0.2124 * Zp - 0.00122 * pow(Zp, 2));
1679 q *= f;
1680 }
1681
1682 return q;
1683}
1684
1685
1686
1687
1697
1699{
1700 // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
1701 // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
1702 // energy after the target (in MeV/nucleon).
1703 // We use the empirical parameterization of Leon et al., At. Dat. and Nucl. Dat. Tab. 69, 217 (1998)
1704 // (see EqbmChargeState(Double_t *t,Double_t*) for details)
1705 //
1706 // By default the range of the function is [5,100] MeV/nucleon.
1707 // Change with TF1::SetRange.
1708
1709 if (!fEqbmChargeState) {
1711 name += " + ";
1712 name += GetNucleus(2)->GetSymbol();
1713 name += " LEON";
1714 fEqbmChargeState = new TF1(name.Data(),
1715 const_cast<KV2Body*>(this), &KV2Body::EqbmChargeState, 5, 100, 0, "KV2Body", "EqbmChargeState");
1716 fEqbmChargeState->SetNpx(1000);
1717 }
1718 return fEqbmChargeState;
1719}
1720
1721
1722
1726
1728{
1729 // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1730 // for solid targets
1731
1732 KVNucleus* proj = GetNucleus(1);
1733 Double_t Zp = proj->GetZ();
1734 proj->SetEnergy(t[0]*proj->GetA());
1735 KVNucleus* targ = GetNucleus(2);
1736 Double_t Zt = targ->GetZ();
1737 Double_t Vp = proj->Beta();
1738 const Double_t V0 = 2.19e+06 / TMath::C();
1739 Double_t x = -0.019 * pow(Zp, -0.52) * Vp / V0;
1740 x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x) / 1.68;
1741 x = pow(x, 1. + 1.8 / Zp);
1742 Double_t q = 0.07 / x + 6. + 0.3 * pow(x, 0.5) + 10.37 * x + pow(x, 4);
1743 q = (Zp * (12 * x + pow(x, 4.))) / q;
1744 return q;
1745}
1746
1747
1757
1759{
1760 // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
1761 // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
1762 // energy after the target (in MeV/nucleon).
1763 // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1764 // This formula is valid for solid targets.
1765 //
1766 // By default the range of the function is [5,100] MeV/nucleon.
1767 // Change with TF1::SetRange.
1768
1769 if (!fEqbmChargeStateShSol) {
1771 name += " + ";
1772 name += GetNucleus(2)->GetSymbol();
1773 name += " SHIWIETZ-SOLID";
1774 fEqbmChargeStateShSol = new TF1(name.Data(),
1775 const_cast<KV2Body*>(this), &KV2Body::eqbm_charge_state_shiwietz_solid, 5, 100, 0, "KV2Body", "eqbm_charge_state_shiwietz_solid");
1777 }
1778 return fEqbmChargeStateShSol;
1779}
1780
1781
1782
1786
1788{
1789 // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1790 // for gas targets
1791
1792 KVNucleus* proj = GetNucleus(1);
1793 Double_t Zp = proj->GetZ();
1794 proj->SetEnergy(t[0]*proj->GetA());
1795 KVNucleus* targ = GetNucleus(2);
1796 Double_t Zt = targ->GetZ();
1797 Double_t Vp = proj->Beta();
1798 const Double_t V0 = 2.19e+06 / TMath::C();
1799 Double_t x = 0.03 - 0.017 * pow(Zp, -0.52) * Vp / V0;
1800 x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x);
1801 x = pow(x, 1. + 0.4 / Zp);
1802 Double_t q = 1428. - 1206.*pow(x, 0.5) + 690 * x + pow(x, 6.);
1803 q = (Zp * (376.*x + pow(x, 6.))) / q;
1804 return q;
1805}
1806
1807
1808
1818
1820{
1821 // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
1822 // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
1823 // energy after the target (in MeV/nucleon).
1824 // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1825 // This formula is valid for gas targets.
1826 //
1827 // By default the range of the function is [5,100] MeV/nucleon.
1828 // Change with TF1::SetRange.
1829
1830 if (!fEqbmChargeStateShGas) {
1832 name += " + ";
1833 name += GetNucleus(2)->GetSymbol();
1834 name += " SHIWIETZ-GAS";
1835 fEqbmChargeStateShGas = new TF1(name.Data(),
1836 const_cast<KV2Body*>(this), &KV2Body::eqbm_charge_state_shiwietz_gas, 5, 100, 0, "KV2Body", "eqbm_charge_state_shiwietz_gas");
1838 }
1839 return fEqbmChargeStateShGas;
1840}
1841
1842
1843
1844
1849
1850TF1* KV2Body::GetXSecRuthLabFunc(Int_t OfNucleus, Double_t theta_min, Double_t theta_max) const
1851{
1852 // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the Lab as a
1853 // function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1854 // By default, theta_min = 1 degree & theta_max = 179 degrees
1855
1856 if (theta_max > GetMaxAngleLab(OfNucleus)) {
1857 theta_max = GetMaxAngleLab(OfNucleus);
1858 Info("GetXSecRuthLabFunc", "Maximum angle set to %lf degrees", theta_max);
1859 }
1860
1861 TString name = "RuthXSec: ";
1862 name += GetNucleus(1)->GetSymbol();
1863 name += " + ";
1864 name += GetNucleus(2)->GetSymbol();
1865 name += " ";
1866 Double_t elab = GetNucleus(1)->GetEnergy() / GetNucleus(1)->GetA();
1867 name += Form("%6.1f AMeV ", elab);
1868 if (OfNucleus == 3) name += "(projectile)";
1869 else name += "(target)";
1870 TF1* fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1871 if (!fXSecRuthLab) {
1872 fXSecRuthLab = new TF1(name.Data(),
1873 const_cast<KV2Body*>(this), &KV2Body::XSecRuthLab, theta_min, theta_max, 1, "KV2Body", "XSecRuthLab");
1874 fXSecRuthLab->SetParameter(0, OfNucleus);
1875 fXSecRuthLab->SetNpx(1000);
1876 }
1877 fXSecRuthLab->SetRange(theta_min, theta_max); //in case TF1 already exists, but new range is required
1878 return fXSecRuthLab;
1879}
1880
1881
1882
1895
1896TF1* KV2Body::GetXSecRuthLabIntegralFunc(Int_t OfNucleus, Double_t theta_min, Double_t theta_max) const
1897{
1898 // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the Lab as a
1899 // function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1900 //
1901 // This function is equal to sin(theta)*dsigma/domega, i.e. it is the integrand
1902 // needed for calculating total cross-sections integrated over solid angle.
1903 //
1904 // WARNING: when integrating this function using TF1::Integral, the result must
1905 // be multiplied by TMath::DegToRad(), because 'x' is in degrees rather than radians,
1906 // e.g. the integrated cross-section in barns is given by
1907 //
1908 // GetXSecRuthLabIntegralFunc(OfNucleus)->Integral(theta_min, theta_max)*TMath::DegToRad()
1909
1910 TString name = "RuthXSecInt: ";
1911 name += GetNucleus(1)->GetSymbol();
1912 name += " + ";
1913 name += GetNucleus(2)->GetSymbol();
1914 name += " ";
1915 Double_t elab = GetNucleus(1)->GetEnergy() / GetNucleus(1)->GetA();
1916 name += Form("%6.1f AMeV ", elab);
1917 if (OfNucleus == 3) name += "(projectile)";
1918 else name += "(target)";
1919
1920 TF1* fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1921 if (!fXSecRuthLab) {
1922 fXSecRuthLab = new TF1(name.Data(),
1923 const_cast<KV2Body*>(this), &KV2Body::XSecRuthLabInt, theta_min, theta_max, 1, "KV2Body", "XSecRuthLabInt");
1924 fXSecRuthLab->SetParameter(0, OfNucleus);
1925 fXSecRuthLab->SetNpx(1000);
1926 }
1927 fXSecRuthLab->SetRange(theta_min, theta_max); //in case TF1 already exists, but new range is required
1928 return fXSecRuthLab;
1929}
1930
1931
1932
1933
1943
1945{
1946 // Find at most two solutions x1 and x2 between xmin and xmax for which fonc->Eval(x) = val
1947 // i.e. we use TF1::GetX for a function which may be (at most) double-valued in range (xmin, xmax)
1948 // This is adapted to the case of the lab angle vs. CM angle of scattered particles.
1949 // if fonc has a maximum between xmin and xmax, and if val < max, we look for x1 between
1950 // xmin and the maximum, and x2 between the maximum and xmax.
1951 // If not (single-valued function) we look for x1 between xmin and xmax.
1952 // This method returns the number of roots found.
1953 // If val > maximum of fonc between xmin and xmax, there is no solution: we return 0.
1954
1955 x1 = x2 = xmin - 1.;
1956 Double_t max = fonc->GetMaximum(xmin, xmax);
1957 Int_t nRoots = 0;
1958 if (val > max) return nRoots;
1959 Double_t maxX = fonc->GetMaximumX(xmin, xmax);
1960 nRoots = 1;
1961 if (TMath::AreEqualAbs(val, max, 1.e-10)) {
1962 // value corresponds to maximum of function
1963 x1 = maxX;
1964 return nRoots;
1965 }
1966 // 2 roots if 'fonc' has a maximum xmin<maxX<xmax and if fonc(xmin)<val && fonc(xmax)<val
1967 if ((maxX < xmax) && (maxX > xmin) && (fonc->Eval(xmin) < val) && (fonc->Eval(xmax) < val)) nRoots = 2;
1968 Double_t xmax1 = (nRoots == 2 ? maxX : xmax);
1969 x1 = fonc->GetX(val, xmin, xmax1);
1970 if (nRoots == 1) return nRoots;
1971 else
1972 x2 = fonc->GetX(val, maxX, xmax);
1973 return nRoots;
1974}
1975
1976
1977
1978
1979
int Int_t
#define f(i)
#define c(i)
#define e(i)
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
#define X(type, name)
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
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 g
char name[80]
float xmin
float * q
float xmax
#define gROOT
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
Double_t XSecRuthLabInt(Double_t *, Double_t *)
Definition KV2Body.cpp:1333
Double_t TETAMIN[5]
defined only for nuclei 3 et 4
Definition KV2Body.h:181
Double_t eqbm_charge_state_shiwietz_gas(Double_t *t, Double_t *)
Definition KV2Body.cpp:1787
virtual ~KV2Body()
Definition KV2Body.cpp:362
Double_t VC[5]
cm velocities
Definition KV2Body.h:177
void SetTarget(const KVNucleus &)
Set target for reaction.
Definition KV2Body.cpp:314
TF1 * GetEqbmChargeStateFunc() const
Definition KV2Body.cpp:1698
void Set4thNucleus()
Definition KV2Body.cpp:422
TF1 * fELabVsThetaLab[5]
Definition KV2Body.h:199
TF1 * GetXSecRuthLabFunc(Int_t OfNucleus=3, Double_t theta_min=1., Double_t theta_max=179.) const
Definition KV2Body.cpp:1850
Double_t K[5]
ratio of c.m. velocity to velocity of nucleus in c.m. v_cm/v_i_cm
Definition KV2Body.h:179
Double_t GetMaxAngleLab(Int_t i) const
Definition KV2Body.cpp:536
std::vector< KVNucleus > fNuclei
nuclei involved in calculation
Definition KV2Body.h:168
void Print(Option_t *opt="") const
Definition KV2Body.cpp:813
TF1 * GetXSecRuthLabIntegralFunc(Int_t OfNucleus=3, Double_t theta_min=1., Double_t theta_max=179.) const
Definition KV2Body.cpp:1896
Double_t WCT
total cm energy
Definition KV2Body.h:174
Double_t GetSphereDureReactionXSec(Double_t r0=1.05)
Definition KV2Body.cpp:1528
Double_t ELabVsThetaLab(Double_t *, Double_t *)
Function calculating lab energy of nucleus par[0] for any lab angle x[0].
Definition KV2Body.cpp:989
TF1 * fThetaLabVsThetaCM[5]
Definition KV2Body.h:197
Double_t GetMinAngleLab(Int_t i) const
Definition KV2Body.cpp:554
Double_t GetThetaLab(Double_t ThetaCM, Int_t OfNucleus) const
Definition KV2Body.h:286
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Definition KV2Body.cpp:339
Double_t fIntPrec
Precision of the TF1::Integral method.
Definition KV2Body.h:207
Double_t GetXSecRuthLab(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Definition KV2Body.cpp:1279
Double_t GetQReaction() const
Calculate Q-value for reaction, including dissipated (excitation) energy.
Definition KV2Body.cpp:480
TF1 * fEqbmChargeStateShSol
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al solid)
Definition KV2Body.h:185
Double_t EqbmChargeState(Double_t *t, Double_t *)
Definition KV2Body.cpp:1635
TF1 * GetThetaLabVsThetaCMFunc(Int_t OfNucleus) const
Definition KV2Body.cpp:1149
void SetOutgoing(const KVNucleus &proj_out)
Definition KV2Body.cpp:397
TF1 * fKoxReactionXSec
function Kox reaction cross-section [barns] vs. E/A projectile
Definition KV2Body.h:183
TF1 * GetELabVsThetaCMFunc(Int_t OfNucleus) const
Definition KV2Body.cpp:1044
Double_t XSecRuthCM(Double_t *, Double_t *)
Definition KV2Body.cpp:1193
TF1 * GetShiwietzEqbmChargeStateFuncForGasTargets() const
Definition KV2Body.cpp:1819
Double_t XSecRuthLab(Double_t *, Double_t *)
Definition KV2Body.cpp:1300
TF1 * GetKoxReactionXSecFunc() const
Definition KV2Body.cpp:1587
Double_t BCM
beta of centre of mass
Definition KV2Body.h:172
TF1 * fEqbmChargeStateShGas
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al gas)
Definition KV2Body.h:186
Double_t GetELab(Double_t ThetaCM, Int_t OfNucleus) const
Definition KV2Body.h:291
Int_t FindRoots(TF1 *, Double_t, Double_t, Double_t, Double_t &, Double_t &) const
Definition KV2Body.cpp:1944
Double_t GetEDiss() const
Definition KV2Body.h:251
static Double_t GetVelocity(Double_t mass, Double_t E)
Definition KV2Body.cpp:380
Double_t GetIntegratedXSecRuthLab(Float_t th1, Float_t th2, Float_t phi1=-1, Float_t phi2=-1, Int_t OfNucleus=3)
Definition KV2Body.cpp:1387
Double_t GetBmaxFromReactionXSec(Double_t ReacXsec)
Definition KV2Body.cpp:1549
TVector3 VCM
velocity of centre of mass
Definition KV2Body.h:171
Double_t KoxReactionXSec(Double_t *, Double_t *)
Definition KV2Body.cpp:1490
Double_t ELabVsThetaCM(Double_t *, Double_t *)
Function calculating lab energy of nucleus par[0] for any CM angle x[0].
Definition KV2Body.cpp:1006
Double_t ThetaLabVsThetaCM(Double_t *, Double_t *)
Definition KV2Body.cpp:1025
Double_t GetCMGamma() const
Definition KV2Body.h:269
Double_t TETAMAX[5]
defined only for nuclei 3 et 4
Definition KV2Body.h:180
Bool_t fSetOutgoing
= kTRUE if SetOutgoing is called before CalculateKinematics
Definition KV2Body.h:204
Double_t WC[5]
cm energy of each nucleus
Definition KV2Body.h:175
Double_t EC[5]
cm energies
Definition KV2Body.h:178
Double_t GetExcitEnergy() const
Definition KV2Body.h:241
Double_t GetMinThetaCMFromThetaLab(Int_t OfNucleus, Double_t theta, Int_t OtherNucleus) const
Definition KV2Body.cpp:1224
Double_t WLT
total lab energy
Definition KV2Body.h:173
Int_t GetVLab(Int_t OfNucleus, Double_t ThetaLab, Int_t AngleNucleus, Double_t &e1, Double_t &e2) const
Definition KV2Body.cpp:1088
TF1 * fXSecRuthLab[5]
Definition KV2Body.h:202
Double_t BassIntBarrier()
Definition KV2Body.cpp:1463
TVector3 GetCMVelocity() const
Return vector velocity of centre of mass of reaction (units: cm/ns)
Definition KV2Body.cpp:572
void CalculateKinematics()
Definition KV2Body.cpp:677
void init()
Default initialisations.
Definition KV2Body.cpp:32
KV2Body()
default ctor
Definition KV2Body.cpp:61
Double_t GetIntegratedXsec(Double_t b1, Double_t b2)
Definition KV2Body.cpp:1568
TF1 * GetXSecRuthCMFunc(Int_t OfNucleus=3, Double_t theta_cm_min=1., Double_t theta_cm_max=179.) const
Definition KV2Body.cpp:1428
Double_t GetLabGrazingAngle(Int_t i=1) const
Definition KV2Body.cpp:603
Double_t GetQGroundStates() const
Calculate Q-value for reaction, assuming all nuclei in ground state.
Definition KV2Body.cpp:499
KVNucleus * GetNucleus(Int_t i) const
Definition KV2Body.cpp:456
Double_t XSecRuthCMVsThetaCM(Double_t *, Double_t *)
Definition KV2Body.cpp:1251
Double_t GetCMEnergy() const
Return available kinetic energy in centre of mass.
Definition KV2Body.cpp:522
TF1 * fEqbmChargeState
function equilibrium charge state of projectile vs. E/A projectile (Leon et al)
Definition KV2Body.h:184
Double_t eqbm_charge_state_shiwietz_solid(Double_t *t, Double_t *)
Definition KV2Body.cpp:1727
Double_t GetXSecRuthCM(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Definition KV2Body.cpp:1172
Double_t fEDiss
dissipated energy, 0 means elastic scattering
Definition KV2Body.h:169
TF1 * fELabVsThetaCM[5]
Definition KV2Body.h:198
Int_t GetThetaCM(Double_t ThetaLab, Int_t OfNucleus, Double_t &t1, Double_t &t2) const
Definition KV2Body.cpp:972
TF1 * GetShiwietzEqbmChargeStateFuncForSolidTargets() const
Definition KV2Body.cpp:1758
TF1 * GetELabVsThetaLabFunc(Int_t OfNucleus) const
Definition KV2Body.cpp:1064
Utility class for kinematical transformations of KVParticle class.
Description of properties and kinematics of atomic nuclei.
Definition KVNucleus.h:126
const Char_t * GetSymbol(Option_t *opt="") const
Definition KVNucleus.cpp:81
Double_t GetExcitEnergy() const
Definition KVNucleus.h:283
Double_t GetMassExcess(Int_t z=-1, Int_t a=-1) const
Bool_t IsDefined() const
Definition KVNucleus.h:202
Int_t GetA() const
Int_t GetZ() const
Return the number of proton / atomic number.
TVector3 GetMomentum() const
Definition KVParticle.h:604
Double_t GetEnergy() const
Definition KVParticle.h:621
static Double_t C()
Double_t GetKE() const
Definition KVParticle.h:614
void SetFrame(const Char_t *frame, const KVFrameTransform &)
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
TVector3 GetV() const
Definition KVParticle.h:671
void SetEnergy(Double_t e)
Definition KVParticle.h:599
Double_t GetMass() const
Definition KVParticle.h:574
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73
void Begin(TString delim) const
Definition KVString.cpp:565
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition KVString.cpp:695
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
virtual void SetNpx(Int_t npx=100)
virtual Double_t GetX(Double_t y, Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
virtual Double_t GetMaximum(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual void SetParameter(const TString &name, Double_t value)
virtual Double_t GetMaximumX(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
Double_t Beta() const
Double_t E() const
virtual void SetTitle(const char *title="")
virtual void Warning(const char *method, const char *msgfmt,...) const
R__ALWAYS_INLINE Bool_t IsZombie() const
virtual void Error(const char *method, const char *msgfmt,...) const
void MakeZombie()
virtual void Info(const char *method, const char *msgfmt,...) const
void Obsolete(const char *method, const char *asOfVers, const char *removedFromVers) const
Double_t Atof() const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
void SetXYZ(Double_t x, Double_t y, Double_t z)
Double_t Mag() const
void compound()
double beta(double x, double y)
T Mag(const SVector< T, D > &rhs)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
Double_t x[n]
double max(double x, double y)
Double_t Min(Double_t a, Double_t b)
constexpr Double_t C()
Double_t ASin(Double_t)
Double_t Exp(Double_t x)
Double_t ATan(Double_t)
constexpr Double_t K()
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
constexpr Double_t PiOver2()
Double_t ATan2(Double_t y, Double_t x)
Double_t Power(Double_t x, Double_t y)
constexpr Double_t E()
constexpr Double_t DegToRad()
Double_t Sqrt(Double_t x)
Double_t Cos(Double_t)
constexpr Double_t Pi()
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Double_t Abs(Double_t d)
Double_t Sin(Double_t)
constexpr Double_t RadToDeg()
Double_t Log10(Double_t x)
v2
v1
auto * th2
auto * th1
TArc a
auto * t1
ClassImp(TPyArg)