REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorGarfieldDriftProcess.cxx
1
16
17#include "TRestDetectorGarfieldDriftProcess.h"
18
19#include <TGeoBBox.h>
20#include <TRandom3.h>
21
23
24#if defined USE_Garfield
25using namespace Garfield;
26#include <ComponentConstant.hh>
27
28#include "TGDMLParse.h"
29#endif
30
31#include <stdio.h>
32
33using namespace std;
34
35const double cmTomm = 10.;
36
37#if defined USE_Garfield
38
39TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess() : fRandom(0), fGfSensor(0) {
40 Initialize();
41}
42
43// __________________________________________________________
44// TODO : Perhaps this constructor should be removed
45// since we will always load the config from TRestRun
46// when we use AddProcess. It would be necessary only if we use the
47// process stand alone but even then we could just call LoadConfig
48// __________________________________________________________
49
50TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess(const char* configFilename)
51 : fRandom(0), fGfSensor(0) {
52 Initialize();
53
54 if (LoadConfigFromFile(configFilename) == -1) LoadDefaultConfig();
55
56 PrintMetadata();
57
58 if (fReadout == nullptr) fReadout = new TRestDetectorReadout(configFilename);
59}
60
61// TRestDetectorGarfieldDriftProcess destructor
62TRestDetectorGarfieldDriftProcess::~TRestDetectorGarfieldDriftProcess() {
63 delete fReadout;
64 delete fGeometry;
65 delete fOutputHitsEvent;
66}
67
68void TRestDetectorGarfieldDriftProcess::LoadDefaultConfig() {
69 SetName("garfieldDriftProcess-Default");
70 SetTitle("Default config");
71
72 cout << "TRestDetectorGarfieldDriftProcess metadata not found. Loading default values" << endl;
73
74 fDriftPotential = 1000;
75 fGasPressure = 10;
76}
77
78void TRestDetectorGarfieldDriftProcess::LoadConfig(const string& configFilename, const string& name) {
79 if (LoadConfigFromFile(configFilename, name)) LoadDefaultConfig();
80
81 if (fDriftPotential == PARAMETER_NOT_FOUND_DBL) {
82 fDriftPotential = 1000; // V
83 }
84}
85
86#endif
87
89 SetSectionName(ClassName());
90 SetLibraryVersion(LIBRARY_VERSION);
91
92 fRandom = new TRandom3(0);
93 fInputHitsEvent = nullptr;
94 fOutputHitsEvent = new TRestDetectorHitsEvent();
95
96#if defined USE_Garfield
97 fReadout = nullptr;
98 fGas = nullptr;
99 fGeometry = nullptr;
100 fPEReduction = 1.;
101 fStopDistance = 2; // default distance from readout to stop drift set to 2mm
102#endif
103}
104
105#if defined USE_Garfield
107 // Function to be executed once at the beginning of process
108 // (before starting the process of the events)
109
110 // Start by calling the InitProcess function of the abstract class.
111 // Comment this if you don't want it.
112 // TRestEventProcess::InitProcess();
113
114 // Getting gas data
115 fGas = GetMetadata<TRestDetectorGas>();
116 if (fGas != nullptr) {
117 if (fGasPressure <= 0)
118 fGasPressure = fGas->GetPressure();
119 else
120 fGas->SetPressure(fGasPressure);
121
122 } else {
123 cout << "REST_WARNING. No TRestDetectorGas found in TRestRun." << endl;
124 }
125
126 // Getting readout data
127 fReadout = GetMetadata<TRestDetectorReadout>();
128 if (fReadout == nullptr) {
129 cout << "REST ERROR: Readout has not been initialized" << endl;
130 exit(-1);
131 }
132
133 // Getting geometry data
134 if (fGDML_Filename == "") {
135 cout << "REST ERROR no geometry GDML file name parameter (gdml_file) given "
136 "in TRestDetectorGarfieldDriftProcess "
137 << endl;
138 exit(-1);
139 }
140
141 TGeoVolume* geovol = nullptr;
142 fGeometry = new TRestDetectorGeometry();
143 fGeometry->InitGfGeometry();
144
145 geovol = TGDMLParse::StartGDML(fGDML_Filename);
146 if (!geovol) {
147 cout << "REST ERROR when TRestDetectorGeometry read GDML file " << fGDML_Filename << endl;
148 exit(-1);
149 }
150
151 fGeometry->SetTopVolume(geovol);
152 fGeometry->CloseGeometry();
153 fGeometry->DefaultColors();
154 // fGeometry->UpdateElements();
155
156 cout << "TRestDetectorGarfieldDriftProcess GetVerboseLevel : "
157 << static_cast<int>(this->GetVerboseLevel()) << endl;
158
159 // analyze GDML geometry to find major elements (gas volume, electrodes,
160 // readout)
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();
167 }
168 cout << "****** itnode " << itnode->GetName() << endl;
169 itnode->PrintCandidates();
170
171 TGeoVolume* itvol = itnode->GetVolume();
173 cout << " * * itvolume " << itvol->GetName() << endl;
174 itvol->Print();
175 }
176 cout << " * * itvolume " << itvol->GetName() << endl;
177 itvol->Print();
178
179 TGeoMedium* itmed = itvol->GetMedium();
181 cout << " * * itmed " << itmed->GetName() << endl;
182
183 // gas volume
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();
190 }
191 cout << " -> gas volume SetMedium itmed " << itmed->GetName() << " fGas "
192 << fGas->GetGasMedium()->GetName() << endl;
193 fGas->GetGasMedium()->PrintGas();
194 }
195
196 // drift electrode (it should be called a cathode, not anode, no ?)
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;
203 }
204
205 // micromegas readout electrode
206 if ((strncmp(itvol->GetName(), "micromegasVol", 13) == 0)) {
207 fGeometry->AddReadoutElecNode(itnode);
209 cout << " -> readout volume " << endl;
210 cout << " -> readout volume " << endl;
211 }
212 }
213
214 // For the moment we set only constant electric field in Garfield
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;
219
220 TGeoNode* readoutnode = fGeometry->GetReadoutElecNode(ii);
221 cout << "ii " << ii << " readoutnode " << readoutnode << " visible " << readoutnode->IsVisible()
222 << endl;
223 if (!readoutnode) continue;
224 cout << "readoutnode visible " << readoutnode->IsVisible() << endl;
225 readoutnode->ls();
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]
232 : 0); // converted to mm
233 for (int jj = 0; jj < fReadout->GetNumberOfReadoutPlanes(); jj++) {
234 TRestDetectorReadoutPlane* readoutplane = fReadout->GetReadoutPlane(jj);
235 planeZpos = readoutplane->GetPosition().Z();
236 cout << " jj " << jj << " matrixZpos " << matrixZpos << " planeZpos " << planeZpos << endl;
237 if (fabs(planeZpos - matrixZpos) > 1)
238 continue; // we search for fReadout entry at same Z position
239 rdPlaneFound = true;
240 planeZvec = readoutplane->GetNormal().z();
241 cout << " planeZvec " << planeZvec << endl;
242 }
243
244 if (rdPlaneFound) {
245 double field = fDriftPotential /
246 (fGeometry->GetDriftElecNode()->GetMatrix()->GetTranslation()[2] - planeZpos);
247 ComponentConstant* cmp = new ComponentConstant();
248 // cmp->SetElectricField(0, 0, fElectricField); // assuming V/cm
249 cmp->SetElectricField(0, 0, field); // assuming V/cm
250 cmp->SetPotential(0, 0, planeZpos / 10., 0); // potential at readout level set to 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;
256 }
257 }
258
259 // create Sensor corresponding to readout plane
260 // only the first component used, so 1 sensor only for the moment
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);
266
267 double xmin = std::numeric_limits<Double_t>::min(), xmax = std::numeric_limits<Double_t>::max(),
268 ymin = std::numeric_limits<Double_t>::min(), ymax = std::numeric_limits<Double_t>::max();
269 TGeoShape* readoutshape = fGeometry->GetReadoutElecNode(0)->GetVolume()->GetShape();
270 TGeoMatrix* readoutmatrix = fGeometry->GetReadoutElecNode(0)->GetMatrix();
271 double xmid = 10. * readoutmatrix->GetTranslation()[0],
272 ymid = 10. * readoutmatrix->GetTranslation()[1];
273 if (readoutshape->InheritsFrom("TGeoBBox")) {
274 // keep drifting electrons 10cm around the readout volume in x and y
275 TGeoBBox* readoutbox = (TGeoBBox*)readoutshape;
276 xmin = xmid - 10. * readoutbox->GetDX() - 100;
277 xmax = xmid + 10. * readoutbox->GetDX() + 100;
278 ymin = ymid - 10. * readoutbox->GetDY() - 100;
279 ymax = ymid + 10. * readoutbox->GetDY() + 100;
280 }
281 double driftelecZpos = fGeometry->GetDriftElecNode()
282 ->GetMatrix()
283 ->GetTranslation()[2]; // hope that all these objects are really there...
284 // drift area defined from bouding box of readout shape
285 fGfSensor->SetArea(xmin / 10., ymin / 10., planeZpos / 10. + fStopDistance, xmax / 10., ymax / 10.,
286 driftelecZpos); // drift stops fStopDistance above the readout plane
287 fGeometry->AddGfSensor(fGfSensor);
288 printf(" GfSensor area x- %lf y- %lf z- %lf x+ %lf y+ %lf z+ %lf\n", xmin / 10., ymin / 10.,
289 planeZpos / 10. + 0.2, xmax / 10., ymax / 10., driftelecZpos);
290
291 fGfDriftMethod = new DRIFT_METHOD();
292 // fGfDriftMethod->EnableDebugging();
293 // fGfDriftMethod->EnableVerbose();
294 fGfDriftMethod->SetSensor(fGfSensor);
295 // fGfDriftMethod->SetTimeSteps(0.2);
296 fGfDriftMethod->SetCollisionSteps(100);
297 }
298
299 if (!fGfDriftMethod || !fGfSensor) {
300 cout << "REST ERROR Garfield not fully initialized in "
301 "TRestDetectorGarfieldDriftProcess::InitProcess "
302 << endl;
303 cout << " fGfDriftMethod " << fGfDriftMethod << " fGfSensor " << fGfSensor << " GetNbOfGfComponent "
304 << fGeometry->GetNbOfGfComponent() << endl;
305 cout << " matrixZpos " << matrixZpos << " planeZpos " << planeZpos << " planeZvec " << planeZvec
306 << endl;
307 exit(-1);
308 }
309
310 // for (register double ix=-3000; ix<3000; ix+=10) {
311 // printf(" ix %lf medium %p is_inside %d", ix,
312 // fGeometry->GetGfGeometry()->GetMedium(ix/10.,5,0),
313 // fGeometry->GetGfGeometry()->IsInside(ix/10.,5,0));
314 // cout << " moduleId " <<
315 // fReadout->GetReadoutPlane(0)->isInsideDriftVolume( ix, 50, 0 );
316 // fGeometry->SetCurrentPoint(ix/10., 5,0);
317 // if (fGeometry->IsOutside()) cout << " is outside " << endl; else cout <<
318 // " is inside" << endl; TGeoNode* cnode =fGeometry ->GetCurrentNode();
319 // TGeoNode* cnode2= fGeometry->FindNode(ix/10., 5,0);
320 // std::string name(cnode->GetMedium()->GetMaterial()->GetName());
321 // register Medium* tmed = fGeometry->GetGfMedium(ix/10., 5,0);
322 // cout << " node " << cnode->GetName() << " material " << name;
323 // cout << " node2 " << cnode2->GetName() << " material " <<
324 // cnode2->GetMedium()->GetMaterial()->GetName(); if (tmed) cout << " medium
325 // " << tmed->GetName(); cout << endl;
326 // }
327}
328
329//------------------------------------------------------------------------------
330
331Int_t TRestDetectorGarfieldDriftProcess::FindModule(Int_t readoutPlane, Double_t x, Double_t y) {
332 // TODO verify this
333 TRestDetectorReadoutPlane* plane = fReadout->GetReadoutPlane(readoutPlane);
334 for (size_t md = 0; md < plane->GetNumberOfModules(); md++) {
335 if ((*plane)[md].IsInside({x, y})) {
336 return md;
337 }
338 }
339
340 return -1;
341}
342
343#endif
344
346#if defined USE_Garfield
347 fInputHitsEvent = (TRestDetectorHitsEvent*)inputEvent;
348
349 double x, y, z, energy;
350 double xi, yi, zi, ti, xf, yf, zf, tf, energyf;
351 int status;
352
354 cout << "Number of hits : " << fInputHitsEvent->GetNumberOfHits() << endl;
355 cout << "--------------------------" << endl;
356 fInputHitsEvent->PrintEvent(20);
357 }
358
359 for (unsigned int hit = 0; hit < fInputHitsEvent->GetNumberOfHits(); hit++) {
360 x = fInputHitsEvent->GetX(hit);
361 y = fInputHitsEvent->GetY(hit);
362 z = fInputHitsEvent->GetZ(hit);
363 energy = fInputHitsEvent->GetEnergy(hit);
364 // if (hit%10==0) printf("hit: x %lf y %lf z %lf Energy %lf\n", x,
365 // y, z, energy); cout << "hit: x " << x << " y " << y <<
366 // " z " << z << " Energy " << energy << endl;
367
368 // assuming energy in keV, W factor in eV
369 Medium* tmed = fGeometry->GetGfMedium(x, y, z);
370 if (!tmed) {
371 cout << "no medium found at x " << x << " y " << y << " z " << z << endl;
372 continue;
373 // } else {
374 // cout << " medium " << tmed->GetName() << " W " <<
375 // tmed->GetW() << " at x " << x << " y " << y << " z " << z <<
376 // endl;
377 }
378 double Wfactor = tmed->GetW();
379 double nb_electrons_real = energy * 1000. / Wfactor / fPEReduction;
380 unsigned int nb_electrons =
381 fRandom->Poisson(nb_electrons_real); // we need an int number of electrons to drift...
382 // if (nb_electrons>0) printf("hit: x %lf y %lf z %lf Energy %lf\n",
383 // x, y, z, energy); if (nb_electrons>0) cout << " : energy " <<
384 // energy << " nb_electrons_real " << nb_electrons_real << "
385 // nb_electrons " << nb_electrons << endl;
386
387 // lets drift all those electrons
388 for (unsigned int iel = 0; iel < nb_electrons; iel++) {
389 // drift an electron
390 fGfDriftMethod->DriftElectron(x / 10., y / 10., z / 10., 0);
391 // fGfDriftMethod->GetEndPoint(xf, yf, zf, tf, status);
392 if (fGfDriftMethod->GetNumberOfElectronEndpoints() == 0) {
393 cout << " TRestDetectorGarfieldDriftProcess::ProcessEvent: no electron end "
394 "point ???"
395 << endl;
396 continue;
397 }
398
399 // create an hit for the drifted electron
400 fGfDriftMethod->GetElectronEndpoint(0, xi, yi, zi, ti, xf, yf, zf, tf,
401 status); // works with AvalancheMC
402 // xi *= 10; yi *= 10; zi *= 10; // cm to mm
403 xf *= cmTomm;
404 yf *= cmTomm;
405 zf *= cmTomm; // cm to mm
406
407 tf /= 1000.; // ns to us
408 energyf = fPEReduction * Wfactor / 1000.; // back to keV
409 fOutputHitsEvent->AddHit(xf, yf, zf, energyf, tf);
410
411 // if (nb_electrons>0) cout << " : Nendpoint " <<
412 // fGfDriftMethod->GetNumberOfElectronEndpoints() << " xi " <<
413 // xi << " yi " << yi << " zi " << zi << " ti " << ti
414 // << " xf " << xf << " yf " << yf << " zf " <<
415 // zf << " tf " << tf << " status " << status ;
416 }
417 }
418
419 if (fOutputHitsEvent->GetNumberOfHits() == 0) {
420 return nullptr;
421 }
422
423 // fSignalEvent->PrintEvent();
424
426 cout << "TRestDetectorElectronDiffusionProcess. Hits added : " << fOutputHitsEvent->GetNumberOfHits()
427 << endl;
428 cout << "TRestDetectorElectronDiffusionProcess. Hits total energy : "
429 << fOutputHitsEvent->GetTotalEnergy() << endl;
430 cout << " fTimedHitsEvent " << fOutputHitsEvent << " class " << fOutputHitsEvent->ClassName() << endl;
431 fOutputHitsEvent->PrintEvent(20);
432 }
433
434 return fOutputHitsEvent;
435#else
436 fInputHitsEvent = (TRestDetectorHitsEvent*)inputEvent;
437 fOutputHitsEvent = fInputHitsEvent;
438 RESTDebug << "nullptr process" << RESTendl;
439 return inputEvent;
440
441#endif
442}
443
444#if defined USE_Garfield
445
447 // Function to be executed once at the end of the process
448 // (after all events have been processed)
449
450 // Start by calling the EndProcess function of the abstract class.
451 // Comment this if you don't want it.
452 // TRestEventProcess::EndProcess();
453}
454
456 fGasPressure = StringToDouble(GetParameter("gasPressure", "-1"));
457 // fElectricField = GetDblParameterWithUnits( "electricField" );
458 fDriftPotential = GetDblParameterWithUnits("driftPotential");
459 fPEReduction = StringToDouble(GetParameter("primaryElectronReduction", "1"));
460 fStopDistance = StringToDouble(GetParameter("stopDistance", "2"));
461
462 // TODO : Still units must be implemented for velocity quantities
463
464 fGDML_Filename = GetParameter("geometryPath", "") + GetParameter("gdml_file", "");
465}
466
467#endif
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.
Definition: TRestEvent.h:38
endl_t RESTendl
Termination flag object for TRestStringOutput.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
@ 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.