Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Heed::HeedMatterDef Class Reference

#include <HeedMatterDef.h>

+ Inheritance diagram for Heed::HeedMatterDef:

Public Member Functions

 HeedMatterDef (void)
 
 HeedMatterDef (EnergyMesh *fenergy_mesh, MatterDef *amatter, AtomPhotoAbsCS *faapacs[], double fW=0.0, double fF=standard_factor_Fano)
 
 HeedMatterDef (EnergyMesh *fenergy_mesh, GasDef *agas, MolecPhotoAbsCS *fampacs[], double fW=0.0, double fF=standard_factor_Fano)
 
 HeedMatterDef (EnergyMesh *fenergy_mesh, const String &gas_notation, MolecPhotoAbsCS *fampacs[], double fW=0.0, double fF=standard_factor_Fano)
 
void replace_epsi12 (const String &file_name)
 
virtual void print (std::ostream &file, int l) const
 
 macro_copy_total (HeedMatterDef)
 

Public Attributes

PassivePtr< MatterDefmatter
 
DynLinArr< PassivePtr< const AtomPhotoAbsCS > > apacs
 
double eldens_cm_3
 
double eldens
 
double xeldens
 
double wpla
 
double radiation_length
 
double Rutherford_const
 
double W
 
double F
 
PassivePtr< EnergyMeshenergy_mesh
 
DynLinArr< double > ACS
 
DynLinArr< double > ICS
 
DynLinArr< double > epsip
 
DynLinArr< double > epsi1
 
DynLinArr< double > epsi2
 
double min_ioniz_pot
 

Detailed Description

Definition at line 37 of file HeedMatterDef.h.

Constructor & Destructor Documentation

◆ HeedMatterDef() [1/4]

Heed::HeedMatterDef::HeedMatterDef ( void  )

Definition at line 13 of file HeedMatterDef.cpp.

14 : eldens_cm_3(0.0),
15 eldens(0.0),
16 xeldens(0.0),
17 wpla(0.0),
20 W(0.0),
21 F(0.0) {
22 ;
23}

◆ HeedMatterDef() [2/4]

Heed::HeedMatterDef::HeedMatterDef ( EnergyMesh fenergy_mesh,
MatterDef amatter,
AtomPhotoAbsCS faapacs[],
double  fW = 0.0,
double  fF = standard_factor_Fano 
)

Definition at line 25 of file HeedMatterDef.cpp.

27 : W(fW), F(fF), energy_mesh(fenergy_mesh) {
28 mfunname("HeedMatterDef::HeedMatterDef(...)");
29 matter.put(amatter);
30 check_econd11(matter->qatom(), <= 0, mcerr);
31 long q = matter->qatom();
32 apacs.put_qel(q);
33 for (long n = 0; n < q; ++n) {
34 apacs[n].put(faapacs[n]);
35 check_econd12(matter->atom(n)->Z(), !=, apacs[n]->get_Z(), mcerr);
36 }
37 check_econd11(F, == 0.0, mcerr);
38 if (W == 0.0) {
39#ifdef CALC_W_USING_CHARGES
40 double mean_I = 0.0;
41 double d = 0.0;
42 for (long n = 0; n < q; ++n) {
43 mean_I +=
44 matter->weight_quan(n) * apacs[n]->get_Z() * apacs[n]->get_I_min();
45 d += matter->weight_quan(n) * apacs[n]->get_Z();
46 }
47 W = coef_I_to_W * mean_I / d;
48#else
49 double mean_I = 0.0;
50 for (long n = 0; n < q; ++n) {
51 mean_I += matter->weight_quan(n) * apacs[n]->get_I_min();
52 }
53 W = coef_I_to_W * mean_I;
54#endif
55 }
56 inite_HeedMatterDef();
57}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:366
#define check_econd12(a, sign, b, stream)
Definition: FunNameStack.h:380
#define mfunname(string)
Definition: FunNameStack.h:67
PassivePtr< MatterDef > matter
Definition: HeedMatterDef.h:39
PassivePtr< EnergyMesh > energy_mesh
Definition: HeedMatterDef.h:50
DynLinArr< PassivePtr< const AtomPhotoAbsCS > > apacs
Definition: HeedMatterDef.h:40
const double coef_I_to_W
Definition: PhotoAbsCS.h:552
#define mcerr
Definition: prstream.h:135

◆ HeedMatterDef() [3/4]

