23#ifndef RestCore_TRestDetectorGas
24#define RestCore_TRestDetectorGas
26#include <TApplication.h>
33#include <TRestMetadata.h>
42#include "TRestDetectorDriftVolume.h"
44#if defined USE_Garfield
45#include <ComponentConstant.hh>
46#include <GeometrySimple.hh>
47#include <MediumMagboltz.hh>
50#include <TrackHeed.hh>
51using namespace Garfield;
56const int RESTGAS_ERROR = -1;
57const int RESTGAS_INTITIALIZED = 0;
58const int RESTGAS_CFG_LOADED = 1;
59const int RESTGAS_GASFILE_LOADED = 2;
65 MediumMagboltz* fGasMedium;
75 Double_t fMaxElectronEnergy;
78 std::vector<TString> fGasComponentName;
80 std::vector<Double_t> fGasComponentFraction;
89 std::vector<Double_t> fEFields;
91 std::vector<Double_t> fBFields;
93 std::vector<Double_t> fAngles;
95 TString fGDMLMaterialRef;
111 Bool_t fTest =
false;
120 void UploadGasToServer(std::string gasFilename);
130 TRestDetectorGas(
const char* configFilename, std::string name =
"",
bool gasGeneration =
false,
150 void CalcGarField(
double Emin,
double Emax,
int n);
152 Int_t
Write(
const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
override;
166 std::cout <<
"REST WARNING. Gas name component n=" << n <<
" requested. But only "
167 <<
GetNofGases() <<
" component(s) in the mixture." << std::endl;
170 return fGasComponentName[n];
176 if (fElectricField == 0) {
177 RESTWarning <<
"TRestDetectorGas::GetDriftVelocity. Warning fElectricField is zero!" <<
RESTendl;
178 RESTWarning <<
" - Use: TRestDetectorGas::SetElectricField( field[V/mm] ) to set the field value"
186 if (fElectricField == 0) {
187 RESTWarning <<
"TRestDetectorGas::GetLongitudinalDiffusion. Warning fElectricField is zero!"
189 RESTWarning <<
" - Use: TRestDetectorGas::SetElectricField( field[V/mm] ) to set the field value"
197 if (fElectricField == 0) {
198 RESTWarning <<
"TRestDetectorGas::GetTransversalDiffusion. Warning fElectricField is zero!"
200 RESTWarning <<
" - Use: TRestDetectorGas::SetElectricField( field[V/mm] ) to set the field value"
206 Double_t GetTownsendCoefficient()
override {
207 if (fElectricField == 0) {
208 RESTWarning <<
"TRestDetectorGas::GetTownsendCoefficient. Warning fElectricField is zero!"
210 RESTWarning <<
" - Use: TRestDetectorGas::SetElectricField( field[V/mm] ) to set the field value"
217 if (fElectricField == 0) {
218 RESTWarning <<
"TRestDetectorGas::GetAttachmentCoefficient. Warning fElectricField is zero!"
220 RESTWarning <<
" - Use: TRestDetectorGas::SetElectricField( field[V/mm] ) to set the field value"
226 void GetGasWorkFunction();
231 std::cout <<
"REST WARNING. Gas fraction for component n=" << n <<
" requested. But only "
232 <<
GetNofGases() <<
" component(s) in the mixture." << std::endl;
236 return fGasComponentFraction[n];
256 inline void PrintGasFileContent() { std::cout <<
fGasFileContent << std::endl; };
std::string ConstructFilename()
Constructs the filename of the pre-generated gas file using the members defined in the RML file.
bool GasFileGenerationEnabled() const
Returns true if the file generation is enabled. False otherwise.
void InitFromConfigFile() override
Loads the gas parameters that define the gas calculation and properties.
void EnableGasGeneration()
void InitFromRootFile() override
Method called after the object is retrieved from root file.
void PrintMetadata() override
Prints the metadata information from the gas.
TString GetGasComponentName(Int_t n)
Returns the gas component n.
Int_t GetNofGases() const
Returns the number of gas elements/compounds present in the gas mixture.
void PrintGasInfo()
Prints the metadata information from the gas.
void PlotTransversalDiffusion(Double_t eMin, Double_t eMax, Int_t nSteps)
It creates a TCanvas where it plots the transversal diffusion as a function of the electric field.
bool GasFileLoaded() const
Returns true if the gas file has been properly loaded. False otherwise.
Double_t GetAttachmentCoefficient(Double_t E)
It returns the attachment coefficient for a given electric field in V/cm.
Double_t GetTransversalDiffusion() override
Returns the transversal diffusion in (cm)^1/2.
Double_t GetMaxElectronEnergy() const
TString GetGasMixture()
Returns a string defining the gas components and fractions.
void SetPressure(Double_t pressure) override
Defines the pressure of the gas.
~TRestDetectorGas()
TRestDetectorGas default destructor.
void PlotTownsendCoefficient(Double_t eMin, Double_t eMax, Int_t nSteps)
It creates a TCanvas where it plots the townsend coefficient as a function of the electric field.
TRestDetectorGas()
TRestDetectorGas default constructor.
Double_t GetDriftVelocity() override
Returns the drift velocity in mm/us.
Double_t GetTownsendCoefficient(Double_t E)
It returns the townsend coefficient for a given electric field in V/cm.
void SetTemperature(Double_t temperature) override
Defines the temperature of the gas.
Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0) override
overwriting the write() method with fStore considered
MediumMagboltz * GetGasMedium() const
Return pointer to Garfield::MediumGas for gas properties.
void LoadGasFile()
It loads a pre-generated gas file corresponding to the gas defined using the RML TRestDetectorGas sec...
void Initialize() override
Defines the metadata section name and initalizes the TRestDetectorGas members.
void AddGasComponent(std::string gasName, Double_t fraction)
Adds a new element/compound to the gas.
std::string FindGasFile(std::string name)
This method tries to find the gas filename given in the argument.
void PlotDriftVelocity(Double_t eMin, Double_t eMax, Int_t nSteps)
It creates a TCanvas where it plots the drift velocity as a function of the electric field.
TString GetGDMLMaterialRef() const
Return reference name of the corresponding material in GDML file.
void PlotLongitudinalDiffusion(Double_t eMin, Double_t eMax, Int_t nSteps)
It creates a TCanvas where it plots the longitudinal diffusion as a function of the electric field.
void GenerateGasFile()
Save a gas file with a structured file name.
Double_t GetGasComponentFraction(Int_t n)
Returns the gas fraction in volume for component n.
void SetMaxElectronEnergy(Double_t energy)
Sets the maximum electron energy to be used in gas generation.
Double_t GetLongitudinalDiffusion() override
Returns the longitudinal diffusion in (cm)^1/2.