REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorReadoutPlane.cxx
1
33
34#include "TRestDetectorReadoutPlane.h"
35
36using namespace std;
37
43
48
53 Int_t nChannels = 0;
54 for (size_t md = 0; md < GetNumberOfModules(); md++) {
55 nChannels += fReadoutModules[md].GetNumberOfChannels();
56 }
57 return nChannels;
58}
59
64void TRestDetectorReadoutPlane::SetNormal(const TVector3& normal) {
65 fNormal = normal.Unit();
66 // prevent user from declaring the zero vector as normal
67 if (TMath::Abs(fNormal.Mag2() - 1.0) > 1E-6) {
68 // only the zero vector will have a magnitude different from 1.0 after normalization
69 RESTError << "TRestDetectorReadoutPlane::SetNormal : normal vector cannot be zero." << RESTendl;
70 exit(1);
71 }
72 UpdateAxes();
73}
74
79 if (height < 0) {
80 fHeight = 0;
81 RESTError << "TRestDetectorReadoutPlane::SetHeight : height cannot be negative." << RESTendl;
82 exit(1);
83 } else {
84 fHeight = height;
85 }
86}
87
92 for (size_t md = 0; md < GetNumberOfModules(); md++) {
93 if (fReadoutModules[md].GetModuleID() == modID) {
94 return &fReadoutModules[md];
95 }
96 }
97
98 cout << "REST ERROR (GetReadoutModuleByID) : Module ID : " << modID << " was not found" << endl;
99 return nullptr;
100}
101
113Double_t TRestDetectorReadoutPlane::GetX(Int_t modID, Int_t chID) {
115
116 TRestDetectorReadoutChannel* rChannel = rModule->GetChannel(chID);
117
118 Double_t x = numeric_limits<Double_t>::quiet_NaN();
119
120 if (rChannel->GetNumberOfPixels() == 1) {
121 Double_t sX = rChannel->GetPixel(0)->GetSize().X();
122 Double_t sY = rChannel->GetPixel(0)->GetSize().Y();
123
124 if (sX > 2 * sY) return x;
125
126 x = rModule->GetPixelCenter(chID, 0).X();
127
128 return x;
129 }
130
131 if (rChannel->GetNumberOfPixels() > 1) {
132 Int_t nPix = rChannel->GetNumberOfPixels();
133
134 // We check the origin of consecutive pixels to check if it goes X or Y
135 // direction. Perhaps more complex readouts need some changes here
136 Double_t x1 = rChannel->GetPixel(0)->GetCenter().X();
137 Double_t x2 = rChannel->GetPixel(nPix - 1)->GetCenter().X();
138
139 Double_t y1 = rChannel->GetPixel(0)->GetCenter().Y();
140 Double_t y2 = rChannel->GetPixel(nPix - 1)->GetCenter().Y();
141
142 /*
143 cout << "Pixel 0. X : " << x1 << " Y : " << y1 endl;
144 cout << "Pixel N. X : " << x2 << " Y : " << y2 endl;
145 */
146
147 Double_t deltaX = abs(x2 - x1);
148 Double_t deltaY = abs(y2 - y1);
149
150 Int_t rotation = (Int_t)(std::round(rModule->GetRotation() * units("degrees")));
151 if (rotation % 90 == 0) {
152 if (rotation / 90 % 2 == 0) // rotation is 0, 180, 360...
153 {
154 if (deltaY > deltaX) x = rModule->GetPixelCenter(chID, 0).X();
155 } else // rotation is 90, 270... We need to invert x and y
156 {
157 if (deltaY < deltaX) x = rModule->GetPixelCenter(chID, 0).X();
158 }
159 } else {
160 // we choose to output x only when deltaY > deltaX under non-90 deg rotation
161 // otherwise it is a y channel and should return nan
162 if (deltaY > deltaX) x = rModule->GetPixelCenter(chID, 0).X();
163 }
164 }
165
166 return x;
167}
168
180Double_t TRestDetectorReadoutPlane::GetY(Int_t modID, Int_t chID) {
182
183 TRestDetectorReadoutChannel* rChannel = rModule->GetChannel(chID);
184
185 Double_t y = numeric_limits<Double_t>::quiet_NaN();
186
187 if (rChannel->GetNumberOfPixels() == 1) {
188 Double_t sX = rChannel->GetPixel(0)->GetSize().X();
189 Double_t sY = rChannel->GetPixel(0)->GetSize().Y();
190
191 if (sY > 2 * sX) return y;
192
193 y = rModule->GetPixelCenter(chID, 0).Y();
194
195 return y;
196 }
197
198 if (rChannel->GetNumberOfPixels() > 1) {
199 Int_t nPix = rChannel->GetNumberOfPixels();
200
201 // We check the origin of consecutive pixels to check if it goes X or Y
202 // direction. Perhaps more complex readouts need some changes here
203 Double_t x1 = rChannel->GetPixel(0)->GetCenter().X();
204 Double_t x2 = rChannel->GetPixel(nPix - 1)->GetCenter().X();
205
206 Double_t y1 = rChannel->GetPixel(0)->GetCenter().Y();
207 Double_t y2 = rChannel->GetPixel(nPix - 1)->GetCenter().Y();
208
209 /*
210 cout << "Pix id : " << rChannel->GetPixel(0)->GetID() << " X1 : " << x1
211 << endl; cout << "Pix id : " << rChannel->GetPixel(1)->GetID() << " X2 :
212 " << x2 << endl; cout << "Pix id : " << rChannel->GetPixel(0)->GetID() <<
213 " Y1 : " << y1 << endl; cout << "Pix id : " <<
214 rChannel->GetPixel(1)->GetID() << " Y2 : " << y2 << endl;
215 */
216
217 Double_t deltaX = abs(x2 - x1);
218 Double_t deltaY = abs(y2 - y1);
219
220 Int_t rotation = (Int_t)std::round(rModule->GetRotation() * units("degrees"));
221 if (rotation % 90 == 0) {
222 if (rotation / 90 % 2 == 0) // rotation is 0, 180, 360...
223 {
224 if (deltaY < deltaX) y = rModule->GetPixelCenter(chID, 0).Y();
225 } else // rotation is 90, 270... We need to invert x and y
226 {
227 if (deltaY > deltaX) y = rModule->GetPixelCenter(chID, 0).Y();
228 }
229 } else {
230 // we choose to output y only when deltaY < deltaX under non-90 deg rotation
231 // otherwise it is a x channel and should return nan
232 if (deltaY < deltaX) y = rModule->GetPixelCenter(chID, 0).Y();
233 }
234 }
235
236 return y;
237}
238
246Int_t TRestDetectorReadoutPlane::FindChannel(Int_t module, const TVector2& position) {
247 const auto relativePosition = position - TVector2{fPosition.X(), fPosition.Y()};
248
249 // TODO : check first if (modX,modY) is inside the module.
250 // If not return error.
251 // FindChannel will take a long time to search for the channel if it is not
252 // there. It will be faster
253
254 return fReadoutModules[module].FindChannel(relativePosition);
255}
256
261Double_t TRestDetectorReadoutPlane::GetDistanceTo(const TVector3& position) const {
262 const TVector3 diff = position - fPosition;
263 return diff.Dot(fNormal);
264}
265
274 TVector3 pos = TVector3(0, 0, z);
275
276 return isZInsideDriftVolume(pos);
277}
278
287 for (size_t m = 0; m < GetNumberOfModules(); m++)
288 if (fReadoutModules[m].IsDaqIDInside(daqId)) return true;
289
290 return false;
291}
292
302Int_t TRestDetectorReadoutPlane::isZInsideDriftVolume(const TVector3& position) {
303 TVector3 posNew = TVector3(position.X() - fPosition.X(), position.Y() - fPosition.Y(), position.Z());
304
305 Double_t distance = GetDistanceTo(posNew);
306
307 if (distance > 0 && distance < fHeight) {
308 return 1;
309 }
310
311 return 0;
312}
313
325Int_t TRestDetectorReadoutPlane::GetModuleIDFromPosition(const TVector3& position) const {
326 Double_t distance = GetDistanceTo(position);
327 if (distance >= 0 && distance <= fHeight) {
328 const TVector2 positionInPlane = GetPositionInPlane(position);
329 for (size_t m = 0; m < GetNumberOfModules(); m++) {
330 auto& module = fReadoutModules[m];
331 if (module.IsInside(positionInPlane)) {
332 return module.GetModuleID();
333 }
334 }
335 }
336
337 return -1;
338}
339
344void TRestDetectorReadoutPlane::Print(Int_t DetailLevel) {
345 if (DetailLevel >= 0) {
346 RESTMetadata << "-- Readout plane : " << GetID() << RESTendl;
347 RESTMetadata << "----------------------------------------------------------------" << RESTendl;
348 RESTMetadata << "-- Position : X = " << fPosition.X() << " mm, "
349 << " Y : " << fPosition.Y() << " mm, Z : " << fPosition.Z() << " mm" << RESTendl;
350 RESTMetadata << "-- Normal vector : X = " << fNormal.X() << " mm, "
351 << " Y : " << fNormal.Y() << " mm, Z : " << fNormal.Z() << " mm" << RESTendl;
352 RESTMetadata << "-- X-axis vector : X = " << fAxisX.X() << " mm, "
353 << " Y : " << fAxisX.Y() << " mm, Z : " << fAxisX.Z() << " mm" << RESTendl;
354 RESTMetadata << "-- Y-axis vector : Y = " << fAxisY.X() << " mm, Y : " << fAxisY.Y()
355 << " mm, Z : " << fAxisY.Z() << " mm" << RESTendl;
356 RESTMetadata << "-- Cathode Position : X = " << GetCathodePosition().X() << " mm, "
357 << " Y : " << GetCathodePosition().Y() << " mm, Z : " << GetCathodePosition().Z()
358 << " mm" << RESTendl;
359 RESTMetadata << "-- Height : " << fHeight << " mm" << RESTendl;
360 RESTMetadata << "-- Charge collection : " << fChargeCollection << RESTendl;
361 RESTMetadata << "-- Total modules : " << GetNumberOfModules() << RESTendl;
362 RESTMetadata << "-- Total channels : " << GetNumberOfChannels() << RESTendl;
363 RESTMetadata << "----------------------------------------------------------------" << RESTendl;
364
365 for (size_t i = 0; i < GetNumberOfModules(); i++) {
366 fReadoutModules[i].Print(DetailLevel - 1);
367 }
368 }
369}
370
375
381 Double_t x[4];
382 Double_t y[4];
383
384 double xmin, xmax, ymin, ymax;
385
386 GetBoundaries(xmin, xmax, ymin, ymax);
387
388 TH2Poly* readoutHistogram = new TH2Poly("ReadoutHistogram", "ReadoutHistogram", xmin, xmax, ymin, ymax);
389
390 for (size_t mdID = 0; mdID < this->GetNumberOfModules(); mdID++) {
392
393 int nChannels = module->GetNumberOfChannels();
394
395 for (int ch = 0; ch < nChannels; ch++) {
396 TRestDetectorReadoutChannel* channel = module->GetChannel(ch);
397 Int_t nPixels = channel->GetNumberOfPixels();
398
399 for (int px = 0; px < nPixels; px++) {
400 for (int v = 0; v < 4; v++) {
401 x[v] = module->GetPixelVertex(ch, px, v).X();
402 y[v] = module->GetPixelVertex(ch, px, v).Y();
403 }
404
405 readoutHistogram->AddBin(4, x, y);
406 }
407 }
408 }
409
410 readoutHistogram->SetStats(false);
411
412 return readoutHistogram;
413}
414
418void TRestDetectorReadoutPlane::GetBoundaries(double& xmin, double& xmax, double& ymin, double& ymax) {
419 Double_t x[4];
420 Double_t y[4];
421
422 xmin = std::numeric_limits<Double_t>::max(), xmax = std::numeric_limits<Double_t>::min(),
423 ymin = std::numeric_limits<Double_t>::max(), ymax = std::numeric_limits<Double_t>::min();
424
425 for (size_t mdID = 0; mdID < this->GetNumberOfModules(); mdID++) {
427
428 for (int v = 0; v < 4; v++) {
429 x[v] = module->GetVertex(v).X();
430 y[v] = module->GetVertex(v).Y();
431
432 if (x[v] < xmin) xmin = x[v];
433 if (y[v] < ymin) ymin = y[v];
434 if (x[v] > xmax) xmax = x[v];
435 if (y[v] > ymax) ymax = y[v];
436 }
437 }
438}
439
440void TRestDetectorReadoutPlane::UpdateAxes() { // idempotent
441 const TVector3 zUnit = {0, 0, 1};
442 fAxisX = {1, 0, 0};
443 fAxisY = {0, 1, 0};
444
445 constexpr double tolerance = 1E-6;
446
447 // Check if fNormal is different from (0,0,1)
448 if ((fNormal - zUnit).Mag2() < tolerance) {
449 // do nothing
450 } else if ((fNormal + zUnit).Mag2() < tolerance) {
451 // normal vector is opposite to (0,0,1), we must also flip the axes
452 fAxisX = {0, -1, 0};
453 fAxisY = {-1, 0, 0};
454 } else {
455 // Calculate the rotation axis by taking the cross product between the original normal and fNormal
456 TVector3 rotationAxis = zUnit.Cross(fNormal);
457
458 // Calculate the rotation angle using the dot product between the original normal and fNormal
459 double rotationAngle = acos(zUnit.Dot(fNormal) / (zUnit.Mag() * fNormal.Mag()));
460
461 // Rotate the axes around the rotation axis by the rotation angle
462 fAxisX.Rotate(rotationAngle, rotationAxis);
463 fAxisY.Rotate(rotationAngle, rotationAxis);
464 }
465
466 // rotate around normal by rotation angle (angle in radians)
467 fAxisX.Rotate(fRotation, fNormal);
468 fAxisY.Rotate(fRotation, fNormal);
469
470 // verify that fNormal, fAxisX and fAxisY are orthogonal and unitary
471 if (TMath::Abs(fNormal.Mag2() - 1.0) > tolerance || TMath::Abs(fAxisX.Mag2() - 1.0) > tolerance ||
472 TMath::Abs(fAxisY.Mag2() - 1.0) > tolerance) {
473 RESTError << "TRestDetectorReadoutPlane::UpdateAxes() : "
474 << "The normal vector, the X-axis vector and the Y-axis vector must be unitary."
475 << RESTendl;
476 exit(1);
477 }
478 if (TMath::Abs(fNormal.Dot(fAxisX)) > tolerance || TMath::Abs(fNormal.Dot(fAxisY)) > tolerance ||
479 TMath::Abs(fAxisX.Dot(fAxisY)) > tolerance) {
480 RESTError << "TRestDetectorReadoutPlane::UpdateAxes() : "
481 << "The normal vector, the X-axis vector and the Y-axis vector must be orthogonal."
482 << RESTendl;
483 exit(1);
484 }
485 // verify that the correct order of axes is being used: X cross Y = normal (and not - normal)
486 if ((fAxisX.Cross(fAxisY) - fNormal).Mag2() > tolerance) {
487 RESTError
488 << "TRestDetectorReadoutPlane::UpdateAxes() : "
489 << "The normal vector is not the cross product between the X-axis vector and the Y-axis vector."
490 << RESTendl;
491 exit(1);
492 }
493}
494
495void TRestDetectorReadoutPlane::SetRotation(Double_t radians) {
496 // sets fRotation modulo 2pi
497 fRotation = TVector2::Phi_0_2pi(radians);
498 UpdateAxes();
499}
500
501TVector2 TRestDetectorReadoutPlane::GetPositionInPlane(const TVector3& point) const {
502 // Given a point in space, returns the position of the point in the plane using the plane's axes
503 // The position is returned in the plane's local coordinates (in mm)
504 return {fAxisX.Dot(point - fPosition),
505 fAxisY.Dot(point - fPosition)}; // dot product between vectors is the projection
506}
507
508TVector3 TRestDetectorReadoutPlane::GetPositionInWorld(const TVector2& point, Double_t height) const {
509 return fPosition + point.X() * fAxisX + point.Y() * fAxisY + height * fNormal;
510}
511
512Double_t TRestDetectorReadoutPlane::GetDistanceToPlane(const TVector3& point) const {
513 return (point - fPosition).Dot(fNormal);
514}
515
516void TRestDetectorReadoutPlane::SetAxisX(const TVector3& axis) {
517 const TVector3 axisInPlane = (axis - axis.Dot(fNormal) * fNormal).Unit();
518 if (axisInPlane.Mag2() < 1E-6) {
519 RESTError << "TRestDetectorReadoutPlane::SetAxisX() : "
520 << "The X-axis vector must not be parallel to the normal vector." << RESTendl;
521 exit(1);
522 }
523
524 // compute the angle between the new axis and the old axis
525 const Double_t angle = fAxisX.Angle(axisInPlane);
526
527 SetRotation(fRotation - angle);
528}
529
530bool TRestDetectorReadoutPlane::IsInside(const TVector3& point) const {
531 const double distance = GetDistanceToPlane(point);
532 if (distance < 0 || distance > fHeight) {
533 // point is outside the volume defined by the plane
534 return false;
535 }
536 const TVector2 pointInPlane = GetPositionInPlane(point);
537 for (size_t moduleIndex = 0; moduleIndex < GetNumberOfModules(); moduleIndex++) {
538 if (fReadoutModules[moduleIndex].IsInside(pointInPlane)) {
539 return true;
540 }
541 }
542 return false;
543}
544
546 cout << "Adding module" << endl;
547 fReadoutModules.emplace_back(module);
548 // if the module has no name or no type, add the one from the plane
549
550 auto& lastModule = fReadoutModules.back();
551 if (lastModule.GetName().empty()) {
552 lastModule.SetName(fName);
553 }
554 if (lastModule.GetType().empty()) {
555 lastModule.SetType(fType);
556 }
557
558 for (size_t channelIndex = 0; channelIndex < lastModule.GetNumberOfChannels(); channelIndex++) {
559 TRestDetectorReadoutChannel* channel = lastModule.GetChannel(channelIndex);
560 if (channel->GetName().empty()) {
561 channel->SetName(lastModule.GetName());
562 }
563 if (channel->GetType().empty()) {
564 channel->SetType(lastModule.GetType());
565 }
566 }
567}
std::string GetType() const
Returns the channel type.
void SetName(const std::string &name)
Sets the channel name.
Int_t GetNumberOfPixels()
Returns the total number of pixels inside the readout channel.
void SetType(const std::string &type)
Sets the channel type.
std::string GetName() const
Returns the channel name.
TRestDetectorReadoutPixel * GetPixel(int n)
Returns a pointer to the pixel n by index.
TVector2 GetPixelCenter(Int_t channel, Int_t pixel)
Returns the center pixel position for a given channel and pixel indexes.
TRestDetectorReadoutChannel * GetChannel(size_t n)
Returns a pointer to a readout channel by index.
Double_t GetRotation() const
Returns the module rotation in degrees.
TVector2 GetPixelVertex(Int_t channel, Int_t pixel, Int_t vertex)
Returns any of the pixel vertex position for a given channel and pixel indexes.
TVector2 GetVertex(int n) const
Returns the coordinates of the specified vertex index n. The physical coordinates relative to the rea...
size_t GetNumberOfChannels() const
Returns the total number of channels defined inside the module.
TVector2 GetCenter() const
Returns the center TVector2 position of the pixel.
TVector2 GetSize()
Returns a TVector2 with the pixel size.
Bool_t isDaqIDInside(Int_t daqId)
This method determines if the daqId given is associated to any of the readout readout channels in any...
Int_t GetNumberOfChannels()
Returns the total number of channels in the readout plane.
Int_t isZInsideDriftVolume(Double_t z)
This method determines if a given position in z is inside the drift volume drifting distance for this...
TRestDetectorReadoutModule * GetModuleByID(Int_t modID)
Returns a pointer to a module using its internal module id.
void SetNormal(const TVector3 &normal)
It updates the value of the normal vector and recalculates the corresponding X and Y axis.
TH2Poly * GetReadoutHistogram()
Creates and returns a TH2Poly object with the readout pixel description.
Double_t fRotation
Rotation (in radians) of the readout plane around the normal vector.
Int_t GetModuleIDFromPosition(const TVector3 &position) const
This method returns the module id where pos is found. The z-coordinate must be found in between the c...
Double_t fChargeCollection
The fraction of charge/energy this readout plane collects from a hit position.
void Print(Int_t DetailLevel=0)
Prints information with details of the readout plane and modules defined inside the readout plane.
TRestDetectorReadoutPlane()
Default TRestDetectorReadoutPlane constructor.
bool IsInside(const TVector3 &point) const
Check if the point is inside any module of the readout plane.
size_t GetNumberOfModules() const
Returns the total number of modules in the readout plane.
void Draw()
Draws the readout plane using GetReadoutHistogram.
Double_t GetX(Int_t modID, Int_t chID)
Returns the X coordinate of a given channel in a given module using their internal module and channel...
Double_t GetDistanceTo(const TVector3 &pos) const
Returns the perpendicular distance to the readout plane from a given position pos.
Int_t FindChannel(Int_t module, const TVector2 &position)
Finds the readout channel index for a given module stored in a given module index stored in the reado...
Double_t fHeight
A length in mm that confers a 3rd dimension to the readout plane and defines a volume.
void GetBoundaries(double &xmin, double &xmax, double &ymin, double &ymax)
Finds the xy boundaries of the readout plane delimited by the modules.
Double_t GetY(Int_t modID, Int_t chID)
Returns the Y coordinate of a given channel in a given module using their internal module and channel...
virtual ~TRestDetectorReadoutPlane()
Default TRestDetectorReadoutPlane destructor.
TVector3 fPosition
The position of the readout plane that defines the origin (0,0) in plane coordinates.
TVector3 GetCathodePosition() const
Returns a TVector3 with the cathode position.
void SetHeight(Double_t d)
Used to define the height of the readout volume with sign crosscheck.
TVector3 fAxisX
A vector contained in the plane that defines the plane X-axis.
TVector3 fAxisY
A vector contained in the plane that defines the plane Y-axis.
TVector3 fNormal
A vector that defines the plane orientation and the side of the active volume.
Int_t GetID() const
Returns an integer with the plane id number.
std::vector< TRestDetectorReadoutModule > fReadoutModules
< A list of TRestDetectorReadoutModule components contained in the readout plane.
void AddModule(const TRestDetectorReadoutModule &module)
Adds a new module to the readout plane.