Heed::HeedMatterDef::HeedMatterDef ( EnergyMesh fenergy_mesh,
GasDef agas,
MolecPhotoAbsCS fampacs[],
double  fW = 0.0,
double  fF = standard_factor_Fano 
)

Definition at line 59 of file HeedMatterDef.cpp.

61 : W(fW), F(fF), energy_mesh(fenergy_mesh) {
62 mfunname("HeedMatterDef::HeedMatterDef(...)");
63 matter.put(agas);
64 check_econd11(agas->qmolec(), <= 0, mcerr);
65 long qat = agas->qatom();
66 apacs.put_qel(qat);
67 const long qmol = agas->qmolec();
68 long nat = 0;
69 for (long nmol = 0; nmol < qmol; ++nmol) {
70 check_econd12(agas->molec(nmol)->tqatom(), !=, fampacs[nmol]->get_qatom(),
71 mcerr);
72 // number of different atoms in mol
73 const long qa = agas->molec(nmol)->qatom();
74 for (long na = 0; na < qa; ++na) {
75 apacs[nat].put(fampacs[nmol]->get_atom(na).getver());
76 check_econd12(apacs[nat]->get_Z(), !=, agas->molec(nmol)->atom(na)->Z(),
77 mcerr);
78 nat++;
79 }
80 }
81 if (F == 0.0) {
82#ifdef CALC_W_USING_CHARGES
83 double u = 0.0;
84 double d = 0.0;
85 for (long n = 0; n < qmol; ++n) {
86 u += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z() *
87 fampacs[n]->get_F();
88 d += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
89 }
90 F = u / d;
91#else
92 for (long n = 0; n < qmol; ++n) {
93 F += agas->weight_quan_molec(n) * fampacs[n]->get_F();
94 }
95#endif
96 }
97
98 if (W == 0.0) {
99#ifdef CALC_W_USING_CHARGES
100 double u = 0.0;
101 double d = 0.0;
102 for (long n = 0; n < qmol; ++n) {
103 u += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z() *
104 fampacs[n]->get_W();
105 d += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
106 }
107 W = u / d;
108#else
109 for (long n = 0; n < qmol; ++n) {
110 W += agas->weight_quan_molec(n) * fampacs[n]->get_W();
111 }
112#endif
113 }
114 inite_HeedMatterDef();
115}

◆ HeedMatterDef() [4/4]

Heed::HeedMatterDef::HeedMatterDef ( EnergyMesh fenergy_mesh,
const String gas_notation,
MolecPhotoAbsCS fampacs[],
double  fW = 0.0,
double  fF = standard_factor_Fano 
)

Definition at line 117 of file HeedMatterDef.cpp.

120 : W(fW), F(fF), energy_mesh(fenergy_mesh) {
121 mfunnamep("HeedMatterDef::HeedMatterDef(...)");
122 MatterDef* amat = MatterDef::get_MatterDef(gas_notation);
123 GasDef* agas = dynamic_cast<GasDef*>(amat);
124 if (agas == NULL) {
125 funnw.ehdr(mcerr);
126 mcerr << "notation supplied as the gas notation is not appear "
127 << "to be related to gas \n";
128 mcerr << "gas_notation=" << gas_notation << '\n';
129 spexit(mcerr);
130 }
131
132 matter.put(agas);
133 check_econd11(agas->qmolec(), <= 0, mcerr);
134 long qat = agas->qatom();
135 apacs.put_qel(qat);
136 long qmol = agas->qmolec();
137 long nat = 0;
138 for (long nmol = 0; nmol < qmol; ++nmol) {
139 check_econd12(agas->molec(nmol)->tqatom(), !=, fampacs[nmol]->get_qatom(),
140 mcerr);
141 long qa = agas->molec(nmol)->qatom(); //quantity of different atoms in mol
142 for (long na = 0; na < qa; ++na) {
143 apacs[nat].put(fampacs[nmol]->get_atom(na).getver());
144 check_econd12(apacs[nat]->get_Z(), !=, agas->molec(nmol)->atom(na)->Z(),
145 mcerr);
146 nat++;
147 }
148 }
149 if (F == 0.0) {
150#ifdef CALC_W_USING_CHARGES
151 double u = 0.0;
152 double d = 0.0;
153 for (long n = 0; n < qmol; ++n) {
154 u += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z() *
155 fampacs[n]->get_F();
156 d += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
157 }
158 F = u / d;
159#else
160 for (long n = 0; n < qmol; ++n) {
161 F += agas->weight_quan_molec(n) * fampacs[n]->get_F();
162 }
163#endif
164 }
165
166 if (W == 0.0) {
167#ifdef CALC_W_USING_CHARGES
168 double u = 0.0;
169 double d = 0.0;
170 for (long n = 0; n < qmol; ++n) {
171 u += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z() *
172 fampacs[n]->get_W();
173 d += agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
174 }
175 W = u / d;
176#else
177 for (long n = 0; n < qmol; ++n) {
178 W += agas->weight_quan_molec(n) * fampacs[n]->get_W();
179 }
180#endif
181 }
182 inite_HeedMatterDef();
183}
#define mfunnamep(string)
Definition: FunNameStack.h:77
#define spexit(stream)
Definition: FunNameStack.h:536
static MatterDef * get_MatterDef(const String &fnotation)
Definition: MatterDef.cpp:132

