KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVedaLossMaterial.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Wed Feb 2 16:13:08 2011
2 //Author: frankland,,,,
3 
4 #include "KVedaLossMaterial.h"
5 #include <TMath.h>
6 #include "KVIonRangeTable.h"
7 #include <TEnv.h>
8 #include <TF1.h>
10 #include "KVedaLoss.h"
11 
13 
16 
17 
20 
22  : KVIonRangeTableMaterial(), fInvRange(ZMAX_VEDALOSS, 1),
23  fEmin(ZMAX_VEDALOSS), fEmax(ZMAX_VEDALOSS), fCoeff(ZMAX_VEDALOSS, std::vector<Double_t>(14))
24 {
25  // Default constructor
26  for (int i = 0; i < ZMAX_VEDALOSS; i++) {
27  fEmin[i] = 0.0;
28  fEmax[i] = 500.0;
29  }
31 }
32 
33 
34 
37 
38 KVedaLossMaterial::KVedaLossMaterial(const KVIonRangeTable* t, const Char_t* name, const Char_t* type, const Char_t* state,
39  Double_t density, Double_t Z, Double_t A, Double_t)
40  : KVIonRangeTableMaterial(t, name, type, state, density, Z, A), fInvRange(ZMAX_VEDALOSS, 1),
41  fEmin(ZMAX_VEDALOSS), fEmax(ZMAX_VEDALOSS), fCoeff(ZMAX_VEDALOSS, std::vector<Double_t>(14))
42 {
43  // create new material
44  fgTable = static_cast<KVedaLoss*>(const_cast<KVIonRangeTable*>(t));
45  for (int i = 0; i < ZMAX_VEDALOSS; i++) {
46  fEmin[i] = 0.0;
47  fEmax[i] = 500.0;
48  }
50 }
51 
52 
53 
56 
58 {
59  // Destructor
60 }
61 
62 
63 
79 
81 {
82  // Read Z- & A-dependent range parameters for material
83  //
84  // For each material we create 4 TF1 objects:
85  // - KVedaLossMaterial:[type]:Range - gives range in g/cm**2 as a function of particle energy
86  // - KVedaLossMaterial:[type]:StoppingPower - gives dE/dx in MeV/(g/cm**2) as a function of particle energy
87  // - KVedaLossMaterial:[type]:EnergyLoss - gives dE as a function of particle energy
88  // - KVedaLossMaterial:[type]:ResidualEnergy - gives energy after material (0 if particle stops)
89  //
90  // The TF1::fNpx parameter for these functions is defined by the environment variables
91  //
92  // - KVedaLoss.Range.Npx: 20 /* also used for StoppingPower */
93  // - KVedaLoss.EnergyLoss.Npx: 50
94  // - KVedaLoss.ResidualEnergy.Npx: 20
95  //
96 
97  char line[132];
98 
99  //look for energy limits to calculation validity
100  if (!fgets(line, 132, fp)) {
101  Warning("ReadRangeTable", "Problem reading energy limits in range table file for %s (%s)",
102  GetName(), GetType());
103  return kFALSE;
104  }
105  else {
106  if (!strncmp(line, "COMPOUND", 8)) {
107  // material is compound. read composition for TGeoMixture.
108  if (fgets(line, 132, fp)) {}
109  int nel = atoi(line); // read number of elements
110  for (int el = 0; el < nel; el++) {
111  if (fgets(line, 132, fp)) {}
112  int z, a, w;
113  sscanf(line, "%d %d %d", &z, &a, &w);
114  AddCompoundElement(z, a, w);
115  }
116  }
117  else if (!strncmp(line, "MIXTURE", 7)) {
118  // material is mixture. read composition for TGeoMixture.
119  if (fgets(line, 132, fp)) {}
120  int nel = atoi(line); // read number of elements
121  for (int el = 0; el < nel; el++) {
122  if (fgets(line, 132, fp)) {}
123  int z, a, nat;
124  float w;
125  sscanf(line, "%d %d %d %f", &z, &a, &nat, &w);
126  AddMixtureElement(z, a, nat, w);
127  }
128  }
129  }
130  if (!fgets(line, 132, fp)) {
131  Warning("ReadRangeTable", "Problem reading energy limits in range table file for %s (%s)",
132  GetName(), GetType());
133  return kFALSE;
134  }
135  else {
136  while (line[0] == 'Z') {
137  Int_t z1, z2;
138  Float_t e1, e2;
139  sscanf(line, "Z = %d,%d %f < E/A < %f MeV", &z1,
140  &z2, &e1, &e2);
141  if (fgets(line, 132, fp)) {}
142  for (int i = z1; i <= z2; i++) {
143  fEmin[i - 1] = e1;
144  fEmax[i - 1] = e2;
145  }
146  }
147  }
148 
149  // get require Npx value from (user-defined) environment variables
150  Int_t my_npx = gEnv->GetValue("KVedaLoss.Range.Npx", 100);
151 
152  fRange = new TF1(Form("KVedaLossMaterial:%s:Range", GetType()), this, &KVedaLossMaterial::RangeFunc,
153  0., 1.e+03, 0, "KVedaLossMaterial", "RangeFunc");
154  fRange->SetNpx(my_npx);
155 
156  fStopping = new TF1(Form("KVedaLossMaterial:%s:StoppingPower", GetType()), this, &KVedaLossMaterial::StoppingFunc,
157  0., 1.e+03, 0, "KVedaLossMaterial", "StoppingFunc");
158  fStopping->SetNpx(my_npx);
159 
160  my_npx = gEnv->GetValue("KVedaLoss.EnergyLoss.Npx", 100);
161  fDeltaE = new TF1(Form("KVedaLossMaterial:%s:EnergyLoss", GetType()), this, &KVedaLossMaterial::DeltaEFunc,
162  0., 1.e+03, 0, "KVedaLossMaterial", "DeltaEFunc");
163  fDeltaE->SetNpx(my_npx);
164 
165  my_npx = gEnv->GetValue("KVedaLoss.ResidualEnergy.Npx", 100);
166  fEres = new TF1(Form("KVedaLossMaterial:%s:ResidualEnergy", GetType()), this, &KVedaLossMaterial::EResFunc,
167  0., 1.e+03, 0, "KVedaLossMaterial", "EResFunc");
168  fEres->SetNpx(my_npx);
169 
170  for (int count = 0; count < ZMAX_VEDALOSS; count++) {
171 
172  if (sscanf(line, "%lf %lf %lf %lf %lf %lf %lf %lf",
173  &fCoeff[count][0], &fCoeff[count][1],
174  &fCoeff[count][2], &fCoeff[count][3],
175  &fCoeff[count][4], &fCoeff[count][5],
176  &fCoeff[count][6], &fCoeff[count][7])
177  != 8) {
178  Error("ReadRangeTable", "problem reading coeffs 0-7 in range table for %s (%s)", GetName(), GetType());
179  return kFALSE;
180  }
181  if (!fgets(line, 132, fp)) {
182  Warning("ReadRangeTable", "file too short ??? %s (%s)", GetName(), GetType());
183  return kFALSE;
184  }
185  else {
186  if (sscanf(line, "%lf %lf %lf %lf %lf %lf",
187  &fCoeff[count][8], &fCoeff[count][9],
188  &fCoeff[count][10], &fCoeff[count][11],
189  &fCoeff[count][12], &fCoeff[count][13])
190  != 6) {
191  Error("ReadRangeTable", "problem reading coeffs 8-13 in range table for %s (%s)", GetName(), GetType());
192  return kFALSE;
193  }
194  }
195  if (fNoLimits) {
196  // if we ignore nominal validity limits on incident energy, we must still use energy limits
197  // such that all range functions increase monotonically in the energy interval
198  GetRangeFunction(fCoeff[count][0], fCoeff[count][1])->SetRange(GetEminValid(fCoeff[count][0], fCoeff[count][1]), VERY_BIG_ENERGY);
199  Double_t emax = fRange->GetMaximumX() - 1;
200  emax /= fCoeff[count][1];
201  Double_t original_emax = fEmax[count];
202  // the new emax is only accepted if it is > than the nominal emax (400 or 250 AMeV),
203  // and at most 1 GeV/nucleon
204  fEmax[count] = TMath::Min(TMath::Max(original_emax, emax), 1000.);
205  // we may further reduce the upper limit to correspond to the minimum of stopping,
206  // if one exists
207  GetStoppingFunction(fCoeff[count][0], fCoeff[count][1])->SetRange(GetEminValid(fCoeff[count][0], fCoeff[count][1]), GetEmaxValid(fCoeff[count][0], fCoeff[count][1]));
208  emax = fStopping->GetMinimumX();
209  emax /= fCoeff[count][1];
210  // again, the new emax is only accepted if it is > than the nominal emax (400 or 250 AMeV),
211  // and at most 1 GeV/nucleon
212  fEmax[count] = TMath::Min(TMath::Max(original_emax, emax), 1000.);
213  //if(fEmax[count]!=original_emax) Info("ReadRangeTable", "Max. incident E for Z=%d ===> E/A = %f", count+1, fEmax[count]);
214  }
215  if (fgets(line, 132, fp)) {}
216  }
217 
218  return kTRUE;
219 }
220 
221 
222 
228 
230 {
231  // Function parameterising the energy loss of charged particles in this material.
232  //
233  // \param E[0] incident energy given in MeV.
234  // \returns energy loss calculated in MeV.
235 
236  return (E[0] - EResFunc(E, nullptr));
237 }
238 
239 
240 
246 
248 {
249  // Function parameterising the residual energy of charged particles in this material.
250  //
251  // \param E[0] incident energy given in MeV.
252  // \returns residual energy calculated in MeV.
253 
254  Double_t R0 = RangeFunc(E, nullptr);
255  if (R0 < thickness) {// if range < thickness, particle stops: Eres=0
256  fRangeOfLastDE = R0;
257  return 0.0;
258  }
259  else
261 
262  // calculate energy after absorber - invert range function to find Eres corresponding to (R0 - thickness)
263  R0 -= thickness;
264 
265  // invert range function to get energy after absorber
267  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[RF_Z])->GetEnergyPerNucleon(R0, riso) * RF_A;
268  }
269  return fRange->GetX(R0);
270 }
271 
272 
273 
274 
279 
281 {
282  // Return function giving range (in \f$g/cm^2\f$) as a function of energy (in MeV) for
283  // charged particles \f$Z,A\f$ in this material.
284  // \param isoAmat If required, the isotopic mass of the material can be given.
285 
286  RF_Z = Z;
287  RF_A = A;
288  // get parameters for this Z
289  par = &fCoeff[Z - 1];
290  // set up polynomial
291  Double_t x1 = TMath::Log(0.1);
292  Double_t x2 = TMath::Log(0.2);
293  ran = 0.0;
294  for (int j = 2; j < 7; j++)
295  ran += (*par)[j + 1] * TMath::Power(x2, (Double_t)(j - 1));
296  ran += (*par)[2];
297  Double_t y2 = ran;
298  ran = 0.0;
299  for (int jj = 2; jj < 7; jj++)
300  ran += (*par)[jj + 1] * TMath::Power(x1, (Double_t)(jj - 1));
301  ran += (*par)[2];
302  Double_t y1 = ran;
303  adm = (y2 - y1) / (x2 - x1);
304  adn = (y1 - adm * x1);
305  riso = RF_A / (*par)[1];
306  if (isoAmat > 0.0) riso *= (isoAmat / fAmat);
307 
308  fRange->SetRange(0., GetEmaxValid(Z, A));
309 
310  /*
311  * New inversion of range-energy curve (if KVedaLoss::fgNewRangeInversion=kTRUE)
312  */
314  if (!fInvRange[Z]) {
316  }
317  }
318 
319  return fRange;
320 }
321 
322 
323 
328 
330 {
331  // Return function giving stopping power (in \f$MeV/(g/cm^2)\f$) as a function of energy (in MeV) for
332  // charged particles \f$Z,A\f$ in this material.
333  // \param isoAmat If required, the isotopic mass of the material can be given.
334 
335  RF_Z = Z;
336  RF_A = A;
337  // get parameters for this Z
338  par = &fCoeff[Z - 1];
339  // set up polynomial
340  Double_t x1 = TMath::Log(0.1);
341  Double_t x2 = TMath::Log(0.2);
342  ran = 0.0;
343  for (int j = 2; j < 7; j++)
344  ran += (*par)[j + 1] * TMath::Power(x2, (Double_t)(j - 1));
345  ran += (*par)[2];
346  Double_t y2 = ran;
347  ran = 0.0;
348  for (int jj = 2; jj < 7; jj++)
349  ran += (*par)[jj + 1] * TMath::Power(x1, (Double_t)(jj - 1));
350  ran += (*par)[2];
351  Double_t y1 = ran;
352  adm = (y2 - y1) / (x2 - x1);
353  adn = (y1 - adm * x1);
354  riso = RF_A / (*par)[1];
355  if (isoAmat > 0.0) riso *= (isoAmat / fAmat);
356  fStopping->SetRange(0., GetEmaxValid(Z, A));
357  return fStopping;
358 }
359 
360 
361 
366 
368 {
369  // Function parameterising the range of charged particles in this material.
370  // \param E[0] energy is given in MeV.
371  // \returns range calculated in units of \f$g/cm^2\f$
372 
373  eps = E[0] / RF_A;
374  dleps = TMath::Log(eps);
375  if (eps < 0.1)
376  ran = adm * dleps + adn;
377  else {
378  DLEP = dleps;
379  ran = (*par)[2] + (*par)[3] * DLEP;
380  ran += (*par)[4] * (DLEP *= dleps);
381  ran += (*par)[5] * (DLEP *= dleps);
382  ran += (*par)[6] * (DLEP *= dleps);
383  ran += (*par)[7] * (DLEP *= dleps);
384  }
385 
386  // range in g/cm**2
387  return riso * TMath::Exp(ran) * KVUnits::mg;
388 }
389 
390 
391 
396 
398 {
399  // Function parameterising the stopping power of charged particles in this material.
400  // \param E[0] energy is given in MeV.
401  // \returns stopping power calculated in units of \f$Mev/(g/cm^2)\f$
402 
403  eps = E[0] / RF_A;
404  dleps = TMath::Log(eps);
405  if (eps < 0.1) {
406  ran = adm * dleps + adn;
407  drande = E[0] / (riso * TMath::Exp(ran) * KVUnits::mg) / adm;
408  return drande;
409  }
410  DLEP = dleps;
411  ran = (*par)[2] + (*par)[3] * DLEP;
412  drande = (*par)[3];
413  for (int i = 4; i < 8; i++) {
414  drande += (i - 2) * (*par)[i] * DLEP;
415  ran += (*par)[i] * (DLEP *= dleps);
416  }
417  // range in g/cm**2
418  return E[0] / (riso * TMath::Exp(ran) * KVUnits::mg) / drande;
419 }
420 
421 
422 
428 
430 {
431  // Return function giving energy loss (in MeV) as a function of incident energy (in MeV) for
432  // charged particles (Z,A) traversing (or not) material
433  // \param e thickness (in \f$g/cm^2\f$) of this material.
434  // \param isoAmat If required, the isotopic mass of the material can be given.
435 
436  GetRangeFunction(Z, A, isoAmat);
437  thickness = e;
438  fDeltaE->SetRange(0., GetEmaxValid(Z, A));
439  return fDeltaE;
440 }
441 
442 
443 
449 
451 {
452  // Return function giving residual energy (in MeV) as a function of incident energy (in MeV) for
453  // charged particles (Z,A) traversing (or not) material
454  // \param e thickness (in \f$g/cm^2\f$) of this material.
455  // \param isoAmat If required, the isotopic mass of the material can be given.
456 
457  GetRangeFunction(Z, A, isoAmat);
458  thickness = e;
459  fEres->SetRange(0., GetEmaxValid(Z, A));
460  return fEres;
461 }
462 
463 
464 
468 
470 {
471  // \returns range (in \f$g/cm^2\f$) of ion (Z,A) with energy E (MeV) in material.
472  // \param isoAmat change default (isotopic) mass of material,
473 
474  if (Z == 0) return 0.0; //only charged particles
475  TF1* f = GetRangeFunction(Z, A, isoAmat);
476  return f->Eval(E);
477 }
478 
479 
480 
486 
488 {
489  // \returns energy lost (in MeV) by ion (Z,A)
490  // \param E incident energy (MeV)
491  // \param e thickness (in \f$g/cm^2\f$)
492  // \param isoAmat change default (isotopic) mass of material,
493 
494  if (Z == 0) return 0.0; //only charged particles
495  TF1* f = GetDeltaEFunction(e, Z, A, isoAmat);
496  return f->Eval(E);
497 }
498 
499 
500 
506 
508  Double_t isoAmat)
509 {
510  // \returns energy lost (in MeV) by ion (Z,A)
511  // \param E incident energy (MeV)
512  // \param e thickness (in \f$g/cm^2\f$)
513  // \param isoAmat change default (isotopic) mass of material,
514 
515  if (Z == 0) return 0.0; //only charged particles
516  TF1* f = GetEResFunction(e, Z, A, isoAmat);
517  return f->Eval(E);
518 }
519 
520 
521 
527 
529 {
530  // Calculate incident energy (in MeV) for ion (Z,A) for which the range is equal to the
531  // given thickness e (in \f$g/cm^2\f$). At this energy the residual energy of the ion is (just) zero,
532  // for all energies above this energy the residual energy is > 0.
533  // \param isoAmat change default (isotopic) mass of material,
534 
535  if (Z == 0) return 0.0; //only charged particles
536  GetRangeFunction(Z, A, isoAmat);
538  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[Z])->GetEnergyPerNucleon(e, riso) * A;
539  }
540  return fRange->GetX(e);
541 }
542 
543 
544 
548 
550 {
551  // Calculates incident energy (in MeV) of an ion (Z,A) with residual energy Eres (MeV) after thickness e (in \f$g/cm^2\f$).
552  // \param isoAmat change default (isotopic) mass of material,
553 
554  if (Z == 0) return 0.0; //only charged particles
555  GetRangeFunction(Z, A, isoAmat);
556  Double_t R0 = fRange->Eval(Eres) + e;
558  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[Z])->GetEnergyPerNucleon(R0, riso) * A;
559  }
560  return fRange->GetX(R0);
561 }
562 
563 
564 
571 
572 void KVedaLossMaterial::GetParameters(Int_t Zion, Int_t& Aion, std::vector<Double_t>& rangepar)
573 {
574  // For the given ion atomic number, give the reference mass used and the six
575  // parameters for the range function fit
576  // \param[in] Zion ion atomic number
577  // \param[out] Aion reference mass used by table
578  // \param[out] rangepar vector containing the six parameters for the range function fit
579 
580  rangepar = std::vector<Double_t>(fCoeff[Zion - 1].begin() + 2, fCoeff[Zion - 1].end());
581  Aion = fCoeff[Zion - 1][1];
582 }
583 
584 
585 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define VERY_BIG_ENERGY
#define ZMAX_VEDALOSS
maximum atomic number included in range tables
#define f(i)
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
float Float_t
const Bool_t kTRUE
R__EXTERN TEnv * gEnv
int type
char * Form(const char *fmt,...)
virtual const Char_t * GetType() const
Definition: KVBase.h:176
Material for use in energy loss & range calculations.
void AddCompoundElement(Int_t Z, Int_t A, Int_t Natoms)
TF1 * fDeltaE
function parameterising energy loss in material
void AddMixtureElement(Int_t Z, Int_t A, Int_t Natoms, Double_t Proportion)
Double_t fAmat
effective mass number of material
TF1 * fEres
function parameterising residual energy after crossing material
TF1 * fRange
function parameterising range of charged particles in material
TF1 * fStopping
function parameterising stopping power of charged particles in material
Double_t fRangeOfLastDE
range corresponding to last calculated DE
Abstract base class for calculation of range & energy loss of charged particles in matter.
Dedicated optimised inversion of range-energy function for KVedaLoss.
Description of material in the KVedaLoss range table.
Double_t EResFunc(Double_t *, Double_t *)
TObjArray fInvRange
KVedaLossInverseRangeFunction objects.
Float_t GetEminValid(Int_t Z, Int_t A) const
virtual Double_t GetEResOfIon(Int_t Z, Int_t A, Double_t E, Double_t e, Double_t isoAmat=0.)
static KVedaLoss * fgTable
Double_t thickness
in g/cm**2
virtual TF1 * GetDeltaEFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0)
virtual Double_t GetDeltaEOfIon(Int_t Z, Int_t A, Double_t E, Double_t e, Double_t isoAmat=0.)
Double_t StoppingFunc(Double_t *, Double_t *)
Bool_t ReadRangeTable(FILE *fp)
Double_t RF_Z
internal variables used by RangeFunc/DeltaEFunc
std::vector< Double_t > fEmin
Z-dependent minimum energy/nucleon for calculation to be valid.
virtual TF1 * GetRangeFunction(Int_t Z, Int_t A, Double_t isoAmat=0)
virtual ~KVedaLossMaterial()
Destructor.
Double_t RangeFunc(Double_t *, Double_t *)
std::vector< Double_t > fEmax
Z-dependent maximum energy/nucleon for calculation to be valid.
virtual Double_t GetRangeOfIon(Int_t Z, Int_t A, Double_t E, Double_t isoAmat=0.)
static Bool_t fNoLimits
if kTRUE, ignore max E limit for validity of calculation
virtual Double_t GetPunchThroughEnergy(Int_t Z, Int_t A, Double_t e, Double_t isoAmat=0.)
void GetParameters(Int_t Zion, Int_t &Aion, std::vector< Double_t > &rangepar)
std::vector< Double_t > * par
Double_t DeltaEFunc(Double_t *, Double_t *)
KVedaLossMaterial()
Default constructor.
virtual TF1 * GetStoppingFunction(Int_t Z, Int_t A, Double_t isoAmat=0)
Float_t GetEmaxValid(Int_t Z, Int_t A) const
std::vector< std::vector< Double_t > > fCoeff
parameters for range tables
virtual Double_t GetEIncFromEResOfIon(Int_t Z, Int_t A, Double_t Eres, Double_t e, Double_t isoAmat=0.)
virtual TF1 * GetEResFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0)
C++ implementation of VEDALOSS stopping power calculation.
Definition: KVedaLoss.h:62
static Bool_t IsUseNewRangeInversion()
Definition: KVedaLoss.h:91
virtual void SetOwner(Bool_t enable=kTRUE)
virtual const char * GetValue(const char *name, const char *dflt) const
virtual Double_t GetMinimumX(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 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 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 const char * GetName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
TLine * line
const long double mg
Definition: KVUnits.h:74
Double_t Min(Double_t a, Double_t b)
Double_t Exp(Double_t x)
Double_t Power(Double_t x, Double_t y)
constexpr Double_t E()
Double_t Log(Double_t x)
Double_t Max(Double_t a, Double_t b)
#define R0(v, w, x, y, z, i)
auto * a