3#if defined(__GPUCOMPILE__) || !defined(G_MEDIUM_H)
5#if !defined(__GPUCOMPILE__) && !defined(G_MEDIUM_H)
10#include "GPUInterface.hh"
35class GARFIELD_CLASS_NAME(Medium) {
54 const std::string& GetName()
const {
return m_name; }
56 virtual bool IsGas()
const {
return false; }
58 virtual bool IsSemiconductor()
const {
return false; }
60 virtual bool IsConductor()
const {
return false; }
63 void SetTemperature(
const double t);
65 double GetTemperature()
const {
return m_temperature; }
67 void SetPressure(
const double p);
69 double GetPressure()
const {
return m_pressure; }
71 void SetDielectricConstant(
const double eps);
73 double GetDielectricConstant()
const {
return m_epsilon; }
76 unsigned int GetNumberOfComponents()
const {
return m_nComponents; }
78 virtual void GetComponent(
const unsigned int i, std::string& label,
81 virtual void SetAtomicNumber(
const double z);
83 virtual double GetAtomicNumber()
const {
return m_z; }
85 virtual void SetAtomicWeight(
const double a);
87 virtual double GetAtomicWeight()
const {
return m_a; }
89 virtual void SetNumberDensity(
const double n);
91 virtual double GetNumberDensity()
const {
return m_density; }
93 virtual void SetMassDensity(
const double rho);
95 virtual double GetMassDensity()
const;
98 virtual void EnableDrift(
const bool on =
true) { m_driftable = on; }
100 virtual void EnablePrimaryIonisation(
const bool on =
true) {
106 __DEVICE__ bool IsDriftable()
const {
return m_driftable; }
108 __DEVICE__ bool IsMicroscopic()
const {
return m_microscopic; }
110#ifndef __GPUCOMPILE__
112 bool IsIonisable()
const {
return m_ionisable; }
115 void SetW(
const double w) { m_w = w; }
117 double GetW()
const {
return m_w; }
119 void SetFanoFactor(
const double f) { m_fano = f; }
121 double GetFanoFactor()
const {
return m_fano; }
124 void PlotVelocity(
const std::string& carriers, TPad* pad);
126 void PlotDiffusion(
const std::string& carriers, TPad* pad);
128 void PlotTownsend(
const std::string& carriers, TPad* pad);
130 void PlotAttachment(
const std::string& carriers, TPad* pad);
132 void PlotAlphaEta(
const std::string& carriers, TPad* pad);
136 virtual bool ElectronVelocity(
const double ex,
const double ey,
137 const double ez,
const double bx,
138 const double by,
const double bz,
double& vx,
139 double& vy,
double& vz);
142 virtual bool ElectronVelocityFluxBulk(
const double ex,
const double ey,
143 const double ez,
const double bx,
144 const double by,
const double bz,
145 double& wv,
double& wr);
147 virtual bool ElectronDiffusion(
const double ex,
const double ey,
148 const double ez,
const double bx,
149 const double by,
const double bz,
double& dl,
154 virtual bool ElectronDiffusion(
const double ex,
const double ey,
155 const double ez,
const double bx,
156 const double by,
const double bz,
159 virtual bool ElectronTownsend(
const double ex,
const double ey,
160 const double ez,
const double bx,
161 const double by,
const double bz,
164 virtual bool ElectronAttachment(
const double ex,
const double ey,
165 const double ez,
const double bx,
166 const double by,
const double bz,
169 virtual bool ElectronTOFIonisation(
const double ex,
const double ey,
170 const double ez,
const double bx,
171 const double by,
const double bz,
174 virtual bool ElectronTOFAttachment(
const double ex,
const double ey,
175 const double ez,
const double bx,
176 const double by,
const double bz,
179 virtual bool ElectronLorentzAngle(
const double ex,
const double ey,
180 const double ez,
const double bx,
181 const double by,
const double bz,
184 virtual double ElectronMobility();
188 virtual double GetElectronEnergy(
const double px,
const double py,
189 const double pz,
double& vx,
double& vy,
190 double& vz,
const int band = 0);
193 virtual void GetElectronMomentum(
const double e,
double& px,
double& py,
194 double& pz,
int& band);
197 virtual double GetElectronNullCollisionRate(
const int band = 0);
202 __device__ cuda_t GetElectronCollisionRate(
const cuda_t e,
const int band);
204 __device__ bool ElectronCollision(
const cuda_t e,
int& type,
int& level,
205 cuda_t& e1, cuda_t& dx, cuda_t& dy,
206 cuda_t& dz, Particle* secondaries_type,
207 cuda_t* secondaries_energy,
208 int& num_secondaries,
int& ndxc,
int& band);
211 virtual double GetElectronCollisionRate(
const double e,
const int band = 0);
216 double distance = 0.;
219 virtual bool ElectronCollision(
const double e,
int& type,
int& level,
220 double& e1,
double& dx,
double& dy,
double& dz,
221 std::vector<Secondary>& secondaries,
225#ifndef __GPUCOMPILE__
228 virtual bool HoleVelocity(
const double ex,
const double ey,
const double ez,
229 const double bx,
const double by,
const double bz,
230 double& vx,
double& vy,
double& vz);
232 virtual bool HoleDiffusion(
const double ex,
const double ey,
const double ez,
233 const double bx,
const double by,
const double bz,
234 double& dl,
double& dt);
236 virtual bool HoleDiffusion(
const double ex,
const double ey,
const double ez,
237 const double bx,
const double by,
const double bz,
240 virtual bool HoleTownsend(
const double ex,
const double ey,
const double ez,
241 const double bx,
const double by,
const double bz,
244 virtual bool HoleAttachment(
const double ex,
const double ey,
const double ez,
245 const double bx,
const double by,
const double bz,
248 virtual double HoleMobility();
252 virtual bool IonVelocity(
const double ex,
const double ey,
const double ez,
253 const double bx,
const double by,
const double bz,
254 double& vx,
double& vy,
double& vz);
255 bool HasIonVelocity()
const {
return !(m_iVel.empty() && m_iMob.empty()); }
257 virtual bool IonDiffusion(
const double ex,
const double ey,
const double ez,
258 const double bx,
const double by,
const double bz,
259 double& dl,
double& dt);
261 virtual bool IonDissociation(
const double ex,
const double ey,
262 const double ez,
const double bx,
263 const double by,
const double bz,
double& diss);
265 virtual double IonMobility();
268 virtual bool NegativeIonVelocity(
const double ex,
const double ey,
269 const double ez,
const double bx,
270 const double by,
const double bz,
double& vx,
271 double& vy,
double& vz);
273 virtual double NegativeIonMobility();
276 void SetFieldGrid(
double emin,
double emax,
const size_t ne,
bool logE,
277 double bmin = 0.,
double bmax = 0.,
const size_t nb = 1,
278 double amin = HalfPi,
double amax = HalfPi,
279 const size_t na = 1);
281 void SetFieldGrid(
const std::vector<double>& efields,
282 const std::vector<double>& bfields,
283 const std::vector<double>& angles);
285 void GetFieldGrid(std::vector<double>& efields, std::vector<double>& bfields,
286 std::vector<double>& angles);
289 bool SetElectronVelocityE(
const size_t ie,
const size_t ib,
const size_t ia,
291 return SetEntry(ie, ib, ia,
"ElectronVelocityE", m_eVelE, v);
294 bool GetElectronVelocityE(
const size_t ie,
const size_t ib,
const size_t ia,
296 return GetEntry(ie, ib, ia,
"ElectronVelocityE", m_eVelE, v);
299 bool SetElectronVelocityExB(
const size_t ie,
const size_t ib,
const size_t ia,
301 return SetEntry(ie, ib, ia,
"ElectronVelocityExB", m_eVelX, v);
304 bool GetElectronVelocityExB(
const size_t ie,
const size_t ib,
const size_t ia,
306 return GetEntry(ie, ib, ia,
"ElectronVelocityExB", m_eVelX, v);
309 bool SetElectronVelocityB(
const size_t ie,
const size_t ib,
const size_t ia,
311 return SetEntry(ie, ib, ia,
"ElectronVelocityB", m_eVelB, v);
314 bool GetElectronVelocityB(
const size_t ie,
const size_t ib,
const size_t ia,
316 return GetEntry(ie, ib, ia,
"ElectronVelocityB", m_eVelB, v);
319 bool SetElectronFluxVelocity(
const size_t ie,
const size_t ib,
320 const size_t ia,
const double v) {
321 return SetEntry(ie, ib, ia,
"ElectronFluxVelocity", m_eVelWv, v);
324 bool GetElectronFluxVelocity(
const size_t ie,
const size_t ib,
325 const size_t ia,
double& v) {
326 return GetEntry(ie, ib, ia,
"ElectronFluxVelocity", m_eVelWv, v);
329 bool SetElectronBulkVelocity(
const size_t ie,
const size_t ib,
330 const size_t ia,
const double v) {
331 return SetEntry(ie, ib, ia,
"ElectronBulkVelocity", m_eVelWr, v);
334 bool GetElectronBulkVelocity(
const size_t ie,
const size_t ib,
335 const size_t ia,
double& v) {
336 return GetEntry(ie, ib, ia,
"ElectronBulkVelocity", m_eVelWr, v);
339 bool SetElectronLongitudinalDiffusion(
const size_t ie,
const size_t ib,
340 const size_t ia,
const double dl) {
341 return SetEntry(ie, ib, ia,
"ElectronLongitudinalDiffusion", m_eDifL, dl);
344 bool GetElectronLongitudinalDiffusion(
const size_t ie,
const size_t ib,
345 const size_t ia,
double& dl) {
346 return GetEntry(ie, ib, ia,
"ElectronLongitudinalDiffusion", m_eDifL, dl);
349 bool SetElectronTransverseDiffusion(
const size_t ie,
const size_t ib,
350 const size_t ia,
const double dt) {
351 return SetEntry(ie, ib, ia,
"ElectronTransverseDiffusion", m_eDifT, dt);
354 bool GetElectronTransverseDiffusion(
const size_t ie,
const size_t ib,
355 const size_t ia,
double& dt) {
356 return GetEntry(ie, ib, ia,
"ElectronTransverseDiffusion", m_eDifT, dt);
359 bool SetElectronTownsend(
const size_t ie,
const size_t ib,
const size_t ia,
360 const double alpha) {
361 return SetEntry(ie, ib, ia,
"ElectronTownsend", m_eAlp, alpha);
364 bool GetElectronTownsend(
const size_t ie,
const size_t ib,
const size_t ia,
366 return GetEntry(ie, ib, ia,
"ElectronTownsend", m_eAlp, alpha);
369 bool SetElectronAttachment(
const size_t ie,
const size_t ib,
const size_t ia,
371 return SetEntry(ie, ib, ia,
"ElectronAttachment", m_eAtt, eta);
374 bool GetElectronAttachment(
const size_t ie,
const size_t ib,
const size_t ia,
376 return GetEntry(ie, ib, ia,
"ElectronAttachment", m_eAtt, eta);
379 bool SetElectronTOFIonisation(
const size_t ie,
const size_t ib,
380 const size_t ia,
const double v) {
381 return SetEntry(ie, ib, ia,
"ElectronTOFIonisation", m_eRIon, v);
384 bool GetElectronTOFIonisation(
const size_t ie,
const size_t ib,
385 const size_t ia,
double& v) {
386 return GetEntry(ie, ib, ia,
"ElectronTOFIonisation", m_eRIon, v);
389 bool SetElectronTOFAttachment(
const size_t ie,
const size_t ib,
390 const size_t ia,
const double v) {
391 return SetEntry(ie, ib, ia,
"ElectronTOFAttachment", m_eRAtt, v);
394 bool GetElectronTOFAttachment(
const size_t ie,
const size_t ib,
395 const size_t ia,
double& v) {
396 return GetEntry(ie, ib, ia,
"ElectronTOFAttachment", m_eRAtt, v);
400 bool SetElectronLorentzAngle(
const size_t ie,
const size_t ib,
401 const size_t ia,
const double lor) {
402 return SetEntry(ie, ib, ia,
"ElectronLorentzAngle", m_eLor, lor);
405 bool GetElectronLorentzAngle(
const size_t ie,
const size_t ib,
406 const size_t ia,
double& lor) {
407 return GetEntry(ie, ib, ia,
"ElectronLorentzAngle", m_eLor, lor);
411 bool SetHoleVelocityE(
const size_t ie,
const size_t ib,
const size_t ia,
413 return SetEntry(ie, ib, ia,
"HoleVelocityE", m_hVelE, v);
416 bool GetHoleVelocityE(
const size_t ie,
const size_t ib,
const size_t ia,
418 return GetEntry(ie, ib, ia,
"HoleVelocityE", m_hVelE, v);
421 bool SetHoleVelocityExB(
const size_t ie,
const size_t ib,
const size_t ia,
423 return SetEntry(ie, ib, ia,
"HoleVelocityExB", m_hVelX, v);
426 bool GetHoleVelocityExB(
const size_t ie,
const size_t ib,
const size_t ia,
428 return GetEntry(ie, ib, ia,
"HoleVelocityExB", m_hVelX, v);
431 bool SetHoleVelocityB(
const size_t ie,
const size_t ib,
const size_t ia,
433 return SetEntry(ie, ib, ia,
"HoleVelocityB", m_hVelB, v);
436 bool GetHoleVelocityB(
const size_t ie,
const size_t ib,
const size_t ia,
438 return GetEntry(ie, ib, ia,
"HoleVelocityB", m_hVelB, v);
442 bool SetHoleLongitudinalDiffusion(
const size_t ie,
const size_t ib,
443 const size_t ia,
const double dl) {
444 return SetEntry(ie, ib, ia,
"HoleLongitudinalDiffusion", m_hDifL, dl);
447 bool GetHoleLongitudinalDiffusion(
const size_t ie,
const size_t ib,
448 const size_t ia,
double& dl) {
449 return GetEntry(ie, ib, ia,
"HoleLongitudinalDiffusion", m_hDifL, dl);
452 bool SetHoleTransverseDiffusion(
const size_t ie,
const size_t ib,
453 const size_t ia,
const double dt) {
454 return SetEntry(ie, ib, ia,
"HoleTransverseDiffusion", m_hDifT, dt);
457 bool GetHoleTransverseDiffusion(
const size_t ie,
const size_t ib,
458 const size_t ia,
double& dt) {
459 return GetEntry(ie, ib, ia,
"HoleTransverseDiffusion", m_hDifT, dt);
462 bool SetHoleTownsend(
const size_t ie,
const size_t ib,
const size_t ia,
463 const double alpha) {
464 return SetEntry(ie, ib, ia,
"HoleTownsend", m_hAlp, alpha);
467 bool GetHoleTownsend(
const size_t ie,
const size_t ib,
const size_t ia,
469 return GetEntry(ie, ib, ia,
"HoleTownsend", m_hAlp, alpha);
472 bool SetHoleAttachment(
const size_t ie,
const size_t ib,
const size_t ia,
474 return SetEntry(ie, ib, ia,
"HoleAttachment", m_hAtt, eta);
477 bool GetHoleAttachment(
const size_t ie,
const size_t ib,
const size_t ia,
479 return GetEntry(ie, ib, ia,
"HoleAttachment", m_hAtt, eta);
486 bool SetIonMobility(
const std::vector<double>& fields,
487 const std::vector<double>& mobilities,
488 const bool negativeIons =
false);
490 bool SetIonMobility(
const size_t ie,
const size_t ib,
const size_t ia,
493 bool GetIonMobility(
const size_t ie,
const size_t ib,
const size_t ia,
495 return GetEntry(ie, ib, ia,
"IonMobility", m_iMob, mu);
499 bool SetIonLongitudinalDiffusion(
const size_t ie,
const size_t ib,
500 const size_t ia,
const double dl) {
501 return SetEntry(ie, ib, ia,
"IonLongitudinalDiffusion", m_iDifL, dl);
504 bool GetIonLongitudinalDiffusion(
const size_t ie,
const size_t ib,
505 const size_t ia,
double& dl) {
506 return GetEntry(ie, ib, ia,
"IonLongitudinalDiffusion", m_iDifL, dl);
509 bool SetIonTransverseDiffusion(
const size_t ie,
const size_t ib,
510 const size_t ia,
const double dt) {
511 return SetEntry(ie, ib, ia,
"IonTransverseDiffusion", m_iDifT, dt);
514 bool GetIonTransverseDiffusion(
const size_t ie,
const size_t ib,
515 const size_t ia,
double& dt) {
516 return GetEntry(ie, ib, ia,
"IonTransverseDiffusion", m_iDifT, dt);
519 bool SetIonDissociation(
const size_t ie,
const size_t ib,
const size_t ia,
521 return SetEntry(ie, ib, ia,
"IonDissociation", m_iDis, diss);
524 bool GetIonDissociation(
const size_t ie,
const size_t ib,
const size_t ia,
526 return GetEntry(ie, ib, ia,
"IonDissociation", m_iDis, diss);
530 bool SetNegativeIonMobility(
const size_t ie,
const size_t ib,
const size_t ia,
533 bool GetNegativeIonMobility(
const size_t ie,
const size_t ib,
const size_t ia,
535 return GetEntry(ie, ib, ia,
"NegativeIonMobility", m_nMob, mu);
539 virtual void ResetTables();
541 void ResetElectronVelocity() {
548 void ResetElectronDiffusion() {
553 void ResetElectronTownsend() { m_eAlp.clear(); }
554 void ResetElectronAttachment() { m_eAtt.clear(); }
555 void ResetElectronTOFRates() {
559 void ResetElectronLorentzAngle() { m_eLor.clear(); }
561 void ResetHoleVelocity() {
566 void ResetHoleDiffusion() {
571 void ResetHoleTownsend() { m_hAlp.clear(); }
572 void ResetHoleAttachment() { m_hAtt.clear(); }
574 void ResetIonMobility() {
578 void ResetIonDiffusion() {
582 void ResetIonDissociation() { m_iDis.clear(); }
583 void ResetNegativeIonMobility() {
588 void VelocityFromMobility(
589 const std::vector<std::vector<std::vector<double> > >& mob,
590 std::vector<std::vector<std::vector<double> > >& vel);
594 void SetExtrapolationMethodVelocity(
const std::string& extrLow,
595 const std::string& extrHigh);
596 void SetExtrapolationMethodDiffusion(
const std::string& extrLow,
597 const std::string& extrHigh);
598 void SetExtrapolationMethodTownsend(
const std::string& extrLow,
599 const std::string& extrHigh);
600 void SetExtrapolationMethodAttachment(
const std::string& extrLow,
601 const std::string& extrHigh);
602 void SetExtrapolationMethodIonMobility(
const std::string& extrLow,
603 const std::string& extrHigh);
604 void SetExtrapolationMethodIonDissociation(
const std::string& extrLow,
605 const std::string& extrHigh);
608 void SetInterpolationMethodVelocity(
const unsigned int intrp);
609 void SetInterpolationMethodDiffusion(
const unsigned int intrp);
610 void SetInterpolationMethodTownsend(
const unsigned int intrp);
611 void SetInterpolationMethodAttachment(
const unsigned int intrp);
612 void SetInterpolationMethodIonMobility(
const unsigned int intrp);
613 void SetInterpolationMethodIonDissociation(
const unsigned int intrp);
616 virtual double ScaleElectricField(
const double e)
const {
return e; }
617 virtual double UnScaleElectricField(
const double e)
const {
return e; }
618 virtual double ScaleVelocity(
const double v)
const {
return v; }
619 virtual double ScaleDiffusion(
const double d)
const {
return d; }
620 virtual double ScaleDiffusionTensor(
const double d)
const {
return d; }
621 virtual double ScaleTownsend(
const double alpha)
const {
return alpha; }
622 virtual double ScaleAttachment(
const double eta)
const {
return eta; }
623 virtual double ScaleLorentzAngle(
const double lor)
const {
return lor; }
624 virtual double ScaleDissociation(
const double diss)
const {
return diss; }
628 virtual bool GetOpticalDataRange(
double& emin,
double& emax,
629 const unsigned int i = 0);
631 virtual bool GetDielectricFunction(
const double e,
double& eps1,
double& eps2,
632 const unsigned int i = 0);
634 virtual bool GetPhotoAbsorptionCrossSection(
const double e,
double& sigma,
635 const unsigned int i = 0);
636 virtual double GetPhotonCollisionRate(
const double e);
637 virtual bool PhotonCollision(
const double e,
int& type,
int& level,
638 double& e1,
double& ctheta,
639 std::vector<Secondary>& secondaries);
642 void EnableDebugging() { m_debug =
true; }
643 void DisableDebugging() { m_debug =
false; }
646 virtual double CreateGPUTransferObject(MediumGPU*& med_gpu);
649 std::string m_className =
"Medium";
651 static int m_idCounter;
654 unsigned int m_nComponents = 1;
656 std::string m_name =
"";
658 double m_temperature = 293.15;
660 double m_pressure = 760.;
662 double m_epsilon = 1.;
668 double m_density = 0.;
675 bool m_driftable =
false;
676 bool m_microscopic =
false;
677 bool m_ionisable =
false;
679#ifndef __GPUCOMPILE__
687 bool m_isChanged =
true;
690 bool m_debug =
false;
693 bool m_tab2d =
false;
696 std::vector<double> m_eFields;
697 std::vector<double> m_bFields;
698 std::vector<double> m_bAngles;
701 std::vector<std::vector<std::vector<double> > > m_eVelE;
702 std::vector<std::vector<std::vector<double> > > m_eVelX;
703 std::vector<std::vector<std::vector<double> > > m_eVelB;
704 std::vector<std::vector<std::vector<double> > > m_eDifL;
705 std::vector<std::vector<std::vector<double> > > m_eDifT;
706 std::vector<std::vector<std::vector<double> > > m_eAlp;
707 std::vector<std::vector<std::vector<double> > > m_eAtt;
708 std::vector<std::vector<std::vector<double> > > m_eLor;
709 std::vector<std::vector<std::vector<double> > > m_eVelWv;
710 std::vector<std::vector<std::vector<double> > > m_eVelWr;
711 std::vector<std::vector<std::vector<double> > > m_eRIon;
712 std::vector<std::vector<std::vector<double> > > m_eRAtt;
714 std::vector<std::vector<std::vector<std::vector<double> > > > m_eDifM;
717 std::vector<std::vector<std::vector<double> > > m_hVelE;
718 std::vector<std::vector<std::vector<double> > > m_hVelX;
719 std::vector<std::vector<std::vector<double> > > m_hVelB;
720 std::vector<std::vector<std::vector<double> > > m_hDifL;
721 std::vector<std::vector<std::vector<double> > > m_hDifT;
722 std::vector<std::vector<std::vector<double> > > m_hAlp;
723 std::vector<std::vector<std::vector<double> > > m_hAtt;
725 std::vector<std::vector<std::vector<std::vector<double> > > > m_hDifM;
728 std::vector<std::vector<std::vector<double> > > m_iMob;
729 std::vector<std::vector<std::vector<double> > > m_iVel;
730 std::vector<std::vector<std::vector<double> > > m_iDifL;
731 std::vector<std::vector<std::vector<double> > > m_iDifT;
732 std::vector<std::vector<std::vector<double> > > m_iDis;
734 std::vector<std::vector<std::vector<double> > > m_nMob;
735 std::vector<std::vector<std::vector<double> > > m_nVel;
738 unsigned int m_eThrAlp = 0;
739 unsigned int m_eThrAtt = 0;
740 unsigned int m_hThrAlp = 0;
741 unsigned int m_hThrAtt = 0;
742 unsigned int m_iThrDis = 0;
745 std::pair<unsigned int, unsigned int> m_extrVel = {0, 1};
746 std::pair<unsigned int, unsigned int> m_extrDif = {0, 1};
747 std::pair<unsigned int, unsigned int> m_extrAlp = {0, 1};
748 std::pair<unsigned int, unsigned int> m_extrAtt = {0, 1};
749 std::pair<unsigned int, unsigned int> m_extrLor = {0, 1};
750 std::pair<unsigned int, unsigned int> m_extrMob = {0, 1};
751 std::pair<unsigned int, unsigned int> m_extrDis = {0, 1};
754 unsigned int m_intpVel = 2;
755 unsigned int m_intpDif = 2;
756 unsigned int m_intpAlp = 2;
757 unsigned int m_intpAtt = 2;
758 unsigned int m_intpLor = 2;
759 unsigned int m_intpMob = 2;
760 unsigned int m_intpDis = 2;
762 bool Velocity(
const double ex,
const double ey,
const double ez,
763 const double bx,
const double by,
const double bz,
764 const std::vector<std::vector<std::vector<double> > >& velE,
765 const std::vector<std::vector<std::vector<double> > >& velB,
766 const std::vector<std::vector<std::vector<double> > >& velX,
767 const double q,
double& vx,
double& vy,
double& vz)
const;
768 bool VelocityFluxBulk(
769 const double ex,
const double ey,
const double ez,
const double bx,
770 const double by,
const double bz,
771 const std::vector<std::vector<std::vector<double> > >& velWv,
772 const std::vector<std::vector<std::vector<double> > >& velWr,
773 double& wv,
double& wr)
const;
774 static void Langevin(
const double ex,
const double ey,
const double ez,
775 double bx,
double by,
double bz,
const double mu,
776 double& vx,
double& vy,
double& vz);
777 static void Langevin(
const double ex,
const double ey,
const double ez,
778 double bx,
double by,
double bz,
const double mu,
779 const double muH,
double& vx,
double& vy,
double& vz);
780 bool Diffusion(
const double ex,
const double ey,
const double ez,
781 const double bx,
const double by,
const double bz,
782 const std::vector<std::vector<std::vector<double> > >& difL,
783 const std::vector<std::vector<std::vector<double> > >& difT,
784 double& dl,
double& dt)
const;
786 const double ex,
const double ey,
const double ez,
const double bx,
787 const double by,
const double bz,
788 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
789 double cov[3][3])
const;
790 bool Alpha(
const double ex,
const double ey,
const double ez,
const double bx,
791 const double by,
const double bz,
792 const std::vector<std::vector<std::vector<double> > >& tab,
793 unsigned int intp,
const unsigned int thr,
794 const std::pair<unsigned int, unsigned int>& extr,
795 double& alpha)
const;
796 double GetAngle(
const double ex,
const double ey,
const double ez,
797 const double bx,
const double by,
const double bz,
798 const double e,
const double b)
const;
799 bool Interpolate(
const double e,
const double b,
const double a,
800 const std::vector<std::vector<std::vector<double> > >& table,
801 double& y,
const unsigned int intp,
802 const std::pair<unsigned int, unsigned int>& extr,
803 const bool logval =
false)
const;
805 double Interpolate1D(
const double e,
const std::vector<double>& table,
806 const std::vector<double>& fields,
807 const unsigned int intpMeth,
808 const std::pair<unsigned int, unsigned int>& extr,
809 const bool logval =
false)
const;
811 bool SetEntry(
const size_t i,
const size_t j,
const size_t k,
812 const std::string& fcn,
813 std::vector<std::vector<std::vector<double> > >& tab,
815 bool GetEntry(
const size_t i,
const size_t j,
const size_t k,
816 const std::string& fcn,
817 const std::vector<std::vector<std::vector<double> > >& tab,
820 void SetExtrapolationMethod(
const std::string& low,
const std::string& high,
821 std::pair<unsigned int, unsigned int>& extr,
822 const std::string& fcn);
823 bool GetExtrapolationIndex(std::string str,
unsigned int& nb)
const;
825 const std::vector<std::vector<std::vector<double> > >& tab)
const;
827 void Clone(std::vector<std::vector<std::vector<double> > >& tab,
828 const std::vector<double>& efields,
829 const std::vector<double>& bfields,
830 const std::vector<double>& angles,
const unsigned int intp,
831 const std::pair<unsigned int, unsigned int>& extr,
832 const double init,
const std::string& label);
833 void Clone(std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
834 const size_t n,
const std::vector<double>& efields,
835 const std::vector<double>& bfields,
836 const std::vector<double>& angles,
const unsigned int intp,
837 const std::pair<unsigned int, unsigned int>& extr,
838 const double init,
const std::string& label);
840 void Init(
const size_t nE,
const size_t nB,
const size_t nA,
841 std::vector<std::vector<std::vector<double> > >& tab,
843 void Init(
const size_t nE,
const size_t nB,
const size_t nA,
const size_t nT,
844 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
851 friend class MediumGas;
852 friend class MediumMagboltz;
855 enum class MediumType { Medium = 0, MediumGas, MediumMagboltz };
857 MediumType m_MediumType{MediumType::Medium};