Member Function Documentation

◆ macro_copy_total()

Heed::HeedMatterDef::macro_copy_total ( HeedMatterDef  )

◆ print()

void Heed::HeedMatterDef::print ( std::ostream &  file,
int  l 
) const
virtual

Definition at line 366 of file HeedMatterDef.cpp.

366 {
367 if (l <= 0) return;
368 Ifile << "HeedMatterDef:\n";
369 indn.n += 2;
370 matter->print(file, 1);
371 if (l >= 2) {
372 long q = matter->qatom();
373 Ifile << "Printing " << q << " photoabsorption cross sections:\n";
374 indn.n += 2;
375 for (long n = 0; n < q; ++n) {
376 apacs[n]->print(file, l - 1);
377 }
378 indn.n -= 2;
379 }
380 Iprintan(file, eldens_cm_3, "1/cm^3");
381 Iprintan(file, eldens, "MeV^3");
382 Iprintan(file, xeldens, "MeV^2/cm");
383 Iprintn(file, wpla);
385 Iprintan(file, Rutherford_const, "1/cm^3");
386 Iprintn(file, W);
387 Iprintn(file, F);
388 Iprintn(file, min_ioniz_pot);
389 Iprintn(file, energy_mesh->get_q());
390 if (l >= 2) {
391 long qe = energy_mesh->get_q();
392 long ne;
393 indn.n += 2;
394 Ifile << " ne energy ACS(Mb) ICS(Mb) ACS(1/MeV^2) "
395 "ICS(1/MeV^2) epsip epsi1 epsi2 "
396 "(1+epsi1)^2+epsi2^2\n";
397 for (ne = 0; ne < qe; ne++) {
398 //double et = pow(energy_mesh->get_ec(ne), 2.0);
399 Ifile << std::setw(3) << ne << ' ' << std::setw(12)
400 << energy_mesh->get_e(ne) << ' ' << std::setw(12) << ACS[ne] << ' '
401 << std::setw(12) << ICS[ne] << ' ' << std::setw(12)
402 << ACS[ne] * C1_MEV2_MBN << ' ' << std::setw(12)
403 << ICS[ne] * C1_MEV2_MBN << ' ' << std::setw(12) << epsip[ne] << ' '
404 << std::setw(12) << epsi1[ne] << ' '
405 // << std::setw(12) << epsip[ne] * et << ' '
406 // << std::setw(12) << epsi1[ne] * et << ' '
407 << std::setw(12) << epsi2[ne] << ' ' << std::setw(12)
408 << pow((1 + epsi1[ne]), 2.0) + pow(epsi2[ne], 2.0) << " \n";
409 }
410 indn.n -= 2;
411 }
412 indn.n -= 2;
413}
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336
DynLinArr< double > epsi2
Definition: HeedMatterDef.h:62
DynLinArr< double > epsip
Definition: HeedMatterDef.h:56
DynLinArr< double > ACS
Definition: HeedMatterDef.h:54
DynLinArr< double > ICS
Definition: HeedMatterDef.h:55
DynLinArr< double > epsi1
Definition: HeedMatterDef.h:61
const double C1_MEV2_MBN
indentation indn
Definition: prstream.cpp:13
#define Iprintan(file, name, addition)
Definition: prstream.h:223
#define Ifile
Definition: prstream.h:207
#define Iprintn(file, name)
Definition: prstream.h:216

◆ replace_epsi12()

void Heed::HeedMatterDef::replace_epsi12 ( const String file_name)

Definition at line 314 of file HeedMatterDef.cpp.

