KaliVeda
Toolkit for HIC analysis
KVedaLoss.cpp
1 //Created by KVClassFactory on Wed Feb 2 15:49:27 2011
2 //Author: frankland,,,,
3 
4 #include "KVedaLoss.h"
5 #include "KVedaLossMaterial.h"
6 #include "KVedaLossRangeFitter.h"
7 #include <TString.h>
8 #include <TSystem.h>
9 #include <TEnv.h>
10 #include <KVSystemDirectory.h>
11 #include <KVSystemFile.h>
12 #include "TGeoMaterial.h"
13 
15 
18 
19 
35 
37 {
38  // Call this static method with yes=kTRUE in order to recalculate the nominal limits
39  // on incident ion energies for which the range tables are valid.
40  //
41  // Normally all range, \f$dE\f$, \f$E_{res}\f$ functions are limited to range \f$0\leq E\leq E_{max}\f$,
42  // where \f$E_{max}\f$ is nominal maximum energy for which range tables are valid
43  // (usually 400MeV/u for \f$Z<3\f$, 250MeV/u for \f$Z>3\f$).
44  //
45  // If higher energies are required, call static method KVedaLoss::SetIgnoreEnergyLimits() **BEFORE ANY MATERIALS ARE CREATED**
46  // in order to recalculate the \f$E_{max}\f$ limits in such a way that:
47  // - range function is always monotonically increasing function of \f$E_{inc}\f$;
48  // - stopping power is concave (i.e. no minimum of stopping power followed by an increase)
49  //
50  // Then, at the most, the new limit will be 1 GeV/nucleon, or
51  // at the least, it will remain at the nominal (400 or 250 MeV/nucleon) level.
53 }
54 
55 
56 
58 
60 {
62 }
63 
64 
65 
68 
70  : KVIonRangeTable("VEDALOSS",
71  "Calculation of range and energy loss of charged particles in matter using VEDALOSS range tables")
72 {
73  // Default constructor
74 
76  if (!CheckMaterialsList()) {
77  Error("KVedaLoss", "Problem reading range tables. Do not use.");
78  }
79 }
80 
81 
82 
85 
87 {
88  // Destructor
89 }
90 
91 
92 
98 
100 {
101  // PRIVATE method - called to initialize fMaterials list of all known materials properties
102  //
103  // Any files in $(WORKING_DIR)/vedaloss/*.dat will also be read, these contain
104  // materials added by the user(s)
105 
106  fMaterials = new KVHashList;
107  fMaterials->SetName("VEDALOSS materials list");
108  fMaterials->SetOwner();
109 
110  auto read_files_in_dir = [ = ](const TString & dir) {
111  bool ok = true;
112  KVSystemDirectory matDir("matDir", dir);
113  TIter nxtfil(matDir.GetListOfFiles());
114  KVSystemFile* fil;
115  while ((fil = (KVSystemFile*)nxtfil())) {
116  if (TString(fil->GetName()).EndsWith(".dat")) ok = ok && ReadMaterials(fil->GetFullPath());
117  }
118  return ok;
119  };
120 
121  bool ok = true;
122  // read all materials in installation directory
123  ok &= read_files_in_dir(GetDATADIRFilePath("data/vedaloss"));
124 
126  ok &= read_files_in_dir(fLocalMaterialsDirectory);
127  }
128  return ok;
129 }
130 
131 
132 
135 
136 Bool_t KVedaLoss::ReadMaterials(const Char_t* DataFilePath) const
137 {
138  // Read and add range tables for materials in file
139 
140  Char_t name[25], gtype[25], state[10];
141  Float_t Amat = 0.;
142  Float_t Dens = 0.;
143  Float_t MoleWt = 0.;
144  Float_t Temp = 19.;
145  Float_t Zmat = 0.;
146  int mat_count = 0;
147 
148  FILE* fp;
149  if (!(fp = fopen(DataFilePath, "r"))) {
150  Error("init_materials()", "Range tables file %s cannot be opened", DataFilePath);
151  return kFALSE;
152  }
153  else {
154  char line[132];
155  while (fgets(line, 132, fp)) { // read lines from file
156 
157  switch (line[0]) {
158 
159  case '/': // ignore comment lines
160  break;
161 
162  case '+': // header lines
163 
164  if (sscanf(line, "+ %s %s %s %f %f %f %f %f",
165  gtype, name, state, &Dens, &Zmat, &Amat,
166  &MoleWt, &Temp)
167  != 8) {
168  Error("init_materials()", "Problem reading file %s", DataFilePath);
169  fclose(fp);
170  return kFALSE;
171  }
172  //found a new material
173  KVedaLossMaterial* tmp_mat = new KVedaLossMaterial(this, name, gtype, state, Dens,
174  Zmat, Amat, MoleWt);
175  if (!tmp_mat->ReadRangeTable(fp)) {
176  delete tmp_mat;
177  return kFALSE;
178  }
179  fMaterials->Add(tmp_mat);
180  tmp_mat->Initialize();
181  ++mat_count;
182  if (tmp_mat->IsGas()) tmp_mat->SetTemperatureAndPressure(19., 1.*KVUnits::atm);
183  break;
184  }
185  }
186  fclose(fp);
187  }
188  return kTRUE;
189 }
190 
191 
192 
203 
205 {
206  // Use the RANGE tables (see KVRangeYanez) to generate a new material composed of a single chemical element.
207  //
208  // \param[in] Z atomic number of element \f$Z\f$
209  // \param[in] A [optional] mass number of isotope \f$A\f$
210  //
211  // If the isotope \a A is not specified, we create a material containing the naturally
212  // occuring isotopes of the given element, weighted according to their natural abundances.
213  //
214  // If the element name is "X", this material will be called "natX", for "naturally-occuring X".
215 
216  std::unique_ptr<KVIonRangeTable> yanez(KVIonRangeTable::GetRangeTable("RANGE"));
217  KVIonRangeTableMaterial* mat = yanez->AddElementalMaterial(Z, A);
218  AddMaterial(mat);
219  return GetMaterial(mat->GetName());
220 }
221 
222 
223 
226 
228 {
229  // If the given material is defined in the RANGE tables, import it into VEDALOSS
230 
231  std::unique_ptr<KVIonRangeTable> yanez(KVIonRangeTable::GetRangeTable("RANGE"));
232  if (yanez->GetMaterial(name)) {
233  AddMaterial(yanez->GetMaterial(name));
234  return kTRUE;
235  }
236  return kFALSE;
237 }
238 
239 
240 
251 
252 KVIonRangeTableMaterial* KVedaLoss::AddCompoundMaterial(const Char_t* name, const Char_t* symbol, Int_t nel, Int_t* Z, Int_t* A, Int_t* nat, Double_t dens) const
253 {
254  // Use the RANGE tables (see KVRangeYanez) to add a compound material with a simple formula composed of different elements
255  //
256  // \param[in] name name for the new compound (no spaces)
257  // \param[in] symbol chemical symbol for compound
258  // \param[in] nelem number of elements in compound
259  // \param[in] z[nelem] atomic numbers of elements
260  // \param[in] a[nelem] mass numbers of elements
261  // \param[in] natoms[nelem] number of atoms of each element
262  // \param[in] density in \f$g/cm^{3}\f$, only required if compound is a solid
263 
264  std::unique_ptr<KVIonRangeTable> yanez(KVIonRangeTable::GetRangeTable("RANGE"));
265  KVIonRangeTableMaterial* mat = yanez->AddCompoundMaterial(name, symbol, nel, Z, A, nat, dens);
266  AddMaterial(mat);
267  return GetMaterial(mat->GetName());
268 }
269 
270 
271 
283 
284 KVIonRangeTableMaterial* KVedaLoss::AddMixedMaterial(const Char_t* name, const Char_t* symbol, Int_t nel, Int_t* Z, Int_t* A, Int_t* nat, Double_t* prop, Double_t dens) const
285 {
286  // Use the RANGE tables (see KVRangeYanez) to add a material which is a mixture of either elements or compounds:
287  //
288  // \param[in] name name for the new mixture (no spaces)
289  // \param[in] symbol chemical symbol for mixture
290  // \param[in] nel number of elements in mixture
291  // \param[in] z[nel] atomic numbers of elements
292  // \param[in] a[nel] mass numbers of elements
293  // \param[in] nat[nel] number of atoms of each element
294  // \param[in] prop[nel] proportion by mass in mixture of element
295  // \param[in] density in \f$g/cm^{3}\f$, if mixture is a solid
296 
297  std::unique_ptr<KVIonRangeTable> yanez(KVIonRangeTable::GetRangeTable("RANGE"));
298  KVIonRangeTableMaterial* mat = yanez->AddMixedMaterial(name, symbol, nel, Z, A, nat, prop, dens);
299  AddMaterial(mat);
300  return GetMaterial(mat->GetName());
301 }
302 
303 
304 
305 
308 
310 {
311  // Returns pointer to material of given name or type.
313  if (!M) {
315  }
316  return M;
317 }
318 
319 
320 
329 
331 {
332  // Add a material (taken from a different range table) to VEDALOSS
333  // This means fitting the ranges for Z=1-100 and writing the parameters in a
334  // file which will be stored in
335  //
336  // $(WORKING_DIR)/VEDALOSS/[name].dat
337  //
338  // which will be read at each initialisation to include the new material
339 
340  KVedaLossRangeFitter vlfit;
341  vlfit.SetMaterial(mat);
342  TString matname = mat->GetName();
343  matname.ReplaceAll(" ", "_"); //no spaces in filename
344  matname += ".dat";
345  // check directory exists & make it if necessary
349  }
350  vlfit.DoFits(Form("%s/%s", fLocalMaterialsDirectory.Data(), matname.Data()));
351  ReadMaterials(Form("%s/%s", fLocalMaterialsDirectory.Data(), matname.Data()));
352 }
353 
354 
355 
357 
359 {
360  printf("KVedaLoss::%s\n%s\n", GetName(), GetTitle());
361  printf("\nEnergy loss & range tables loaded for %d materials:\n\n", fMaterials->GetEntries());
362  fMaterials->ls();
363 }
364 
365 
366 
371 
373 {
374  // Create and fill a list of all materials for which range tables exist.
375  // Each entry is a TNamed with the name and type (title) of the material.
376  // User's responsibility to delete list after use (it owns its objects).
377 
378  if (CheckMaterialsList()) {
379  TObjArray* list = new TObjArray(fMaterials->GetEntries());
380  list->SetOwner(kTRUE);
381  TIter next(fMaterials);
382  KVedaLossMaterial* mat;
383  while ((mat = (KVedaLossMaterial*)next())) {
384  list->Add(new TNamed(mat->GetName(), mat->GetType()));
385  }
386  return list;
387  }
388  return 0;
389 }
390 
391 
int Int_t
bool Bool_t
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h prop
char name[80]
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
virtual const Char_t * GetType() const
Definition: KVBase.h:177
static const Char_t * GetDATADIRFilePath(const Char_t *namefile="")
Definition: KVBase.cpp:74
static const Char_t * GetWORKDIRFilePath(const Char_t *namefile="")
Definition: KVBase.cpp:127
Extended version of ROOT THashList.
Definition: KVHashList.h:29
Material for use in energy loss & range calculations.
void SetTemperatureAndPressure(Double_t T, Double_t P)
Abstract base class for calculation of range & energy loss of charged particles in matter.
KVIonRangeTableMaterial * GetMaterial(const Char_t *material) const
Returns pointer to material of given name or type.
static KVIonRangeTable * GetRangeTable(const Char_t *name)
Generates an instance of the KVIonRangeTable plugin class corresponding to given name.
virtual void SetOwner(Bool_t enable=kTRUE)
virtual TObject * FindObjectByType(const Char_t *) const
virtual void Add(TObject *obj)
virtual TObject * FindObject(const char *name) const
Extension of ROOT TSystemDirectory class, handling browsing directories on disk.
virtual TList * GetListOfFiles() const
Extended ROOT TSystemFile with added info on file size etc.
Definition: KVSystemFile.h:18
const Char_t * GetFullPath() const
Definition: KVSystemFile.h:49
Description of material in the KVedaLoss range table.
static void SetNoLimits(Bool_t on=kTRUE)
Bool_t ReadRangeTable(FILE *fp)
static Bool_t CheckIon(Int_t Z)
Fit a range table using the VEDALOSS functional.
void SetMaterial(KVIonRangeTableMaterial *m)
Sets range table to fit. Also finds material with closest Z in VEDALOSS library.
void DoFits(TString output_file, Int_t Zmin=1, Int_t Zmax=100)
Perform fits to range tables for elements from Zmin to Zmax.
C++ implementation of VEDALOSS stopping power calculation.
Definition: KVedaLoss.h:63
KVIonRangeTableMaterial * AddCompoundMaterial(const Char_t *, const Char_t *, Int_t, Int_t *, Int_t *, Int_t *, Double_t=-1.0) const
Definition: KVedaLoss.cpp:252
TObjArray * GetListOfMaterials()
Definition: KVedaLoss.cpp:372
Bool_t init_materials() const
Definition: KVedaLoss.cpp:99
KVIonRangeTableMaterial * GetMaterialWithNameOrType(const Char_t *material) const
Returns pointer to material of given name or type.
Definition: KVedaLoss.cpp:309
KVedaLoss()
Default constructor.
Definition: KVedaLoss.cpp:69
virtual ~KVedaLoss()
Destructor.
Definition: KVedaLoss.cpp:86
static KVHashList * fMaterials
static list of all known materials
Definition: KVedaLoss.h:64
KVIonRangeTableMaterial * AddElementalMaterial(Int_t Z, Int_t A=0) const
Definition: KVedaLoss.cpp:204
TString fLocalMaterialsDirectory
Definition: KVedaLoss.h:65
void AddMaterial(KVIonRangeTableMaterial *) const
Definition: KVedaLoss.cpp:330
Bool_t ReadMaterials(const Char_t *path) const
Read and add range tables for materials in file.
Definition: KVedaLoss.cpp:136
Bool_t CheckIon(Int_t Z, Int_t A) const
Definition: KVedaLoss.cpp:59
KVIonRangeTableMaterial * AddMixedMaterial(const Char_t *, const Char_t *, Int_t, Int_t *, Int_t *, Int_t *, Double_t *, Double_t=-1.0) const
Definition: KVedaLoss.cpp:284
Bool_t AddRANGEMaterial(const Char_t *name) const
If the given material is defined in the RANGE tables, import it into VEDALOSS.
Definition: KVedaLoss.cpp:227
static Bool_t fgNewRangeInversion
static flag for using new KVedaLossInverseRangeFunction
Definition: KVedaLoss.h:74
static void SetIgnoreEnergyLimits(Bool_t yes=kTRUE)
Definition: KVedaLoss.cpp:36
void Print(Option_t *="") const
Definition: KVedaLoss.cpp:358
Bool_t CheckMaterialsList() const
Definition: KVedaLoss.h:68
void ls(Option_t *option="") const override
void SetName(const char *name)
virtual Int_t GetEntries() const
virtual void SetOwner(Bool_t enable=kTRUE)
const char * GetName() const override
const char * GetTitle() const override
void Add(TObject *obj) override
virtual void Error(const char *method, const char *msgfmt,...) const
const char * Data() const
TString & ReplaceAll(const char *s1, const char *s2)
virtual int Chmod(const char *file, UInt_t mode)
virtual int mkdir(const char *name, Bool_t recursive=kFALSE)
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
TLine * line
const long double atm
Definition: KVUnits.h:79
ClassImp(TPyArg)