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