314 {
315 mfunnamep("void HeedMatterDef::replace_epsi12(const String& file_name)");
316
317#ifdef USE_STLSTRING
318 std::ifstream file(file_name.c_str());
319#else
320 std::ifstream file(hist_file_name);
321#endif
322 if (!file) {
323 funnw.ehdr(mcerr);
324 mcerr << "cannot open file " << file_name << std::endl;
325 spexit(mcerr);
326 } else {
327 mcout << "file " << file_name << " is opened" << std::endl;
328 }
329 long qe = 0; // number of points in input mesh
330 file >> qe;
331 check_econd11(qe, <= 2, mcerr);
332
333 DynLinArr<double> ener(qe);
334 DynLinArr<double> eps1(qe);
335 DynLinArr<double> eps2(qe);
336
337 for (long ne = 0; ne < qe; ++ne) {
338 file >> ener[ne] >> eps1[ne] >> eps2[ne];
339 check_econd11(eps2[ne], < 0.0, mcerr);
340 if (ne > 0) {
341 check_econd12(ener[ne], <, ener[ne - 1], mcerr);
342 }
343 }
344
345 PointCoorMesh<double, DynLinArr<double> > pcmd(qe, &(ener));
346 double emin = ener[0] - 0.5 * (ener[1] - ener[0]);
347 double emax = ener[qe - 1] + 0.5 * (ener[qe - 1] - ener[qe - 2]);
348
349 qe = energy_mesh->get_q();
350 for (long ne = 0; ne < qe; ++ne) {
351 double ec = energy_mesh->get_ec(ne);
354 pcmd, // dimension q
355 eps1, // dimension q-1
356 ec, 0, 1, emin, 1, emax);
359 pcmd, // dimension q
360 eps2, // dimension q-1
361 ec, 1, 1, emin, 1, emax);
362 // Iprint3n(mcout, ec, epsi1[ne], epsi2[ne]);
363 }
364}
#define mcout
Definition: prstream.h:133
T t_value_straight_point_ar(const M &mesh, const D &y, T x, int s_ban_neg, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
Definition: tline.h:2511

Member Data Documentation

◆ ACS

DynLinArr<double> Heed::HeedMatterDef::ACS

Definition at line 54 of file HeedMatterDef.h.

Referenced by print().

◆ apacs

◆ eldens

double Heed::HeedMatterDef::eldens

Definition at line 43 of file HeedMatterDef.h.

Referenced by print().

◆ eldens_cm_3

double Heed::HeedMatterDef::eldens_cm_3

Definition at line 42 of file HeedMatterDef.h.

Referenced by print().

◆ energy_mesh

PassivePtr<EnergyMesh> Heed::HeedMatterDef::energy_mesh

◆ epsi1

DynLinArr<double> Heed::HeedMatterDef::epsi1

Definition at line 61 of file HeedMatterDef.h.

Referenced by print(), and replace_epsi12().

◆ epsi2

DynLinArr<double> Heed::HeedMatterDef::epsi2

Definition at line 62 of file HeedMatterDef.h.

Referenced by print(), and replace_epsi12().

◆ epsip

DynLinArr<double> Heed::HeedMatterDef::epsip

Definition at line 56 of file HeedMatterDef.h.

Referenced by print().

◆ F

double Heed::HeedMatterDef::F

Definition at line 49 of file HeedMatterDef.h.

Referenced by Garfield::TrackHeed::GetFanoFactor(), HeedMatterDef(), and print().

◆ ICS

DynLinArr<double> Heed::HeedMatterDef::ICS

Definition at line 55 of file HeedMatterDef.h.

Referenced by print().

◆ matter

◆ min_ioniz_pot

double Heed::HeedMatterDef::min_ioniz_pot

Definition at line 63 of file HeedMatterDef.h.

Referenced by print().

◆ radiation_length

double Heed::HeedMatterDef::radiation_length

Definition at line 46 of file HeedMatterDef.h.

Referenced by print().

◆ Rutherford_const

double Heed::HeedMatterDef::Rutherford_const

Definition at line 47 of file HeedMatterDef.h.

Referenced by print().

◆ W

double Heed::HeedMatterDef::W

Definition at line 48 of file HeedMatterDef.h.

Referenced by Garfield::TrackHeed::GetW(), HeedMatterDef(), and print().

◆ wpla

double Heed::HeedMatterDef::wpla

Definition at line 45 of file HeedMatterDef.h.

Referenced by print().

◆ xeldens

double Heed::HeedMatterDef::xeldens

Definition at line 44 of file HeedMatterDef.h.

Referenced by print().


The documentation for this class was generated from the following files: