KaliVeda
Toolkit for HIC analysis
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 
23 using namespace std;
24 
26 
27 
28 
29 
32 void 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 
61 KV2Body::KV2Body(): fNuclei(5), fEDiss(0)
62 {
63  //default ctor
64  init();
65 }
66 
67 
68 
88 
89 KV2Body::KV2Body(const Char_t* systemname) : fNuclei(5), fEDiss(0)
90 {
91  //Set up calculation defining entrance channel following this prescription:
92  //
93  // [Projectile_Symbol]+[Target_Symbol]@[Incident_Energy]MeV/A
94  // [Projectile_Symbol]+[Target_Symbol]@[Incident_Energy]MeV/u
95  // [Projectile_Symbol]+[Target_Symbol]@[Incident_Energy]MeV
96  // [Projectile_Symbol]+[Target_Symbol]
97  //
98  //Any spaces will be ignored.
99  //
100  //Valid examples:
101  //
102  // 129Xe+119Sn@50.0MeV/A
103  // 58Ni + 64Ni @ 32 MeV/u
104  // 12C + 1H @ 132 MeV
105  // U+U@5MeV/A
106  // Ta+Zn
107  //
108  //If the format is not respected, this object will be a zombie (IsZombie() returns kTRUE)
109 
110  init();
111  KVString sener = "";
112  Double_t ee = 0;
113  bool beam_in_mev = false;
114  KVString syst(systemname);
115  syst.ReplaceAll(" ", "");
116  KVString scouple(syst);
117 
118  if (syst.Contains("@")) {
119  syst.Begin("@");
120  scouple = syst.Next();
121  sener = syst.Next();
122  }
123  if (sener != "") {
124  if (sener.Contains("MeV/A")) {
125  sener.ReplaceAll("MeV/A", "");
126  ee = sener.Atof();
127  }
128  else if (sener.Contains("MeV/u")) {
129  sener.ReplaceAll("MeV/u", "");
130  ee = sener.Atof();
131  }
132  else if (sener.Contains("MeV")) {
133  sener.ReplaceAll("MeV", "");
134  ee = sener.Atof();
135  beam_in_mev = true;
136  }
137  else {
138  MakeZombie();
139  }
140  }
141  scouple.Begin("+");
142  KVNucleus nnuc1(scouple.Next());
143  if(beam_in_mev)
144  nnuc1.SetEnergy(ee);
145  else
146  nnuc1.SetEnergy(ee*nnuc1.GetA());
147  if (nnuc1.IsZombie()) MakeZombie();
148  else {
149  fNuclei[1] = nnuc1;
150  KVNucleus nnuc2(scouple.Next());
151  if (nnuc2.IsZombie()) MakeZombie();
152  else {
153  fNuclei[2] = nnuc2;
154  }
155  }
156  if (!IsZombie()) {
157  SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), ee));
158  }
159 }
160 
161 
162 
165 
166 KV2Body::KV2Body(KVNucleus*, KVNucleus* cib, KVNucleus* proj_out, Double_t): fNuclei(5), fEDiss(0)
167 {
168  // Deprecated. Do not use.
169 
170  Obsolete("KV2Body(KVNucleus*, KVNucleus*, KVNucleus*, Double_t)", "1.11", "1.12");
171 
172  if (!cib) {
173  Info("KV2Body", "Use constructor KV2Body(const KVNucleus& compound) to simulate binary decay of compound");
174  }
175  else if (!proj_out) {
176  Info("KV2Body", "Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ) to define entrance channel");
177  }
178  else {
179  Info("KV2Body", "Use constructor KV2Body(const KVNucleus& proj, const KVNucleus& targ, const KVNucleus& outgoing) to define entrance & exit channels");
180  }
181  init();
182 }
183 
184 
185 
209 
210 KV2Body::KV2Body(const KVNucleus& compound, double Exx)
211  : fNuclei(5), fEDiss(-Exx)
212 {
213  // Set up kinematics of binary decay of compound nucleus.
214  //
215  // The excitation energy of the CN (in MeV) is either given:
216  //
217  // - in the KVNucleus object itself (i.e. compound.GetExcitEnergy() gives the E* of CN)
218  // - or in the variable Exx: in this case any E* in the KVNucleus will be ignored
219  // - or by calling SetExcitEnergy or SetDissipatedEnergy after this constructor with the NEGATIVE E*
220  //
221  //Usage:
222  //~~~~~~~~~~{.cpp}
223  // KVNucleus CN("204Pb");
224  // CN.SetExcitEnergy(1.5*CN.GetA());
225  // KV2Body CNdecay(CN);
226  //
227  // /* or: KV2Body CNdecay("204Pb", 1.5*204);
228  // * or: KV2Body CNdecay("204Pb");
229  // * CNdecay.SetExcitEnergy(-1.5*204); */
230  //
231  // // calculate decay
232  // KVNucleus alpha("4He", 67.351/4.);
233  // CNdecay.SetOutgoing(alpha);
234  //~~~~~~~~~~
235 
236  init();
237  fNuclei[1] = compound;
238  if (fEDiss == 0) fEDiss = compound.GetExcitEnergy();
239 }
240 
241 
242 
255 
256 KV2Body::KV2Body(const KVNucleus& proj, const KVNucleus& targ, double Ediss):
257  fNuclei(5), fEDiss(Ediss)
258 {
259  //Set up calculation of basic binary reaction for given projectile and target.
260  //
261  //By default the dissipated energy Ediss is zero (elastic reaction).
262  //
263  //Usage:
264  //~~~~~~~~~~{.cpp}
265  // KV2Body reaction(KVNucleus("129Xe",50), "natSn");
266  //
267  // // or with C++11 (ROOT6):
268  // KV2Body reaction({"129Xe",50}, "natSn");
269  //~~~~~~~~~~
270 
271  init();
272  fNuclei[1] = proj;
273  fNuclei[2] = targ;
274  SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
275 }
276 
277 
278 
295 
296 KV2Body::KV2Body(const KVNucleus& proj, const KVNucleus& targ, const KVNucleus& outgoing, double Ediss):
297  fNuclei(5), fEDiss(Ediss)
298 {
299  //Set up calculation of basic binary reaction for given projectile and target with definition
300  //of exit channel (outgoing projectile- or target-like fragment).
301  //
302  //By default the dissipated energy is zero (elastic reaction).
303  //
304  //Any excitation energy of the outgoing fragment will be taken into account, e.g.
305  //for quasi-elastic scattering leaving the target in an excited state.
306  //
307  //Usage:
308  //~~~~~~~~~~{.cpp}
309  // KV2Body reaction(KVNucleus("129Xe",50), "natSn", KVNucleus("24Mg",35));
310  //
311  // // or with C++11 (ROOT6):
312  // KV2Body reaction({"129Xe",50}, "natSn", {"24Mg",35});
313  //~~~~~~~~~~
314 
315  init();
316  fNuclei[1] = proj;
317  fNuclei[2] = targ;
318  if (outgoing.GetExcitEnergy() > 0) fEDiss += outgoing.GetExcitEnergy();
319  SetOutgoing(outgoing);
320  SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
321 }
322 
323 
324 
327 
328 void KV2Body::SetTarget(const KVNucleus& targ)
329 {
330  //Set target for reaction.
331  fNuclei[2] = targ;
333 }
334 
335 
336 
339 
341 {
342  //Set target for reaction
343 
344  fNuclei[2].SetZandA(z, a);
346 }
347 
348 
349 
352 
354 {
355  //Set projectile for reaction.
356  fNuclei[1] = proj;
358 }
359 
360 
361 
364 
366 {
367  //Set projectile for reaction
368  fNuclei[1].SetZandA(z, a);
370 }
371 
372 
373 
375 
376 KV2Body::~KV2Body()
377 {
380  for (int i = 0; i < 5; i++) {
381  if (fThetaLabVsThetaCM[i]) delete fThetaLabVsThetaCM[i];
382  if (fELabVsThetaCM[i]) delete fELabVsThetaCM[i];
383  if (fELabVsThetaLab[i]) delete fELabVsThetaLab[i];
384  }
385 }
386 
387 
388 
389 
393 
395 {
396  //Static function, calculates relativistic velocity (in cm/ns) from rest mass and total
397  //energy (i.e. KE + mass) 'E' for any particle
398 
399  Double_t p = TMath::Sqrt(E * E - mass * mass);
400  return KVParticle::C() * p / E;
401 }
402 
403 
404 
405 
410 
411 void KV2Body::SetOutgoing(const KVNucleus& proj_out)
412 {
413  // Set outgoing projectile-like nucleus properties.
414  //
415  // The properties of the outgoing target-like nucleus will be deduced from mass, charge and momentum/energy conservation.
416 
417  SetTitle(Form("%s + %s %.1f MeV/A", fNuclei[1].GetSymbol(), fNuclei[2].GetSymbol(), fNuclei[1].GetAMeV()));
419  fNuclei[3] = proj_out;
420  Set4thNucleus();
421  // all nuclei should now have E* set to zero in order to use ground state masses in kinematics
422  for (unsigned i = 1; i <= 4; ++i) if (fNuclei[i].IsDefined()) fNuclei[i].SetExcitEnergy(0);
423 }
424 
425 
426 
427 
435 
437 {
438  // Private method, used to deduce 4th nucleus (target-like) from projectile, target
439  // and outgoing projectile using conservation of mass, momentum and energy.
440  //
441  // if the outgoing nucleus set by the user is equal to the compound nucleus
442  // formed by projectile and target, but the excitation energy of the CN
443  // was not set by the user, we calculate it here.
444 
445  KVNucleus sum;
446  if (GetNucleus(target))
448  else sum = *GetNucleus(projectile);
450  if (!tmp4.IsDefined()) {
451  // nucleus 3 is the CN proj+targ
452  // has E* been defined ?
453  if (fEDiss == 0) fEDiss = sum.GetExcitEnergy();
454  }
455  fNuclei[4] = tmp4;
456 }
457 
458 
459 
460 
467 
469 {
470  //Return pointer to nucleus
471  // \param[in] i nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
472  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
473  //
474  // Will return `nullptr` if any nucleus is undefined
475 
476  if (i >= projectile && i <= target_like) {
477  if (fNuclei[i].IsDefined()) return (KVNucleus*) &fNuclei[i];
478  return nullptr;
479  }
480  Warning("GetNucleus(nucleus_of_interest i)", "Index i out of bounds, i=%d", i);
481  return nullptr;
482 }
483 
484 
485 
486 
489 
491 {
492  //Calculate Q-value for reaction, including dissipated (excitation) energy
493 
494  if (!fSetOutgoing) {
495  Warning("GetQReaction", "Parameters for outgoing nuclei not set");
496  return 0.0;
497  }
499 
500  return QR;
501 }
502 
503 
504 
505 
508 
510 {
511  //Calculate Q-value for reaction, assuming all nuclei in ground state
512 
513  if (!fSetOutgoing) {
514  Warning("GetQGroundStates",
515  "Parameters for outgoing nuclei not set");
516  return 0.0;
517  }
518  Double_t QGG =
523  return QGG;
524 }
525 
526 
527 
528 
531 
533 {
534  //Return available kinetic energy in centre of mass
535 
536  return WCT - (GetNucleus(projectile)->GetMass() +
537  (GetNucleus(target) ? GetNucleus(target)->GetMass() : 0.));
538 }
539 
540 
541 
545 
547 {
548  //Returns maximum scattering angle in lab
549  // \param[in] i either nucleus_of_interest::projectile_like or nucleus_of_interest::target_like
550 
551  if (i < projectile_like || i > target_like) {
552  Warning("GetMaxAngleLab(nucleus_of_interest i)",
553  "Returns maximum scattering angle in lab for projectile_like or target_like nuclei");
554  return 0.;
555  }
556  return TETAMAX[i];
557 }
558 
559 
560 
564 
566 {
567  //Returns minimum scattering angle in lab
568  // \param[in] i either nucleus_of_interest::projectile_like or nucleus_of_interest::target_like
569  if (i < projectile_like || i > target_like) {
570  Warning("GetMinAngleLab(Int_t i)",
571  "Returns minimum scattering angle in lab for projectile_like or target_like nuclei");
572  return 0.;
573  }
574  return TETAMIN[i];
575 }
576 
577 
578 
579 
582 
584 {
585  //Return vector velocity of centre of mass of reaction (units: cm/ns)
586 
587  return VCM;
588 }
589 
590 
591 
592 
597 
599 {
600  //Return velocity of nucleus in the centre of mass of the reaction
601  // \param[in] i nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
602  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
603  return VC[i];
604 }
605 
606 
607 
608 
612 
614 {
615  //Calculate laboratory grazing angle.
616  // \param[in] i nucleus_of_interest::projectile or nucleus_of_interest::target
617 
618  if (i < projectile || i > target) {
619  Warning("GetLabGrazingAngle(nucleus_of_interest i)",
620  "i should be projectile or target");
621  return 0.0;
622  }
623  if (!GetNucleus(target)) {
624  Warning("GetLabGrazingAngle(nucleus_of_interest i)",
625  "No target defined for reaction");
626  return 0.0;
627  }
628  Double_t RP = 1.28 * TMath::Power(GetNucleus(projectile)->GetA(), 1. / 3.)
629  - 0.76 + 0.8 * TMath::Power(GetNucleus(projectile)->GetA(), -1. / 3.);
630  Double_t RT = 1.28 * TMath::Power(GetNucleus(target)->GetA(), 1. / 3.)
631  - 0.76 + 0.8 * TMath::Power(GetNucleus(target)->GetA(), -1. / 3.);
632  Double_t CP = RP * (1. - TMath::Power(1. / RP, 2.));
633  Double_t CT = RT * (1. - TMath::Power(1. / RT, 2.));
634  Double_t RINT = CP + CT + 4.49 - (CT + CP) / 6.35;
635  Double_t BAR =
636  1.44 * GetNucleus(projectile)->GetZ() * GetNucleus(target)->GetZ() / RINT;
637  if (GetCMEnergy() < BAR) {
638  //below Coulomb barrier
639  switch (i) {
640  case projectile:
641  return 180.;
642  break;
643  case target:
644  return 90.;
645  break;
646  case projectile_like:
647  case target_like:
648  break;
649  }
650  }
651  Double_t X = 2. * RINT * GetCMEnergy();
652  Double_t Y =
653  X / (GetNucleus(projectile)->GetZ() * GetNucleus(target)->GetZ() * 1.44) - 1.;
654  Double_t U = 1. / Y;
655  Double_t GR = 2. * TMath::ASin(U);
656  Double_t GPART = TMath::Pi() - GR;
657  Double_t ARAZ =
658  TMath::Sin(GR) / (TMath::Cos(GR) +
659  GetNucleus(projectile)->GetA() / (1.0 *
660  GetNucleus(target)->GetA()));
661  Double_t GRAZ = TMath::ATan(ARAZ);
662  if (GRAZ <= 0.)
663  GRAZ += TMath::Pi();
664 
665  switch (i) {
666  case projectile: //projectile
667  return (GRAZ * TMath::RadToDeg());
668  break;
669  case target: //target
670  return (GPART * TMath::RadToDeg() / 2.);
671  break;
672  case projectile_like:
673  case target_like:
674  default:
675  break;
676  }
677  return -99.;
678 }
679 
680 
681 
682 
692 
694 {
695  // Called to make kinematics calculation for all nuclei, use once properties of
696  // entrance-channel (and, if necessary, exit-channel) nuclei have been defined.
697  //
698  // If a compound decay is to be calculated (only 1 nucleus in entrance channel),
699  // outgoing (decay) product must be defined before calling this method.
700  //
701  // For a 2-body entrance channel, if no exit-channel nuclei are defined,
702  // we assume elastic scattering (i.e. identical outgoing nuclei)
703 
705  KVNucleus* Nuc2 = GetNucleus(target);
706 
707  // compound decay ?
708  if (!Nuc2 && !fSetOutgoing) {
709  Error("CalculateKinematics", "Set outgoing (decay) product first");
710  return;
711  }
712 
713  // call SetOutgoing if not already done
714  if (!fSetOutgoing) SetOutgoing(*Nuc1);
715 
716  // set everything to zero
717  WLT = WCT = BCM = 0.;
718  for (int i = 0; i < 5; i++) {
719  WC[i] = 0.;
720  VC[i] = 0.;
721  EC[i] = 0.;
722  K[i] = 0.;
723  TETAMAX[i] = 0.;
724  TETAMIN[i] = 0.;
725  }
726  VCM.SetXYZ(0., 0., 0.);
727 
728  //total energy (T + m) of entrance channel
729  WLT = Nuc1->E();
730  //velocity of CM
731  VCM = Nuc1->GetMomentum();
732  if (Nuc2) {
733  WLT += Nuc2->E();
734  VCM += Nuc2->GetMomentum();
735  }
736 
737  //beta of CM
738  BCM = VCM.Mag() / WLT;
739 
740  VCM *= (1. / WLT) * KVParticle::C();
741  Nuc1->SetFrame("CM", KVFrameTransform(VCM, kFALSE));
742  if (Nuc2)
743  Nuc2->SetFrame("CM", KVFrameTransform(VCM, kFALSE));
744 
745  //total energy in CM
746  WCT = (GetCMGamma() > 0.0 ? WLT / GetCMGamma() : 0.);
747  if (WCT == 0.0)
748  return;
749 
750  //total energy proj in CM
751  WC[1] = Nuc1->GetFrame("CM")->E();
752  //kinetic energy proj in CM
753  EC[1] = Nuc1->GetFrame("CM")->GetKE();
754  VC[1] = Nuc1->GetFrame("CM")->GetVelocity().Mag();
755  K[1] = VCM.Mag() / VC[1];
756  if (Nuc2) {
757  WC[2] = Nuc2->GetFrame("CM")->E();
758  EC[2] = Nuc2->GetFrame("CM")->GetKE();
759  VC[2] = Nuc2->GetFrame("CM")->GetVelocity().Mag();
760  K[2] = VCM.Mag() / VC[2];
761  }
762 
765  //total cm energy of nucleus 3 (quasiproj)
766  WC[3] = (WCT - GetEDiss()) / 2. + (AM3 - AM4) * (AM3 + AM4)
767  / (2. * (WCT - GetEDiss()));
768  VC[3] = GetVelocity(AM3, WC[3]);
769  //cm kinetic energy
770  EC[3] = WC[3] - AM3;
771  if (VC[3] > 0.) K[3] = VCM.Mag() / VC[3];
772 
773  Double_t T3MAX = 0.;
774  if (TMath::AreEqualAbs(K[3], 1., 1.e-16))
775  T3MAX = TMath::PiOver2();
776  else if (K[3] < 1.)
777  T3MAX = TMath::Pi();
778  else
779  T3MAX = TMath::ATan((1. / TMath::Sqrt(K[3] * K[3] - 1.)) / GetCMGamma());
780  TETAMAX[3] = T3MAX * TMath::RadToDeg();
781 
782  if (!GetNucleus(target_like))
783  return; //only valid for binary channels
784  //total cm energy of nucleus 4 (quasitarg)
785  WC[4] = WCT - GetEDiss() - WC[3];
786  VC[4] = GetVelocity(AM4, WC[4]);
787  //cm kinetic energy
788  EC[4] = WC[4] - AM4;
789  if (VC[4] > 0.) K[4] = VCM.Mag() / VC[4];
790 
791  Double_t T4MAX = 0.;
792  if (TMath::AreEqualAbs(K[4], 1., 1.e-16))
793  T4MAX = TMath::PiOver2();
794  else if (K[4] < 1.)
795  T4MAX = TMath::Pi();
796  else
797  T4MAX = TMath::ATan((1. / TMath::Sqrt(K[4] * K[4] - 1.)) / GetCMGamma());
798  if (TMath::Abs(GetQReaction()) < 1.E-08 && K[4] < 1.)
799  T4MAX = TMath::PiOver2();
800  TETAMAX[4] = T4MAX * TMath::RadToDeg();
801 }
802 
803 
804 
814 
816 {
817  // \return symbolic representation of the reaction as \f$A(X,Y)B\f$
818  // where
819  // + \f$A\f$ is the target nucleus;
820  // + \f$X\f$ is the projectile nucleus;
821  // + \f$Y\f$ is the outgoing projectile-like nucleus;
822  // + \f$B\f$ is the outgoing target-like nucleus.
823  //
824  // If dissipated/excitation energy is non-zero, we add \f$E_{X}=...\f$ at the end.
825 
826  TString eqn;
828  eqn.Form("%s(%s,%s)%s",
829  GetNucleus(target)->GetLatexSymbol(), GetNucleus(projectile)->GetLatexSymbol(),
830  GetNucleus(projectile_like)->GetLatexSymbol(), GetNucleus(target_like)->GetLatexSymbol());
831  else if (GetNucleus(projectile_like) && !GetNucleus(target_like)) // fusion
832  eqn.Form("%s(%s,%s)",
833  GetNucleus(target)->GetLatexSymbol(), GetNucleus(projectile)->GetLatexSymbol(),
834  GetNucleus(projectile_like)->GetLatexSymbol());
835  else if (!GetNucleus(target) && GetNucleus(target_like)) // fission
836  eqn.Form("(%s,%s)%s",
837  GetNucleus(projectile)->GetLatexSymbol(),
838  GetNucleus(projectile_like)->GetLatexSymbol(), GetNucleus(target_like)->GetLatexSymbol());
839  else if (!GetNucleus(projectile_like) && !GetNucleus(target_like)) // elastic scattering
840  eqn.Form("%s(%s,%s)%s",
841  GetNucleus(target)->GetLatexSymbol(), GetNucleus(projectile)->GetLatexSymbol(),
842  GetNucleus(projectile)->GetLatexSymbol(), GetNucleus(target)->GetLatexSymbol());
843  if (TMath::Abs(GetEDiss()) > 0)
844  eqn.Append(Form(" E_{X}=%.1f MeV", GetEDiss()));
845  return eqn;
846 }
847 
848 
849 
850 
855 
857 {
858  //Returns kinetic energy of nucleus in the centre of mass of the reaction
859  // \param[in] i nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
860  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
861  return EC[i];
862 }
863 
864 
865 
866 
873 
874 void KV2Body::Print(Option_t* opt) const
875 {
876  // Print out characteristics of reaction.
877  //
878  // If a two-body exit channel has been defined, you can use the following options:
879  // opt = "ruth" : list Rutherford scattering cross-sections as a function of angle in laboratory frame
880  // opt = "lab" : list energies and angles in laboratory frame
881 
882  cout << " ***** REACTION " << GetNucleus(projectile)->GetSymbol();
883  if (GetNucleus(target))
884  cout << " + " << GetNucleus(target)->GetSymbol();
886  cout << " ---> " << GetNucleus(projectile_like)->GetSymbol();
887  }
888  if (GetNucleus(target_like))
889  cout << " + " << GetNucleus(target_like)->GetSymbol();
890  cout << " ******" << endl;
891 
892  cout << " E.LAB = " << GetNucleus(projectile)->GetEnergy() << " MEV";
894  cout << " QGG = " << GetQGroundStates() << " MEV";
895  }
896  cout << endl;
897  cout << " E.EXC = " << GetExcitEnergy() << " MEV";
899  cout << " ==> Q-REACTION = " << GetQReaction() << " MEV";
900  }
901  cout << endl;
902  cout << " AVAILABLE ENERGY IN C.M. : ECM = " << GetCMEnergy() <<
903  " MEV (" << GetCMEnergy() / (GetNucleus(projectile)->GetA() +
905  GetA() : 0.)) << " MEV/A)" << endl;
906  cout << " PROJECTILE VELOCITY IN LAB " << GetNucleus(projectile)->GetV().Mag()
907  << " CM/NS ( " << GetNucleus(projectile)->Beta() << " * C )" << endl;
908  cout << " VELOCITY OF C.M. " << GetCMVelocity().
909  Mag() << " CM/NS" << endl;
910 
911  for (int i=1; i<=4; ++i) {
913  cout << " ENERGY - VELOCITY OF NUCLEUS " << i << " IN CM : " <<
914  GetCMEnergy((nucleus_of_interest)i) << " MEV " << GetCMVelocity((nucleus_of_interest)i) << " CM/NS (K=" << K[i] << ")" <<
915  endl;
916  }
918  cout << " MAXIMUM SCATTERING ANGLE IN LABORATORY" << endl;
919  cout << " THETA #projectile_like# " << GetMaxAngleLab(projectile_like) << " DEG." <<
920  endl;
921  if (GetNucleus(target_like))
922  cout << " THETA #target_like# " << GetMaxAngleLab(target_like) << " DEG." <<
923  endl;
924  }
925  if (GetNucleus(target)) {
926  cout << " GRAZING ANGLE IN LABORATORY : PROJECTILE " <<
927  GetLabGrazingAngle(projectile) << " DEG." << endl;
928  cout << " GRAZING ANGLE IN LABORATORY : TARGET " <<
929  GetLabGrazingAngle(target) << " DEG." << endl;
930  }
931 
933  Int_t nsteps = 15;
934 
935  if (!strcmp(opt, "lab") || !strcmp(opt, "LAB")
936  || !strcmp(opt, "Lab")) {
937  //laboratory angular distribution
938  cout << endl <<
939  " LABORATORY ANGULAR DISTRIBUTION" << endl
940  << endl;
941  cout <<
942  "Th-PLF [deg.] E-PLF [MeV] V-PLF [cm/ns] Th-TLF [deg.] E-TLF [MeV] Th-PLF [deg.]"
943  << endl;
944  cout <<
945  " (LAB) (LAB) (LAB) (LAB) (LAB) (C.M)"
946  << endl;
947  Double_t dtheta = GetMaxAngleLab(projectile_like) / nsteps;
948  for (int step = 0; step <= nsteps; step++) {
949  Double_t theta = dtheta * step;
950  auto elabP = GetELab(projectile_like, theta);
951  if(!elabP)
952  continue;
953  auto tcmP = GetThetaCM(theta, projectile_like);
954  auto vP = GetVLab(projectile_like, theta);
955  auto tlT = GetThetaLab(target_like, theta, projectile_like);
956  auto elabT = GetELab(target_like, theta, projectile_like);
957  if(elabP->second())
958  printf
959  (" %6.2f %7.2f /%7.2f %5.2f /%5.2f %6.2f /%6.2f %7.2f /%7.2f %6.2f/%6.2f\n",
960  theta, *elabP->first(), *elabP->second(), *vP->first(), *vP->second(),
961  *tlT->first(), *tlT->second(), *elabT->first(), *elabT->second(), *tcmP->first(), *tcmP->second());
962  else
963  printf
964  (" %6.2f %7.2f %5.2f %6.2f %7.2f %6.2f\n",
965  theta, *elabP->first(), *vP->first(), *tlT->first(), *elabT->first(), *tcmP->first());
966  }
967  }
968  if (GetNucleus(target) && (!strcmp(opt, "ruth") || !strcmp(opt, "RUTH")
969  || !strcmp(opt, "Ruth"))) {
970  //laboratory angular distribution with Rutherford x-section
971  cout << endl <<
972  " RUTHERFORD LABORATORY ANGULAR DISTRIBUTION" <<
973  endl << endl;
974  cout <<
975  " TETA 3 TETA 3 EN.DIF.3 SIG.RUT. SIG.RUT."
976  << endl;
977  cout <<
978  " (LAB) (CM) (LAB) (LAB) (b/sr) (CM) (b/sr)"
979  << endl;
980  //go up to grazing angle
981  Double_t dtheta = GetLabGrazingAngle(projectile) / nsteps;
982  for (int step = 0; step < nsteps; step++) {
983  Double_t theta = dtheta * (step + 1);
984  auto elabP = GetELab(projectile_like, theta);
985  if(!elabP)
986  continue;
987  auto tcm = GetThetaCM(theta, projectile_like);
988  printf
989  (" %6.2f %6.2f %7.2f %10.4g %10.4g\n",
990  theta, *tcm->first(), *elabP->first(),
991  GetXSecRuthLab(theta), GetXSecRuthCM(theta));
992  }
993  }
994  }
995 }
996 
997 
998 
1003 
1004 KV2Body::KinematicSolutions_v KV2Body::GetELab(nucleus_of_interest OfNucleus, Double_t ThetaLab, std::optional<nucleus_of_interest> AngleNucleus) const
1005 {
1006  // \returns laboratory kinetic energy(ies) of nucleus OfNucleus as a function of the lab angle ThetaLab
1007  // \param[in] AngleNucleus if given, ThetaLab is AngleNucleus' lab angle
1008 
1009  // calculate CM angle(s)
1010  if(!AngleNucleus)
1011  AngleNucleus = OfNucleus;
1012  auto AN_theta_cm = GetThetaCM(ThetaLab, *AngleNucleus);
1013  if(AN_theta_cm)
1014  {
1015  KinematicSolutions ON_theta_cm;
1016  int i=0;
1017  AN_theta_cm->iterate_solutions([&](double x){ON_theta_cm.solutions[i++] = GetThetaCM(OfNucleus, x, *AngleNucleus);});
1018  // calculate lab energie(s) corresponding to CM angle(s)
1019  KinematicSolutions energies;
1020  i=0;
1021  ON_theta_cm.iterate_solutions([&](double x){energies.solutions[i++] = GetELab(x, OfNucleus);});
1022  return energies;
1023  }
1024  return {};
1025 }
1026 
1027 
1028 
1044 
1046 {
1047  // \returns CM angle(s) of nucleus OfNucleus for a given lab angle, if solution(s) exist
1048  //
1049  // \note A one-to-one correspondence between \f$\theta_{cm}\f$ and \f$\theta_{lab}\f$
1050  // only exists for:
1051  // + in direct kinematics, the projectile-like nucleus for elastic or inelastic scattering;
1052  // + in symmetric kinematics, both nuclei (only for elastic scattering);
1053  // + in any kinematics, the target-like nucleus for elastic scattering.
1054  //
1055  // In all other cases, there are two \f$\theta_{cm}\f$ values for each \f$\theta_{lab}\f$,
1056  // up to the maximum lab scattering angle which corresponds to a single \f$\theta_{cm}\f$.
1057  //
1058  // If ThetaLab>max lab scattering angle for nucleus, no solution exists.
1059  //
1060  // In all cases, the solution with the smallest \f$\theta_{cm}\f$ corresponds to the largest \f$E_{lab}\f$.
1061 
1062  auto fonc = GetThetaLabVsThetaCMFunc(OfNucleus);
1063  return FindRoots(fonc,ThetaLab);
1064 }
1065 
1066 
1067 
1071 
1072 std::optional<double> KV2Body::GetThetaCM(Double_t ThetaLab, nucleus_of_interest OfNucleus, kinematic_solution kine_sol) const
1073 {
1074  // \return \f$\theta_{cm}\f$ corresponding to \f$\theta_{lab}\f$ for given nucleus of interest and kinematic solution,
1075  // if defined.
1076 
1077  auto kine = GetThetaCM(ThetaLab, OfNucleus);
1078  if(kine)
1079  {
1080  if(kine_sol == kinematic_solution::low_E_branch)
1081  return kine->second();
1082  return kine->first();
1083  }
1084  return {};
1085 }
1086 
1087 
1088 
1093 
1095 {
1096  // Function calculating lab energy of nucleus par[0] for any lab angle x[0]
1097  //
1098  // par[1] is choice of kinematic solution
1099 
1100  auto elab = GetELab((nucleus_of_interest)par[0], x[0], (nucleus_of_interest)par[0]);
1101  if(elab && elab->second() && ((kinematic_solution)par[1])==kinematic_solution::low_E_branch)
1102  return *elab->second();
1103  return elab ? *elab->first() : 0;
1104 }
1105 
1106 
1107 
1110 
1112 {
1113  // Function calculating lab energy of nucleus par[0] for any CM angle x[0]
1114 
1115  nucleus_of_interest OfNucleus = (nucleus_of_interest)par[0];
1116  Double_t TCM = x[0] * TMath::DegToRad();
1117  Double_t WL =
1118  WC[OfNucleus] * GetCMGamma() * (1. +
1119  BCM * (VC[OfNucleus] / KVParticle::C()) *
1120  TMath::Cos(TCM));
1121  return (WL - GetNucleus(OfNucleus)->GetMass());
1122 }
1123 
1124 
1125 
1129 
1131 {
1132  // Calculate Lab angle of nucleus as function of CM angle x[0]
1133  // par[0] = index of nucleus = 1, 2, 3, 4
1134 
1135  Int_t OfNucleus = (Int_t)par[0];
1136  Double_t sin_theta, cos_theta;
1137  if(x[0]==0)
1138  {
1139  return 0;
1140  }
1141  else if(x[0]==180)
1142  {
1143  if(TMath::AreEqualAbs(K[OfNucleus],1.0,1.e-10)){
1144  return 90.;
1145  }
1146  else if(K[OfNucleus]<1){
1147  return 180.;
1148  }
1149  return 0.;
1150  }
1151 
1152  auto ThetaCM = x[0]*TMath::DegToRad();
1153  Double_t ThetaL = TMath::ATan2(TMath::Sin(ThetaCM), (K[OfNucleus] + TMath::Cos(ThetaCM)) * GetCMGamma()) * TMath::RadToDeg();
1154 
1155  if (ThetaL < 0.) ThetaL += 180.;
1156  return ThetaL;
1157 }
1158 
1159 
1160 
1165 
1167 {
1168  // Return TF1 giving lab energy of nucleus as function of CM angle
1169  // \param[in] OfNucleus nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
1170  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
1171 
1172  if (!fELabVsThetaCM[OfNucleus]) {
1173  fELabVsThetaCM[OfNucleus] = new TF1(Form("KV2Body:ELabVsThetaCM:%d", OfNucleus),
1174  const_cast<KV2Body*>(this), &KV2Body::ELabVsThetaCM, 0, 180, 1, "KV2Body", "ELabVsThetaCM");
1175  fELabVsThetaCM[OfNucleus]->SetNpx(1000);
1176  fELabVsThetaCM[OfNucleus]->SetParameter(0, OfNucleus);
1177  }
1178  return fELabVsThetaCM[OfNucleus];
1179 }
1180 
1181 
1182 
1197 
1199 {
1200  // Return TF1 giving lab energy of nucleus as function of its lab angle
1201  // \param[in] OfNucleus nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
1202  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
1203  // \param[in] kine_sol in case of multiple solutions, choose either high- or low-energy branch
1204  //
1205  // \note A one-to-one correspondence between \f$\theta_{cm}\f$ and \f$\theta_{lab}\f$
1206  // only exists for:
1207  // + in direct kinematics, the projectile-like nucleus for elastic or inelastic scattering;
1208  // + in symmetric kinematics, both nuclei (only for elastic scattering);
1209  // + in any kinematics, the target-like nucleus for elastic scattering.
1210  // In all other cases, there are two \f$\theta_{cm}\f$ values for each \f$\theta_{lab}\f$,
1211  // up to the maximum lab scattering angle which corresponds to a single \f$\theta_{cm}\f$.
1212  // In all cases, the solution with the smallest \f$\theta_{cm}\f$ corresponds to the largest \f$E_{lab}\f$.
1213  if (!fELabVsThetaLab[OfNucleus]) {
1214  fELabVsThetaLab[OfNucleus] = new TF1(Form("KV2Body:ELabVsThetaLab:%d", OfNucleus),
1215  const_cast<KV2Body*>(this), &KV2Body::ELabVsThetaLab, 0, 180, 2, "KV2Body", "ELabVsThetaLab");
1216  fELabVsThetaLab[OfNucleus]->SetNpx(1000);
1217  fELabVsThetaLab[OfNucleus]->SetParameter(0, OfNucleus);
1218  }
1219  fELabVsThetaLab[OfNucleus]->SetParameter(1, kine_sol);
1220  return fELabVsThetaLab[OfNucleus];
1221 }
1222 
1223 
1224 
1225 
1229 
1230 KV2Body::KinematicSolutions_v KV2Body::GetVLab(nucleus_of_interest OfNucleus, Double_t ThetaLab, std::optional<nucleus_of_interest> AngleNucleus) const
1231 {
1232  // \returns vector of laboratory velocity(ies) of nucleus OfNucleus as a function of the lab angle ThetaLab
1233  // \param[in] AngleNucleus if given, ThetaLab is AngleNucleus' lab angle
1234 
1235  auto elab = GetELab(OfNucleus, ThetaLab, AngleNucleus);
1236  if(!elab) // no solutions
1237  return {};
1238 
1239  KinematicSolutions velo;
1240  int i=0;
1241  elab->iterate_solutions([&](double e){ velo.solutions[i++] = GetVelocity(GetNucleus(OfNucleus)->GetMass(), e + GetNucleus(OfNucleus)->GetMass());} );
1242  return velo;
1243 }
1244 
1245 
1246 
1247 
1250 
1252 {
1253  // \returns laboratory angle(s) of nucleus OfNucleus as a function of the lab angle ThetaLab corresponding to nucleus AngleNucleus
1254 
1255  if (ThetaLab > TETAMAX[AngleNucleus]) return {};
1256  if (!(TMath::Abs(OfNucleus - AngleNucleus) % 2)) {
1257  // same nucleus!
1258  return {{ThetaLab,{}}};
1259  }
1260  // calculate CM angle(s)
1261  auto TCM = GetThetaCM(ThetaLab, AngleNucleus);
1262  if(!TCM) // no solutions
1263  return {};
1264 
1265  KinematicSolutions tcm;
1266  int i = 0;
1267  TCM->iterate_solutions([&](double T){ tcm.solutions[i++] = GetThetaCM(OfNucleus, T, AngleNucleus); });
1268  // calculate lab angle(s) of corresponding CM angle(s)
1269  i = 0;
1270  KinematicSolutions tlab;
1271  tcm.iterate_solutions([&](double t){ tlab.solutions[i++] = GetThetaLab(t, OfNucleus); });
1272  return tlab;
1273 }
1274 
1275 
1276 
1281 
1283 {
1284  // \return TF1 giving lab angle of nucleus as function of CM angle
1285  // \param[in] OfNucleus nucleus_of_interest::projectile or nucleus_of_interest::target (entrance channel),
1286  // or nucleus_of_interest::projectile_like or nucleus_of_interest::target_like (exit channel)
1287 
1288  if (!fThetaLabVsThetaCM[OfNucleus]) {
1289  fThetaLabVsThetaCM[OfNucleus] = new TF1(Form("KV2Body:ThetaLabVsThetaCM:%d", OfNucleus),
1290  const_cast<KV2Body*>(this), &KV2Body::ThetaLabVsThetaCM, 0, 180, 1, "KV2Body", "ThetaLabVsThetaCM");
1291  fThetaLabVsThetaCM[OfNucleus]->SetNpx(1000);
1292  fThetaLabVsThetaCM[OfNucleus]->SetParameter(0, OfNucleus);
1293  }
1294  return fThetaLabVsThetaCM[OfNucleus];
1295 }
1296 
1297 
1298 
1318 
1320 {
1321  // Fill and return TGraph with complete correlation between lab energy [MeV] and angle [deg]
1322  //
1323  // \param ofNuc index of outgoing nucleus for which to draw correlation:
1324  // + ofNuc=3: projectile-like
1325  // + ofNuc=4: target-like
1326  // \param npoints number of points in graph [default=50]
1327  //
1328  // \returns a KVDrawable<TGraph> object which can either be used directly as a TGraph with the `->` operator:
1329  //~~~{.cpp}
1330  //KV2Body toto;
1331  //toto.GraphELabVsThetaLab()->Draw("al");
1332  //~~~
1333  // or with a chain of instructions to define the colour, style, etc. of the lines, points:
1334  //~~~{.cpp}
1335  //toto.GraphELabVsThetaLab().LineColor(kRed).LineWidth(2).LineStyle(kDashed)->Draw("al");
1336  //~~~
1337  //
1338  // See \ref Example_KV2Body.C for a full example of use.
1339 
1340  auto gr = new TGraph;
1341  gr->SetTitle(Form("%s: %s", GetNucleus(OfNucleus)->GetLatexSymbol(), GetReactionEquation().Data()));
1342  double dth = 180. / (npoints - 1.);
1343  int i = 0;
1344  while (i < npoints) {
1345  auto thcm = i * dth;
1346  auto thlab = GetThetaLabVsThetaCMFunc(OfNucleus)->Eval(thcm);
1347  auto elab = GetELabVsThetaCMFunc(OfNucleus)->Eval(thcm);
1348  gr->AddPoint(thlab, elab);
1349  ++i;
1350  }
1351  return gr;
1352 }
1353 
1354 
1355 
1359 
1361 {
1362  // Calculate CM Rutherford cross-section (b/sr) in the CM as a function of
1363  // scattering angle in the CM frame for nucleus par[0]
1364 
1365  Double_t PB =
1366  1.44 * GetNucleus(projectile)->GetZ() * GetNucleus(target)->GetZ() /
1367  GetCMEnergy();
1368  Double_t TCM = x[0] * TMath::DegToRad();
1369  auto OfNucleus = (nucleus_of_interest)par[0];
1370  if (OfNucleus == target || OfNucleus == target_like) TCM = TMath::Pi() - TCM; // get projectile CM angle from target CM angle
1371  Double_t D = 1. / (16. * (TMath::Power(TMath::Sin(TCM / 2.), 4.)));
1372  return ((TMath::Power(PB, 2.)) * D / 100.);
1373 }
1374 
1375 
1376 
1381 
1383 {
1384  // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the CM as a
1385  // function of projectile_like or target_like CM scattering angle
1386  // By default, theta_min = 1 degree & theta_max = 179 degrees
1387 
1388  TString name = "RuthXSecCM: ";
1390  name += " + ";
1392  name += " ";
1394  name += Form("%6.1f AMeV ", elab);
1395  if (OfNucleus == projectile_like) name += "(projectile)";
1396  else name += "(target)";
1397  TF1* fXSecRuthCM = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
1398  if (!fXSecRuthCM) {
1399  fXSecRuthCM = new TF1(name.Data(),
1400  const_cast<KV2Body*>(this), &KV2Body::XSecRuthCMVsThetaCM,
1401  theta_cm_min, theta_cm_max, 1, "KV2Body", "XSecRuthCMVsThetaCM");
1402  fXSecRuthCM->SetParameter(0, OfNucleus);
1403  fXSecRuthCM->SetNpx(1000);
1404  }
1405  fXSecRuthCM->SetRange(theta_cm_min, theta_cm_max); //in case TF1 already exists, but new range is required
1406  return fXSecRuthCM;
1407 }
1408 
1409 
1410 
1421 
1423 {
1424  //Calculate Rutherford cross-section (b/sr) in the CM as a function of projectile_like or target_like lab scattering angle
1425  //
1426  // x[0]: lab angle
1427  // par[0]: nucleus_of_interest
1428  // par[1]: kinematic_solution
1429  //
1430  // If no kinematic solution exists (angle > max lab angle), return 0
1431  //
1432  // If low energy branch is requested but does not exist, return 0
1433 
1434  auto tcm = GetThetaCM(x[0], (nucleus_of_interest)par[0]);
1435  if(tcm)
1436  {
1437  if((kinematic_solution)par[1]==low_E_branch)
1438  {
1439  if(tcm->second())
1440  {
1441  auto theta = *tcm->second();
1442  return XSecRuthCMVsThetaCM(&theta, par);
1443  }
1444  return 0;
1445  }
1446  else
1447  {
1448  auto theta = *tcm->first();
1449  return XSecRuthCMVsThetaCM(&theta, par);
1450  }
1451  }
1452  return 0;
1453 }
1454 
1455 
1456 
1462 
1464 {
1465  //Calculate Rutherford cross-section (b/sr) in the CM as a
1466  //function of projectile_like or target_like lab scattering angle
1467  //
1468  // For multi-valued kinematics, we return by default the high energy branch
1469 
1470  std::array<double,2> par = { (double)OfNucleus, (double)kineSol };
1471  return const_cast<KV2Body*>(this)->XSecRuthCM(&ThetaLab, par.data());
1472 }
1473 
1474 
1475 
1487 
1489 {
1490  //Calculate Rutherford cross-section (b/sr) in the Lab as a
1491  //function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
1492  //
1493  // x[0]: lab angle
1494  // par[0]: nucleus_of_interest
1495  // par[1]: kinematic_solution
1496  //
1497  // Checked with LISE++ calculations: Jacobian is correct for both branches
1498  // (1+K*Cos(Theta_cm) can become negative for low E branch (Theta_cm>90),
1499  // correct cross-section is absolute value of DSIDTB*RLC)
1500 
1501  Double_t DSIDTB = XSecRuthCM(x, par);
1502  auto T3CM_v = GetThetaCM(x[0], (nucleus_of_interest)par[0], (kinematic_solution)par[1]);
1503  if(!T3CM_v)
1504  return 0;
1505  auto T3CM = *T3CM_v;
1506 
1507  Double_t T3L = GetThetaLab(T3CM, (nucleus_of_interest)par[0]) * TMath::DegToRad();
1508  T3CM *= TMath::DegToRad();
1509  Double_t RLC = (TMath::Power(TMath::Sin(T3CM), 3.)) /
1510  ((TMath::Power(TMath::Sin(T3L), 3.)) * GetCMGamma() *
1511  (1. + K[(int)par[0]] * TMath::Cos(T3CM)));
1512  DSIDTB *= RLC;
1513  return TMath::Abs(DSIDTB);
1514 }
1515 
1516 
1517 
1523 
1525 {
1526  //Calculate Rutherford cross-section (b/sr) in the Lab as a
1527  //function of projectile_like or target_like lab scattering angle
1528  //
1529  // For multi-valued kinematics, we return by default the high energy branch
1530 
1531  std::vector<double> par = { (double)OfNucleus, (double)kineSol };
1532  return const_cast<KV2Body*>(this)->XSecRuthLab(&ThetaLab, par.data());
1533 }
1534 
1535 
1536 
1545 
1547 {
1548  //Rutherford cross-section (b/sr) function in the Lab as a
1549  //function of projectile_like or target_like lab scattering angle x[0]
1550  //including 'sin theta' factor needed for integrating over solid angles.
1551  //
1552  // x[0]: lab angle
1553  // par[0]: nucleus_of_interest
1554  // par[1]: kinematic_solution
1555 
1556  return XSecRuthLab(x, par) * TMath::Sin(x[0] * TMath::DegToRad());
1557 }
1558 
1559 
1560 //Double_t KV2Body::GetIntegratedXSecRuthLab(KVTelescope* tel, Int_t OfNucleus)
1561 //{
1562 // //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1563 // //polar and azimuthal angular range of the given KVTelescope.
1564 // // if (OfNucleus==3) => X-section for scattered projectile
1565 // // if (OfNucleus==4) => X-section for scattered target
1566 // //
1567 // //The returned value is in barns
1568 //
1569 // return GetIntegratedXSecRuthLab(tel->GetThetaMin(), tel->GetThetaMax(), tel->GetPhiMin(), tel->GetPhiMax(), OfNucleus);
1570 //}
1571 //
1573 //Double_t KV2Body::GetIntegratedXSecRuthLab(KVDetector* det, Int_t OfNucleus)
1574 //{
1575 // //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1576 // //polar and azimuthal angular range of the given KVDetector. These will be taken
1577 // //from the parent KVTelescope of the detector.
1578 // // if (OfNucleus==3) => X-section for scattered projectile
1579 // // if (OfNucleus==4) => X-section for scattered target
1580 // //
1581 // //The returned value is in barns
1582 //
1583 // KVTelescope* tel = (KVTelescope*)det->GetParentStructure("TELESCOPE");
1584 // if (!det) {
1585 // Error("GetIntegratedXSecRuthLab(KVDetector*,Int_t)",
1586 // "Detector has no parent telescope: it has not been positioned in a multidetector geometry");
1587 // return 0;
1588 // }
1589 // return GetIntegratedXSecRuthLab(tel, OfNucleus);
1590 //}
1591 
1592 
1603 
1605 {
1606  //Calculate Integrated Rutherford cross-section (barns) in the Lab using
1607  //polar and azimuthal angular range expressed in degree
1608  // if (OfNucleus==3) This angular range is considered to be the scattered projectile one
1609  // if (OfNucleus==4) This angular range is considered to be the scattered target one
1610  //
1611  //If phi1 ou phi2 ==-1 the azimuthal width is set to 2pi
1612  //Else if phi1=phi2 the azimuthal width is set to 1 ie the integral is only on theta
1613  //
1614  //The returned value is in barns
1615 
1616  if (th2 < th1) return 0;
1617  Double_t dphi = 0;
1618  //azimuthal width expressed in rad
1619  if (phi1 == -1 || phi2 == -1) dphi = 2 * TMath::Pi();
1620  else if (phi1 == phi2) dphi = 1;
1621  else {
1622  dphi = phi2 - phi1;
1623  dphi *= TMath::DegToRad();
1624  }
1625  Double_t theta_min = 1.;
1626  Double_t theta_max = 179.;
1627  if (th1 < theta_min) theta_min = th1;
1628  if (th2 > theta_max) theta_max = th2;
1629 
1630 #if ROOT_VERSION_CODE > ROOT_VERSION(5,99,01)
1631  return GetXSecRuthLabIntegralFunc(OfNucleus, theta_min, theta_max)->Integral(th1, th2, fIntPrec) * TMath::DegToRad() * dphi;
1632 #else
1633  const Double_t* para = 0;
1634  return GetXSecRuthLabIntegralFunc(OfNucleus, theta_min, theta_max)->Integral(th1, th2, para, fIntPrec) * TMath::DegToRad() * dphi;
1635 #endif
1636 }
1637 
1638 
1639 
1640 
1641 
1646 
1648 {
1649  // calculate Bass interaction barrier B_int
1650  //
1651  // r0 = 1.07 fm
1652 
1653  const Double_t r0 = 1.07;
1654  const Double_t e2 = 1.44;
1655  Double_t A1third = pow(GetNucleus(projectile)->GetA(), 1. / 3.);
1656  Double_t A2third = pow(GetNucleus(target)->GetA(), 1. / 3.);
1657  Double_t R12 = r0 * (A1third + A2third);
1658 
1659  Double_t Bint = GetNucleus(projectile)->GetZ() * GetNucleus(target)->GetZ() * e2 / (R12 + 2.7)
1660  - 2.9 * A1third * A2third / (A1third + A2third);
1661 
1662  return Bint;
1663 }
1664 
1665 
1666 
1667 
1673 
1675 {
1676  // calculate Kox reaction X-section (in barns) for a given lab energy of projectile (in MeV/nucleon)
1677  //
1678  // c parameter fitted with Landau function vs. Log10(E/A) (see Fig. 12 of PRC Kox 87)
1679  // by imposing c=0.1 at E/A=10 MeV (goes to 0 as E/A->0)
1680 
1681  const Double_t r0 = 1.05;
1682  const Double_t rc = 1.3;
1683  const Double_t e2 = 1.44;
1684  const Double_t a = 1.9;
1685 
1686  KVNucleus* proj = GetNucleus(projectile);
1687  KVNucleus* targ = GetNucleus(target);
1688  proj->SetEnergy(eproj[0]*proj->GetA());
1690  Double_t ECM = GetCMEnergy();
1691 
1692  Double_t c = 11.3315 * TMath::Landau(TMath::Log10(eproj[0]), 2.47498, 0.554937);
1693  Double_t A1third = pow(proj->GetA(), 1. / 3.);
1694  Double_t A2third = pow(targ->GetA(), 1. / 3.);
1695  Double_t BC = proj->GetZ() * targ->GetZ() * e2 / (rc * (A1third + A2third));
1696  Double_t EFac = 1 - BC / ECM;
1697 
1698  Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1699  pow((A1third + A2third + a * A1third * A2third / (A1third + A2third) - c), 2) *
1700  EFac;
1701 
1702  return Xsec / 100.;
1703 }
1704 
1705 
1706 
1707 
1711 
1713 {
1714  // calculate Reaction Cross Section with the "Sphere Dure"
1715  // approximation
1716 
1717  Double_t A1third = pow(GetNucleus(projectile)->GetA(), 1. / 3.);
1718  Double_t A2third = pow(GetNucleus(target)->GetA(), 1. / 3.);
1719 
1720  Double_t Xsec = TMath::Pi() * pow(r0, 2) *
1721  pow(A1third + A2third, 2);
1722 
1723  return Xsec / 100.;
1724 }
1725 
1726 
1727 
1728 
1732 
1734 {
1735 
1736  // Deduce the maximum impact parameter (in fm) from a given reaction cross section (in barn)
1737  // in the approximation Xsec = \int(0,bmax) 2Pi b db
1738 
1739  return 10 * TMath::Sqrt(ReacXsec / TMath::Pi());
1740 
1741 }
1742 
1743 
1744 
1745 
1746 
1747 
1751 
1753 {
1754 
1755  // Integrate the cross section between two impact parameter (in fm)
1756  // and give the result in barn
1757 
1758  return TMath::Pi() * (TMath::Power(b2, 2.) - TMath::Power(b1, 2.)) / 100;
1759 
1760 }
1761 
1762 
1763 
1764 
1770 
1772 {
1773  // Return pointer to TF1 with Kox reaction X-section in barns as a
1774  // function of projectile lab energy (in Mev/nucleon) for this reaction.
1775  // By default the range of the function is [20,100] MeV/nucleon.
1776  // Change with TF1::SetRange.
1777 
1778  if (!fKoxReactionXSec) {
1780  name += " + ";
1782  fKoxReactionXSec = new TF1(name.Data(),
1783  const_cast<KV2Body*>(this), &KV2Body::KoxReactionXSec, 20, 100, 0, "KV2Body", "KoxReactionXSec");
1784  fKoxReactionXSec->SetNpx(1000);
1785  }
1786  return fKoxReactionXSec;
1787 }
1788 
1789 
1790 
1791 
1818 
1820 {
1821  // Calculate the mean charge state of the projectile after passage through the target,
1822  // assuming that the equilibrium charge state distribution is achieved*.
1823  // t[0] = energy of projectile after the target (in MeV/nucleon)
1824  //
1825  // We use the empirical parameterization of Leon et al., At. Dat. and Nucl. Dat. Tab. 69, 217 (1998)
1826  // developed for heavy ions in the GANIL energy range (it represents a fit to data measured using
1827  // GANIL beams).
1828  //
1829  // *N.B. Concerning equilibrium charge state distributions, it is not easy to know whether, for a given
1830  // combination of projectile, projectile energy, target, and target thickness, the equilibrium
1831  // distribution is reached or not. Here are some comments from the paper cited above which
1832  // may give some guidelines:
1833  //
1834  // "The energies available at the GANIL accelerator range from 24 to 95 MeV/u. Within this energy range,
1835  // the equilibrium charge state is reached only for fairly thick targets (~1 mg/cm2 for C foils)."
1836  //
1837  // "Mean Charge State as a Function of the Target Thickness
1838  // A typical example of the variation of the mean charge as a function of the foil thickness is shown ... It is seen
1839  // that the mean charge initially increases due to the ionization process. Then, the equilibrium state is reached at
1840  // a certain thickness, the so-called equilibrium thickness, due to the equilibration of electron loss and
1841  // capture processes. Finally, for foils thicker than the equilibrium thickness, the mean charge decreases due to the
1842  // slowing down of the ions in matter leading to higher capture cross sections."
1843  //
1844  // It should be noted that, according to the data published in this and other papers, the equilibrium thickness
1845  // decreases with increasing atomic number of the target, and increases with increasing energy of the projectile.
1846 
1847  KVNucleus* proj = GetNucleus(projectile);
1848  Double_t Zp = proj->GetZ();
1849  proj->SetEnergy(t[0]*proj->GetA());
1850  Double_t beta = proj->Beta();
1851  Double_t vp = beta * KVParticle::C();
1852  Double_t Zt = GetNucleus(target)->GetZ();
1853 
1854  Double_t q = Zp * (1. - TMath::Exp(-83.275 * beta / pow(Zp, 0.477)));
1855 
1856  // correction for target Z
1857  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);
1858  q *= g;
1859 
1860  if (Zp > 54) {
1861  // f(Zp) - correction for projectiles with Z>54
1862  Double_t f = 1. - TMath::Exp(-12.905 + 0.2124 * Zp - 0.00122 * pow(Zp, 2));
1863  q *= f;
1864  }
1865 
1866  return q;
1867 }
1868 
1869 
1870 
1871 
1881 
1883 {
1884  // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
1885  // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
1886  // energy after the target (in MeV/nucleon).
1887  // We use the empirical parameterization of Leon et al., At. Dat. and Nucl. Dat. Tab. 69, 217 (1998)
1888  // (see EqbmChargeState(Double_t *t,Double_t*) for details)
1889  //
1890  // By default the range of the function is [5,100] MeV/nucleon.
1891  // Change with TF1::SetRange.
1892 
1893  if (!fEqbmChargeState) {
1895  name += " + ";
1897  name += " LEON";
1898  fEqbmChargeState = new TF1(name.Data(),
1899  const_cast<KV2Body*>(this), &KV2Body::EqbmChargeState, 5, 100, 0, "KV2Body", "EqbmChargeState");
1900  fEqbmChargeState->SetNpx(1000);
1901  }
1902  return fEqbmChargeState;
1903 }
1904 
1905 
1906 
1910 
1912 {
1913  // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1914  // for solid targets
1915 
1916  KVNucleus* proj = GetNucleus(projectile);
1917  Double_t Zp = proj->GetZ();
1918  proj->SetEnergy(t[0]*proj->GetA());
1919  KVNucleus* targ = GetNucleus(target);
1920  Double_t Zt = targ->GetZ();
1921  Double_t Vp = proj->Beta();
1922  const Double_t V0 = 2.19e+06 / TMath::C();
1923  Double_t x = -0.019 * pow(Zp, -0.52) * Vp / V0;
1924  x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x) / 1.68;
1925  x = pow(x, 1. + 1.8 / Zp);
1926  Double_t q = 0.07 / x + 6. + 0.3 * pow(x, 0.5) + 10.37 * x + pow(x, 4);
1927  q = (Zp * (12 * x + pow(x, 4.))) / q;
1928  return q;
1929 }
1930 
1931 
1941 
1943 {
1944  // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
1945  // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
1946  // energy after the target (in MeV/nucleon).
1947  // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1948  // This formula is valid for solid targets.
1949  //
1950  // By default the range of the function is [5,100] MeV/nucleon.
1951  // Change with TF1::SetRange.
1952 
1953  if (!fEqbmChargeStateShSol) {
1955  name += " + ";
1957  name += " SHIWIETZ-SOLID";
1958  fEqbmChargeStateShSol = new TF1(name.Data(),
1959  const_cast<KV2Body*>(this), &KV2Body::eqbm_charge_state_shiwietz_solid, 5, 100, 0, "KV2Body", "eqbm_charge_state_shiwietz_solid");
1961  }
1962  return fEqbmChargeStateShSol;
1963 }
1964 
1965 
1966 
1970 
1972 {
1973  // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
1974  // for gas targets
1975 
1976  KVNucleus* proj = GetNucleus(projectile);
1977  Double_t Zp = proj->GetZ();
1978  proj->SetEnergy(t[0]*proj->GetA());
1979  KVNucleus* targ = GetNucleus(target);
1980  Double_t Zt = targ->GetZ();
1981  Double_t Vp = proj->Beta();
1982  const Double_t V0 = 2.19e+06 / TMath::C();
1983  Double_t x = 0.03 - 0.017 * pow(Zp, -0.52) * Vp / V0;
1984  x = Vp / V0 * pow(Zp, -0.52) * pow(Zt, x);
1985  x = pow(x, 1. + 0.4 / Zp);
1986  Double_t q = 1428. - 1206.*pow(x, 0.5) + 690 * x + pow(x, 6.);
1987  q = (Zp * (376.*x + pow(x, 6.))) / q;
1988  return q;
1989 }
1990 
1991 
1992 
2002 
2004 {
2005  // Return pointer to TF1 giving mean charge state of the projectile after passage through the target,
2006  // assuming that the equilibrium charge state distribution is achieved, as a function of projectile
2007  // energy after the target (in MeV/nucleon).
2008  // G. Shiwietz et al Nucl. Instr. and Meth. in Phys. Res. B 175-177 (2001) 125-131
2009  // This formula is valid for gas targets.
2010  //
2011  // By default the range of the function is [5,100] MeV/nucleon.
2012  // Change with TF1::SetRange.
2013 
2014  if (!fEqbmChargeStateShGas) {
2016  name += " + ";
2018  name += " SHIWIETZ-GAS";
2019  fEqbmChargeStateShGas = new TF1(name.Data(),
2020  const_cast<KV2Body*>(this), &KV2Body::eqbm_charge_state_shiwietz_gas, 5, 100, 0, "KV2Body", "eqbm_charge_state_shiwietz_gas");
2022  }
2023  return fEqbmChargeStateShGas;
2024 }
2025 
2026 
2027 
2028 
2034 
2036 {
2037  // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the Lab as a
2038  // function of projectile_like or target_like lab scattering angle
2039  //
2040  // By default, theta_min = 1 degree & theta_max = 179 degrees
2041 
2042  if (theta_max > GetMaxAngleLab(OfNucleus)) {
2043  theta_max = GetMaxAngleLab(OfNucleus);
2044  Info("GetXSecRuthLabFunc", "Maximum angle set to %lf degrees", theta_max);
2045  }
2046 
2047  TString name = "RuthXSec: ";
2049  name += " + ";
2051  name += " ";
2053  name += Form("%6.1f AMeV ", elab);
2054  if (OfNucleus == projectile_like) name += "(projectile)";
2055  else name += "(target)";
2056  if(kineSol == low_E_branch) name+= "(lowE)";
2057  TF1* fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
2058  if (!fXSecRuthLab) {
2059  fXSecRuthLab = new TF1(name.Data(),
2060  const_cast<KV2Body*>(this), &KV2Body::XSecRuthLab, theta_min, theta_max, 2, "KV2Body", "XSecRuthLab");
2061  fXSecRuthLab->SetParameter(0, OfNucleus);
2062  fXSecRuthLab->SetParameter(1, kineSol);
2063  fXSecRuthLab->SetNpx(1000);
2064  }
2065  fXSecRuthLab->SetRange(theta_min, theta_max); //in case TF1 already exists, but new range is required
2066  return fXSecRuthLab;
2067 }
2068 
2069 
2070 
2083 
2085 {
2086  // Return pointer to TF1 giving Rutherford cross-section (b/sr) in the Lab as a
2087  // function of projectile (OfNucleus=3) or target (OfNucleus=4) lab scattering angle
2088  //
2089  // This function is equal to sin(theta)*dsigma/domega, i.e. it is the integrand
2090  // needed for calculating total cross-sections integrated over solid angle.
2091  //
2092  // WARNING: when integrating this function using TF1::Integral, the result must
2093  // be multiplied by TMath::DegToRad(), because 'x' is in degrees rather than radians,
2094  // e.g. the integrated cross-section in barns is given by
2095  //
2096  // GetXSecRuthLabIntegralFunc(OfNucleus)->Integral(theta_min, theta_max)*TMath::DegToRad()
2097 
2098  TString name = "RuthXSecInt: ";
2100  name += " + ";
2102  name += " ";
2104  name += Form("%6.1f AMeV ", elab);
2105  if (OfNucleus == projectile_like) name += "(projectile)";
2106  else name += "(target)";
2107 
2108  TF1* fXSecRuthLab = (TF1*)gROOT->GetListOfFunctions()->FindObject(name.Data());
2109  if (!fXSecRuthLab) {
2110  fXSecRuthLab = new TF1(name.Data(),
2111  const_cast<KV2Body*>(this), &KV2Body::XSecRuthLabInt, theta_min, theta_max, 2, "KV2Body", "XSecRuthLabInt");
2112  fXSecRuthLab->SetParameter(0, OfNucleus);
2113  fXSecRuthLab->SetParameter(1, kineSol);
2114  fXSecRuthLab->SetNpx(1000);
2115  }
2116  fXSecRuthLab->SetRange(theta_min, theta_max); //in case TF1 already exists, but new range is required
2117  return fXSecRuthLab;
2118 }
2119 
2120 
2121 
2122 
2149 
2151 {
2152  // Find at most two solutions x1 and x2 between xmin and xmax for which fonc->Eval(x) = val
2153  // i.e. we use TF1::GetX for a function which may be (at most) double-valued in range (xmin, xmax)
2154  //
2155  // This is adapted to the case of the lab angle vs. CM angle of scattered particles:
2156  // - if fonc has a maximum between xmin and xmax, and if val < max, we look for x1 between
2157  // xmin and the maximum, and x2 between the maximum and xmax.
2158  // - if not (single-valued function) we look for x1 between xmin and xmax.
2159  //
2160  // If val > maximum of fonc between xmin and xmax, there is no solution
2161  //
2162  // \returns a KinematicSolutions_v which may or may not contain kinematic solutions
2163  //
2164  // Testing and using returned value:
2165  //~~~{.cpp}
2166  // KV2Body kin;
2167  // auto roots = kin.FindRoots(kin.GetThetaLabVsThetaCMFunc(KV2Body::projectile_like), 2.6);
2168  // if(roots)
2169  // {
2170  // if(roots->second())
2171  // std::cout << "Two solutions found: second one is " << *roots->second() << std::endl;
2172  // std::cout << "First solution (always exists) is " << *roots->first() << std::endl;
2173  // }
2174  // else
2175  // std::cout << "No kinematic solutions exist\n";
2176  //~~~
2177 
2178  constexpr Double_t xmin = 0.0;
2179  constexpr Double_t xmax = 180.0;
2180 
2181  // special case: val=0 i.e. lab angle zero degrees
2182  // theta_cm=0 always gives theta_lab=0
2183  // if theta_cm=180 also gives theta_lab=0, we add second solution
2184  if(val==xmin)
2185  {
2186  if(fonc->Eval(xmax)==0)
2187  return {{xmin,xmax}};
2188  return {{xmin,{}}};
2189  }
2190  Double_t max = fonc->GetMaximum(xmin, xmax);
2191  // no solution
2192  if (val > max) return {};
2193  Double_t maxX = fonc->GetMaximumX(xmin, xmax);
2194  int nRoots = 1;
2195  if (TMath::AreEqualAbs(val, max, 1.e-5)) {
2196  // value corresponds to maximum of function
2197  return {{maxX,{}}};
2198  }
2199  // 2 roots if 'fonc' is zero at 180 deg.
2200  if (fonc->Eval(xmax) == 0) nRoots = 2;
2201  Double_t xmax1 = (nRoots == 2 ? maxX : xmax);
2202  auto x1 = fonc->GetX(val, xmin, xmax1);
2203  if (nRoots == 1)
2204  return {{x1,{}}};
2205  auto x2 = fonc->GetX(val, maxX, xmax);
2206  return {{x1,x2}};
2207 }
2208 
2209 
2210 
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:178
Double_t XSecRuthLabInt(Double_t *, Double_t *)
Definition: KV2Body.cpp:1546
Double_t TETAMIN[5]
defined only for nuclei 3 et 4
Definition: KV2Body.h:193
Double_t eqbm_charge_state_shiwietz_gas(Double_t *t, Double_t *)
Definition: KV2Body.cpp:1971
Double_t VC[5]
cm velocities
Definition: KV2Body.h:189
void SetTarget(const KVNucleus &)
Set target for reaction.
Definition: KV2Body.cpp:328
TString GetReactionEquation() const
Definition: KV2Body.cpp:815
void Set4thNucleus()
Definition: KV2Body.cpp:436
TF1 * fELabVsThetaLab[5]
Definition: KV2Body.h:211
Double_t K[5]
ratio of c.m. velocity to velocity of nucleus in c.m. v_cm/v_i_cm
Definition: KV2Body.h:191
std::vector< KVNucleus > fNuclei
nuclei involved in calculation
Definition: KV2Body.h:180
Double_t WCT
total cm energy
Definition: KV2Body.h:186
KVDrawable< TF1 > GetELabVsThetaLabFunc(nucleus_of_interest OfNucleus, kinematic_solution kine_sol=kinematic_solution::high_E_branch) const
Definition: KV2Body.cpp:1198
Double_t GetSphereDureReactionXSec(Double_t r0=1.05)
Definition: KV2Body.cpp:1712
Double_t ELabVsThetaLab(Double_t *, Double_t *)
Definition: KV2Body.cpp:1094
TF1 * fThetaLabVsThetaCM[5]
Definition: KV2Body.h:209
Double_t GetMaxAngleLab(nucleus_of_interest i) const
Definition: KV2Body.cpp:546
KinematicSolutions_v FindRoots(TF1 *, Double_t) const
Definition: KV2Body.cpp:2150
KinematicSolutions_v GetThetaCM(Double_t ThetaLab, nucleus_of_interest OfNucleus) const
Definition: KV2Body.cpp:1045
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Definition: KV2Body.cpp:353
Double_t fIntPrec
Precision of the TF1::Integral method.
Definition: KV2Body.h:218
KVDrawable< TF1 > GetXSecRuthCMFunc(nucleus_of_interest OfNucleus=projectile_like, Double_t theta_cm_min=1., Double_t theta_cm_max=179.) const
Definition: KV2Body.cpp:1382
Double_t GetXSecRuthLab(Double_t ThetaLab, nucleus_of_interest OfNucleus=projectile_like, kinematic_solution kineSol=high_E_branch) const
Definition: KV2Body.cpp:1524
Double_t GetQReaction() const
Calculate Q-value for reaction, including dissipated (excitation) energy.
Definition: KV2Body.cpp:490
std::optional< KinematicSolutions > KinematicSolutions_v
Definition: KV2Body.h:259
TF1 * fEqbmChargeStateShSol
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al solid)
Definition: KV2Body.h:197
Double_t EqbmChargeState(Double_t *t, Double_t *)
Definition: KV2Body.cpp:1819
KVDrawable< TF1 > GetShiwietzEqbmChargeStateFuncForSolidTargets() const
Definition: KV2Body.cpp:1942
void SetOutgoing(const KVNucleus &proj_out)
Definition: KV2Body.cpp:411
TF1 * fKoxReactionXSec
function Kox reaction cross-section [barns] vs. E/A projectile
Definition: KV2Body.h:195
Double_t GetELab(Double_t ThetaCM, nucleus_of_interest OfNucleus) const
Definition: KV2Body.h:358
Double_t XSecRuthCM(Double_t *, Double_t *)
Definition: KV2Body.cpp:1422
Double_t XSecRuthLab(Double_t *, Double_t *)
Definition: KV2Body.cpp:1488
Double_t BCM
beta of centre of mass
Definition: KV2Body.h:184
TF1 * fEqbmChargeStateShGas
function equilibrium charge state of projectile vs. E/A projectile (Shiwietz et al gas)
Definition: KV2Body.h:198
KVDrawable< TGraph > GetGraphELabVsThetaLab(nucleus_of_interest OfNucleus=projectile_like, int npoints=50) const
Definition: KV2Body.cpp:1319
Double_t GetThetaLab(Double_t ThetaCM, nucleus_of_interest OfNucleus) const
Definition: KV2Body.h:351
KVDrawable< TF1 > GetELabVsThetaCMFunc(nucleus_of_interest OfNucleus) const
Definition: KV2Body.cpp:1166
Double_t GetEDiss() const
Definition: KV2Body.h:315
KinematicSolutions_v GetVLab(nucleus_of_interest OfNucleus, Double_t ThetaLab, std::optional< nucleus_of_interest > AngleNucleus={}) const
Definition: KV2Body.cpp:1230
static Double_t GetVelocity(Double_t mass, Double_t E)
Definition: KV2Body.cpp:394
KVDrawable< TF1 > GetEqbmChargeStateFunc() const
Definition: KV2Body.cpp:1882
KVDrawable< TF1 > GetXSecRuthLabFunc(nucleus_of_interest OfNucleus=projectile_like, kinematic_solution kineSol=high_E_branch, Double_t theta_min=1., Double_t theta_max=179.) const
Definition: KV2Body.cpp:2035
Double_t GetBmaxFromReactionXSec(Double_t ReacXsec)
Definition: KV2Body.cpp:1733
Double_t GetLabGrazingAngle(nucleus_of_interest i=projectile) const
Definition: KV2Body.cpp:613
TVector3 VCM
velocity of centre of mass
Definition: KV2Body.h:183
KVDrawable< TF1 > GetShiwietzEqbmChargeStateFuncForGasTargets() const
Definition: KV2Body.cpp:2003
Double_t KoxReactionXSec(Double_t *, Double_t *)
Definition: KV2Body.cpp:1674
nucleus_of_interest
Definition: KV2Body.h:263
@ target_like
Definition: KV2Body.h:267
@ target
Definition: KV2Body.h:265
@ projectile_like
Definition: KV2Body.h:266
@ projectile
Definition: KV2Body.h:264
Double_t ELabVsThetaCM(Double_t *, Double_t *)
Function calculating lab energy of nucleus par[0] for any CM angle x[0].
Definition: KV2Body.cpp:1111
Double_t ThetaLabVsThetaCM(Double_t *, Double_t *)
Definition: KV2Body.cpp:1130
Double_t GetCMGamma() const
Definition: KV2Body.h:334
Double_t GetMinAngleLab(nucleus_of_interest i) const
Definition: KV2Body.cpp:565
Double_t TETAMAX[5]
defined only for nuclei 3 et 4
Definition: KV2Body.h:192
Bool_t fSetOutgoing
= kTRUE if SetOutgoing is called before CalculateKinematics
Definition: KV2Body.h:216
Double_t WC[5]
cm energy of each nucleus
Definition: KV2Body.h:187
Double_t EC[5]
cm energies
Definition: KV2Body.h:190
Double_t GetExcitEnergy() const
Definition: KV2Body.h:304
Double_t WLT
total lab energy
Definition: KV2Body.h:185
TF1 * fXSecRuthLab[5]
Definition: KV2Body.h:214
Double_t BassIntBarrier()
Definition: KV2Body.cpp:1647
TVector3 GetCMVelocity() const
Return vector velocity of centre of mass of reaction (units: cm/ns)
Definition: KV2Body.cpp:583
void CalculateKinematics()
Definition: KV2Body.cpp:693
void Print(Option_t *opt="") const override
Definition: KV2Body.cpp:874
void init()
Default initialisations.
Definition: KV2Body.cpp:32
KVDrawable< TF1 > GetKoxReactionXSecFunc() const
Definition: KV2Body.cpp:1771
KV2Body()
default ctor
Definition: KV2Body.cpp:61
Double_t GetIntegratedXsec(Double_t b1, Double_t b2)
Definition: KV2Body.cpp:1752
Double_t GetQGroundStates() const
Calculate Q-value for reaction, assuming all nuclei in ground state.
Definition: KV2Body.cpp:509
Double_t XSecRuthCMVsThetaCM(Double_t *, Double_t *)
Definition: KV2Body.cpp:1360
KVDrawable< TF1 > GetThetaLabVsThetaCMFunc(nucleus_of_interest OfNucleus) const
Definition: KV2Body.cpp:1282
Double_t GetXSecRuthCM(Double_t ThetaLab, nucleus_of_interest OfNucleus=projectile_like, kinematic_solution kineSol=high_E_branch) const
Definition: KV2Body.cpp:1463
Double_t GetCMEnergy() const
Return available kinetic energy in centre of mass.
Definition: KV2Body.cpp:532
kinematic_solution
Definition: KV2Body.h:269
@ low_E_branch
Definition: KV2Body.h:271
TF1 * fEqbmChargeState
function equilibrium charge state of projectile vs. E/A projectile (Leon et al)
Definition: KV2Body.h:196
Double_t eqbm_charge_state_shiwietz_solid(Double_t *t, Double_t *)
Definition: KV2Body.cpp:1911
KVNucleus * GetNucleus(nucleus_of_interest i) const
Definition: KV2Body.cpp:468
Double_t fEDiss
dissipated energy, 0 means elastic scattering
Definition: KV2Body.h:181
TF1 * fELabVsThetaCM[5]
Definition: KV2Body.h:210
KVDrawable< TF1 > GetXSecRuthLabIntegralFunc(nucleus_of_interest OfNucleus=projectile_like, Double_t theta_min=1., Double_t theta_max=179., kinematic_solution kineSol=high_E_branch) const
Definition: KV2Body.cpp:2084
Double_t GetIntegratedXSecRuthLab(Float_t th1, Float_t th2, Float_t phi1=-1, Float_t phi2=-1, nucleus_of_interest OfNucleus=projectile_like)
Definition: KV2Body.cpp:1604
void Error(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1709
void Warning(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1696
Simple wrapper for objects which can be drawn (graphs, histograms)
Definition: KVDrawable.h:29
Utility class for kinematical transformations of KVParticle class.
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
const Char_t * GetSymbol(Option_t *opt="") const
Definition: KVNucleus.cpp:71
Double_t GetAMeV() const
Definition: KVNucleus.cpp:1378
Double_t GetExcitEnergy() const
Definition: KVNucleus.h:279
Double_t GetMassExcess(Int_t z=-1, Int_t a=-1) const
Definition: KVNucleus.cpp:880
Bool_t IsDefined() const
Definition: KVNucleus.h:198
Int_t GetA() const
Definition: KVNucleus.cpp:796
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:767
TVector3 GetMomentum() const
Definition: KVParticle.h:607
Double_t GetEnergy() const
Definition: KVParticle.h:624
static Double_t C()
Definition: KVParticle.cpp:117
Double_t GetKE() const
Definition: KVParticle.h:617
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:775
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:897
TVector3 GetV() const
Definition: KVParticle.h:674
void SetEnergy(Double_t e)
Definition: KVParticle.h:602
Double_t GetMass() const
Definition: KVParticle.h:577
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 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
virtual void AddPoint(Double_t x, Double_t y)
void SetTitle(const char *title="") override
Double_t Beta() const
Double_t E() const
virtual void SetTitle(const char *title="")
R__ALWAYS_INLINE Bool_t IsZombie() 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
TString & Append(char c, Ssiz_t rep=1)
void Form(const char *fmt,...)
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 T0 &x, const RVec< T1 > &v)
Double_t x[n]
TGraphErrors * gr
double T(double x)
double max(double x, double y)
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)
Maybe single- or double-valued.
Definition: KV2Body.h:234
std::array< std::optional< double >, 2 > solutions
Definition: KV2Body.h:235
void iterate_solutions(std::function< void(double)> f)
Definition: KV2Body.h:236
auto * th2
auto * th1
TArc a
ClassImp(TPyArg)