17#include "TRestDetectorGarfieldDriftProcess.h"
24#if defined USE_Garfield
25using namespace Garfield;
26#include <ComponentConstant.hh>
28#include "TGDMLParse.h"
35const double cmTomm = 10.;
37#if defined USE_Garfield
39TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess() : fRandom(0), fGfSensor(0) {
50TRestDetectorGarfieldDriftProcess::TRestDetectorGarfieldDriftProcess(
const char* configFilename)
51 : fRandom(0), fGfSensor(0) {
54 if (LoadConfigFromFile(configFilename) == -1) LoadDefaultConfig();
62TRestDetectorGarfieldDriftProcess::~TRestDetectorGarfieldDriftProcess() {
65 delete fOutputHitsEvent;
68void TRestDetectorGarfieldDriftProcess::LoadDefaultConfig() {
69 SetName(
"garfieldDriftProcess-Default");
70 SetTitle(
"Default config");
72 cout <<
"TRestDetectorGarfieldDriftProcess metadata not found. Loading default values" << endl;
74 fDriftPotential = 1000;
78void 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 = 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")) {
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;
281 double driftelecZpos = fGeometry->GetDriftElecNode()
283 ->GetTranslation()[2];
285 fGfSensor->SetArea(xmin / 10., ymin / 10., planeZpos / 10. + fStopDistance, xmax / 10., ymax / 10.,
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);
291 fGfDriftMethod =
new DRIFT_METHOD();
294 fGfDriftMethod->SetSensor(fGfSensor);
296 fGfDriftMethod->SetCollisionSteps(100);
299 if (!fGfDriftMethod || !fGfSensor) {
300 cout <<
"REST ERROR Garfield not fully initialized in "
301 "TRestDetectorGarfieldDriftProcess::InitProcess "
303 cout <<
" fGfDriftMethod " << fGfDriftMethod <<
" fGfSensor " << fGfSensor <<
" GetNbOfGfComponent "
304 << fGeometry->GetNbOfGfComponent() << endl;
305 cout <<
" matrixZpos " << matrixZpos <<
" planeZpos " << planeZpos <<
" planeZvec " << planeZvec
331Int_t TRestDetectorGarfieldDriftProcess::FindModule(Int_t readoutPlane, Double_t x, Double_t y) {
335 if ((*plane)[md].IsInside({x, y})) {
346#if defined USE_Garfield
349 double x, y, z, energy;
350 double xi, yi, zi, ti, xf, yf, zf, tf, energyf;
354 cout <<
"Number of hits : " << fInputHitsEvent->GetNumberOfHits() << endl;
355 cout <<
"--------------------------" << endl;
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);
369 Medium* tmed = fGeometry->GetGfMedium(x, y, z);
371 cout <<
"no medium found at x " << x <<
" y " << y <<
" z " << z << endl;
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);
388 for (
unsigned int iel = 0; iel < nb_electrons; iel++) {
390 fGfDriftMethod->DriftElectron(x / 10., y / 10., z / 10., 0);
392 if (fGfDriftMethod->GetNumberOfElectronEndpoints() == 0) {
393 cout <<
" TRestDetectorGarfieldDriftProcess::ProcessEvent: no electron end "
400 fGfDriftMethod->GetElectronEndpoint(0, xi, yi, zi, ti, xf, yf, zf, tf,
408 energyf = fPEReduction * Wfactor / 1000.;
409 fOutputHitsEvent->
AddHit(xf, yf, zf, energyf, tf);
419 if (fOutputHitsEvent->GetNumberOfHits() == 0) {
426 cout <<
"TRestDetectorElectronDiffusionProcess. Hits added : " << fOutputHitsEvent->GetNumberOfHits()
428 cout <<
"TRestDetectorElectronDiffusionProcess. Hits total energy : "
429 << fOutputHitsEvent->GetTotalEnergy() << endl;
430 cout <<
" fTimedHitsEvent " << fOutputHitsEvent <<
" class " << fOutputHitsEvent->ClassName() << endl;
434 return fOutputHitsEvent;
437 fOutputHitsEvent = fInputHitsEvent;
438 RESTDebug <<
"nullptr process" <<
RESTendl;
444#if defined USE_Garfield
458 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.