17 #include "TRestDetectorGarfieldDriftProcess.h"
24 #if defined USE_Garfield
25 using namespace Garfield;
26 #include <ComponentConstant.hh>
28 #include "TGDMLParse.h"
35 const double cmTomm = 10.;
37 #if defined USE_Garfield
39 TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess() : fRandom(0), fGfSensor(0) {
50 TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess(
const char* configFilename)
51 : fRandom(0), fGfSensor(0) {
54 if (LoadConfigFromFile(configFilename) == -1) LoadDefaultConfig();
62 TRestDetectorGarfieldDriftProcess::~TRestDetectorGarfieldDriftProcess() {
65 delete fOutputHitsEvent;
68 void TRestDetectorGarfieldDriftProcess::LoadDefaultConfig() {
69 SetName(
"garfieldDriftProcess-Default");
70 SetTitle(
"Default config");
72 cout <<
"TRestDetectorGarfieldDriftProcess metadata not found. Loading default values" << endl;
74 fDriftPotential = 1000;
78 void TRestDetectorGarfieldDriftProcess::LoadConfig(
const string& configFilename,
const string& name) {
81 if (fDriftPotential == PARAMETER_NOT_FOUND_DBL) {
82 fDriftPotential = 1000;
92 fRandom =
new TRandom3(0);
93 fInputHitsEvent =
nullptr;
96 #if defined USE_Garfield
105 #if defined USE_Garfield
115 fGas = GetMetadata<TRestDetectorGas>();
116 if (fGas !=
nullptr) {
117 if (fGasPressure <= 0)
118 fGasPressure = fGas->GetPressure();
120 fGas->SetPressure(fGasPressure);
123 cout <<
"REST_WARNING. No TRestDetectorGas found in TRestRun." << endl;
127 fReadout = GetMetadata<TRestDetectorReadout>();
128 if (fReadout ==
nullptr) {
129 cout <<
"REST ERROR: Readout has not been initialized" << endl;
134 if (fGDML_Filename ==
"") {
135 cout <<
"REST ERROR no geometry GDML file name parameter (gdml_file) given "
136 "in TRestDetectorGarfieldDriftProcess "
141 TGeoVolume* geovol =
nullptr;
143 fGeometry->InitGfGeometry();
145 geovol = TGDMLParse::StartGDML(fGDML_Filename);
147 cout <<
"REST ERROR when TRestDetectorGeometry read GDML file " << fGDML_Filename << endl;
151 fGeometry->SetTopVolume(geovol);
152 fGeometry->CloseGeometry();
153 fGeometry->DefaultColors();
156 cout <<
"TRestDetectorGarfieldDriftProcess GetVerboseLevel : "
161 TObjArray* thenodes = geovol->GetNodes();
162 for (TIter it = thenodes->begin(); it != thenodes->end(); ++it) {
163 TGeoNode* itnode = (TGeoNode*)(*it);
165 cout <<
"****** itnode " << itnode->GetName() << endl;
166 itnode->PrintCandidates();
168 cout <<
"****** itnode " << itnode->GetName() << endl;
169 itnode->PrintCandidates();
171 TGeoVolume* itvol = itnode->GetVolume();
173 cout <<
" * * itvolume " << itvol->GetName() << endl;
176 cout <<
" * * itvolume " << itvol->GetName() << endl;
179 TGeoMedium* itmed = itvol->GetMedium();
181 cout <<
" * * itmed " << itmed->GetName() << endl;
184 if (fGas->GetGDMLMaterialRef() == itmed->GetName()) {
185 fGeometry->SetGfGeoMedium(itmed->GetName(), fGas);
187 cout <<
" -> gas volume SetMedium itmed " << itmed->GetName() <<
" fGas "
188 << fGas->GetGasMedium()->GetName() << endl;
189 fGas->GetGasMedium()->PrintGas();
191 cout <<
" -> gas volume SetMedium itmed " << itmed->GetName() <<
" fGas "
192 << fGas->GetGasMedium()->GetName() << endl;
193 fGas->GetGasMedium()->PrintGas();
197 if ((strncmp(itvol->GetName(),
"anodeVol", 8) == 0) ||
198 (strncmp(itvol->GetName(),
"cathodeVol", 10) == 0)) {
199 fGeometry->SetDriftElecNode(itnode);
201 cout <<
" -> cathode volume " << endl;
202 cout <<
" -> cathode volume " << endl;
206 if ((strncmp(itvol->GetName(),
"micromegasVol", 13) == 0)) {
207 fGeometry->AddReadoutElecNode(itnode);
209 cout <<
" -> readout volume " << endl;
210 cout <<
" -> readout volume " << endl;
215 double matrixZpos = 0, planeZpos = 0, planeZvec = 0;
216 cout <<
"fGeometry " << fGeometry <<
" nb node " << fGeometry->GetNReadoutElecNodes() << endl;
217 for (
int ii = 0; ii < fGeometry->GetNReadoutElecNodes(); ii++) {
218 bool rdPlaneFound =
false;
220 TGeoNode* readoutnode = fGeometry->GetReadoutElecNode(ii);
221 cout <<
"ii " << ii <<
" readoutnode " << readoutnode <<
" visible " << readoutnode->IsVisible()
223 if (!readoutnode)
continue;
224 cout <<
"readoutnode visible " << readoutnode->IsVisible() << endl;
226 readoutnode->PrintCandidates();
227 TGeoMatrix* readoutmatrix = readoutnode->GetMatrix();
228 if (!readoutmatrix)
continue;
229 cout <<
"readoutmatrix " << endl;
230 readoutmatrix->Print();
231 matrixZpos = (readoutmatrix->IsTranslation() ? 10 * readoutmatrix->GetTranslation()[2]
233 for (
int jj = 0; jj < fReadout->GetNumberOfReadoutPlanes(); jj++) {
236 cout <<
" jj " << jj <<
" matrixZpos " << matrixZpos <<
" planeZpos " << planeZpos << endl;
237 if (fabs(planeZpos - matrixZpos) > 1)
240 planeZvec = readoutplane->
GetNormal().z();
241 cout <<
" planeZvec " << planeZvec << endl;
245 double field = fDriftPotential /
246 (fGeometry->GetDriftElecNode()->GetMatrix()->GetTranslation()[2] - planeZpos);
247 ComponentConstant* cmp =
new ComponentConstant();
249 cmp->SetElectricField(0, 0, field);
250 cmp->SetPotential(0, 0, planeZpos / 10., 0);
251 cmp->SetWeightingField(0, 0, planeZvec,
"m");
252 fGeometry->AddGfComponent(cmp);
253 cout <<
"add component field " << field <<
" planeZpos " << planeZpos <<
" planeZvec "
254 << planeZvec <<
" drift pos "
255 << fGeometry->GetDriftElecNode()->GetMatrix()->GetTranslation()[2] << endl;
261 if (fGeometry->GetNbOfGfComponent() > 0) {
262 fGfSensor =
new Sensor();
263 fGfSensor->AddComponent(fGeometry->GetGfComponent(0));
264 fGfSensor->AddElectrode(fGeometry->GetGfComponent(0),
"m");
265 fGfSensor->SetTimeWindow(0., 0.1, 500);
267 double xmin = -1e9, xmax = 1e9, ymin = -1e9, ymax = 1e9;
268 TGeoShape* readoutshape = fGeometry->GetReadoutElecNode(0)->GetVolume()->GetShape();
269 TGeoMatrix* readoutmatrix = fGeometry->GetReadoutElecNode(0)->GetMatrix();
270 double xmid = 10. * readoutmatrix->GetTranslation()[0],
271 ymid = 10. * readoutmatrix->GetTranslation()[1];
272 if (readoutshape->InheritsFrom(
"TGeoBBox")) {
274 TGeoBBox* readoutbox = (TGeoBBox*)readoutshape;
275 xmin = xmid - 10. * readoutbox->GetDX() - 100;
276 xmax = xmid + 10. * readoutbox->GetDX() + 100;
277 ymin = ymid - 10. * readoutbox->GetDY() - 100;
278 ymax = ymid + 10. * readoutbox->GetDY() + 100;
280 double driftelecZpos = fGeometry->GetDriftElecNode()
282 ->GetTranslation()[2];
284 fGfSensor->SetArea(xmin / 10., ymin / 10., planeZpos / 10. + fStopDistance, xmax / 10., ymax / 10.,
286 fGeometry->AddGfSensor(fGfSensor);
287 printf(
" GfSensor area x- %lf y- %lf z- %lf x+ %lf y+ %lf z+ %lf\n", xmin / 10., ymin / 10.,
288 planeZpos / 10. + 0.2, xmax / 10., ymax / 10., driftelecZpos);
290 fGfDriftMethod =
new DRIFT_METHOD();
293 fGfDriftMethod->SetSensor(fGfSensor);
295 fGfDriftMethod->SetCollisionSteps(100);
298 if (!fGfDriftMethod || !fGfSensor) {
299 cout <<
"REST ERROR Garfield not fully initialized in "
300 "TRestDetectorGarfieldDriftProcess::InitProcess "
302 cout <<
" fGfDriftMethod " << fGfDriftMethod <<
" fGfSensor " << fGfSensor <<
" GetNbOfGfComponent "
303 << fGeometry->GetNbOfGfComponent() << endl;
304 cout <<
" matrixZpos " << matrixZpos <<
" planeZpos " << planeZpos <<
" planeZvec " << planeZvec
330 Int_t TRestDetectorGarfieldDriftProcess::FindModule(Int_t readoutPlane, Double_t x, Double_t y) {
334 if ((*plane)[md].IsInside({x, y})) {
345 #if defined USE_Garfield
348 double x, y, z, energy;
349 double xi, yi, zi, ti, xf, yf, zf, tf, energyf;
353 cout <<
"Number of hits : " << fInputHitsEvent->GetNumberOfHits() << endl;
354 cout <<
"--------------------------" << endl;
358 for (
unsigned int hit = 0; hit < fInputHitsEvent->GetNumberOfHits(); hit++) {
359 x = fInputHitsEvent->
GetX(hit);
360 y = fInputHitsEvent->
GetY(hit);
361 z = fInputHitsEvent->
GetZ(hit);
362 energy = fInputHitsEvent->GetEnergy(hit);
368 Medium* tmed = fGeometry->GetGfMedium(x, y, z);
370 cout <<
"no medium found at x " << x <<
" y " << y <<
" z " << z << endl;
377 double Wfactor = tmed->GetW();
378 double nb_electrons_real = energy * 1000. / Wfactor / fPEReduction;
379 unsigned int nb_electrons =
380 fRandom->Poisson(nb_electrons_real);
387 for (
unsigned int iel = 0; iel < nb_electrons; iel++) {
389 fGfDriftMethod->DriftElectron(x / 10., y / 10., z / 10., 0);
391 if (fGfDriftMethod->GetNumberOfElectronEndpoints() == 0) {
392 cout <<
" TRestDetectorGarfieldDriftProcess::ProcessEvent: no electron end "
399 fGfDriftMethod->GetElectronEndpoint(0, xi, yi, zi, ti, xf, yf, zf, tf,
407 energyf = fPEReduction * Wfactor / 1000.;
408 fOutputHitsEvent->
AddHit(xf, yf, zf, energyf, tf);
418 if (fOutputHitsEvent->GetNumberOfHits() == 0)
return nullptr;
423 cout <<
"TRestDetectorElectronDiffusionProcess. Hits added : " << fOutputHitsEvent->GetNumberOfHits()
425 cout <<
"TRestDetectorElectronDiffusionProcess. Hits total energy : "
426 << fOutputHitsEvent->GetTotalEnergy() << endl;
427 cout <<
" fTimedHitsEvent " << fOutputHitsEvent <<
" class " << fOutputHitsEvent->ClassName() << endl;
431 return fOutputHitsEvent;
434 fOutputHitsEvent = fInputHitsEvent;
435 RESTDebug <<
"nullptr process" <<
RESTendl;
441 #if defined USE_Garfield
455 fDriftPotential = GetDblParameterWithUnits(
"driftPotential");
void Initialize() override
Making default settings.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
virtual void PrintEvent() const
Double_t GetX(int n) const
Returns the X-coordinate of hit entry n in mm.
Double_t GetY(int n) const
Returns the Y-coordinate of hit entry n in mm.
Double_t GetZ(int n) const
Returns the Z-coordinate of hit entry n in mm.
void AddHit(Double_t x, Double_t y, Double_t z, Double_t en, Double_t t=0, REST_HitType type=XYZ)
Adds a new hit to this event.
TVector3 GetNormal() const
Returns a TVector3 with a vector normal to the readout plane.
TVector3 GetPosition() const
Returns a TVector3 with the readout plane position.
size_t GetNumberOfModules() const
Returns the total number of modules in the readout plane.
A metadata class to generate/store a readout description.
virtual void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
virtual void EndProcess()
To be executed at the end of the run (outside event loop)
virtual void InitProcess()
To be executed at the beginning of the run (outside event loop)
A base class for any REST event.
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
Double_t StringToDouble(std::string in)
Gets a double from a string.