12#include "TRestDetectorHits3DReconstructionProcess.h"
28 bool iszombie()
const {
return e == 0; }
29 void add(
const regionhit& h) {
30 x = (x * e + h.x * h.e) / (e + h.e);
31 y = (y * e + h.y * h.e) / (e + h.e);
34 static double dist(
const regionhit& h1,
const regionhit& h2) {
35 return sqrt((h1.x - h2.x) * (h1.x - h2.x) + (h1.y - h2.y) * (h1.y - h2.y));
47 void Print() { cout <<
"x1: " << x1 <<
", x2: " << x2 <<
", y1: " << y1 <<
", y2: " << y2 << endl; }
49 bool iszombie()
const {
return e == 0; }
51 static bool IsRectCross(
const stripehit& s1,
const stripehit& s2) {
52 bool ret = min(s1.x1, s1.x2) <= max(s2.x1, s2.x2) && min(s2.x1, s2.x2) <= max(s1.x1, s1.x2) &&
53 min(s1.y1, s1.y2) <= max(s2.y1, s2.y2) && min(s2.y1, s2.y2) <= max(s1.y1, s1.y2);
57 static bool IsLineSegmentCross(
const stripehit& s1,
const stripehit& s2) {
58 if (((s2.x1 - s1.x1) * (s2.y1 - s2.y2) - (s2.y1 - s1.y1) * (s2.x1 - s2.x2)) *
59 ((s2.x1 - s1.x2) * (s2.y1 - s2.y2) - (s2.y1 - s1.y2) * (s2.x1 - s2.x2)) <
61 ((s1.x1 - s2.x1) * (s1.y1 - s1.y2) - (s1.y1 - s2.y1) * (s1.x1 - s1.x2)) *
62 ((s1.x1 - s2.x2) * (s1.y1 - s1.y2) - (s1.y1 - s2.y2) * (s1.x1 - s1.x2)) <
69 static bool GetCrossPoint(
const stripehit& s1,
const stripehit& s2,
double& x,
double& y) {
70 if (IsRectCross(s1, s2)) {
71 if (IsLineSegmentCross(s1, s2)) {
72 long tmpLeft, tmpRight;
73 tmpLeft = (s2.x2 - s2.x1) * (s1.y1 - s1.y2) - (s1.x2 - s1.x1) * (s2.y1 - s2.y2);
74 tmpRight = (s1.y1 - s2.y1) * (s1.x2 - s1.x1) * (s2.x2 - s2.x1) +
75 s2.x1 * (s2.y2 - s2.y1) * (s1.x2 - s1.x1) -
76 s1.x1 * (s1.y2 - s1.y1) * (s2.x2 - s2.x1);
78 x = (int)((
double)tmpRight / (double)tmpLeft);
80 tmpLeft = (s1.x1 - s1.x2) * (s2.y2 - s2.y1) - (s1.y2 - s1.y1) * (s2.x1 - s2.x2);
81 tmpRight = s1.y2 * (s1.x1 - s1.x2) * (s2.y2 - s2.y1) +
82 (s2.x2 - s1.x2) * (s2.y2 - s2.y1) * (s1.y1 - s1.y2) -
83 s2.y2 * (s2.x1 - s2.x2) * (s1.y2 - s1.y1);
84 y = (int)((
double)tmpRight / (double)tmpLeft);
91 static double dist(
const stripehit& h1,
const stripehit& h2) {
92 if (IsRectCross(h1, h2) && IsLineSegmentCross(h1, h2)) {
97 double d1 = abs(((x - h2.x1) * (h2.y2 - h2.y1) - (h2.x2 - h2.x1) * (y - h2.y1)) /
98 sqrt((h2.x2 - h2.x1) * (h2.x2 - h2.x1) + (h2.y2 - h2.y1) * (h2.y2 - h2.y1)));
102 double d2 = abs(((x - h2.x1) * (h2.y2 - h2.y1) - (h2.x2 - h2.x1) * (y - h2.y1)) /
103 sqrt((h2.x2 - h2.x1) * (h2.x2 - h2.x1) + (h2.y2 - h2.y1) * (h2.y2 - h2.y1)));
107 double d3 = abs(((x - h1.x1) * (h1.y2 - h1.y1) - (h1.x2 - h1.x1) * (y - h1.y1)) /
108 sqrt((h1.x2 - h1.x1) * (h1.x2 - h1.x1) + (h1.y2 - h1.y1) * (h1.y2 - h1.y1)));
112 double d4 = abs(((x - h1.x1) * (h1.y2 - h1.y1) - (h1.x2 - h1.x1) * (y - h1.y1)) /
113 sqrt((h1.x2 - h1.x1) * (h1.x2 - h1.x1) + (h1.y2 - h1.y1) * (h1.y2 - h1.y1)));
115 return min({d1, d2, d3, d4});
119 friend regionhit operator&(
const stripehit& s1,
const stripehit& s2) {
120 if (s1.type * s2.type == XZ * YZ) {
121 if (!s1.iszombie() && !s2.iszombie()) {
124 if (stripehit::GetCrossPoint(s1, s2, x, y)) {
137 type = hits->GetType(i);
139 y1 = hits->GetY(i) - YWidth;
140 y2 = hits->GetY(i) + YWidth;
143 }
else if (type == YZ) {
146 x1 = hits->GetX(i) - XWidth;
147 x2 = hits->GetX(i) + XWidth;
151 e = hits->GetEnergy(i);
154 stripehit(
double x1,
double x2,
double y1,
double y2,
double e = 0, REST_HitType type = unknown) {
166 TRestDetectorHits3DReconstructionProcess::TRestDetectorHits3DReconstructionProcess() {
170TRestDetectorHits3DReconstructionProcess::~TRestDetectorHits3DReconstructionProcess() {
171 delete fOutputHitsEvent;
178 fInputHitsEvent =
nullptr;
193 htemp =
new TH2D(
"hhits",
"hhits", 100, -PlaneMaxX, PlaneMaxX, 100, -PlaneMaxY, PlaneMaxY);
197 fCompareProc = GetFriendLive(
"TRestGeant4ToDetectorHitsProcess");
205 fInputHitsEvent->Sort();
208 double totalambiguity = 0;
210 for (
unsigned int n = 0; n < fInputHitsEvent->GetNumberOfHits(); n++) {
211 double z = fInputHitsEvent->
GetZ(n);
212 if (z == lastz)
continue;
217 map<int, double> ids;
219 for (
int i = n - 1; i >= 0; i--) {
220 if (abs(fInputHitsEvent->
GetZ(i) - z) < fZRange * 2) {
221 ids[i] = TMath::Gaus(abs(fInputHitsEvent->
GetZ(i) - z), 0, fZRange);
226 for (
unsigned int i = n + 1; i < fInputHitsEvent->GetNumberOfHits(); i++) {
227 if (abs(fInputHitsEvent->
GetZ(i) - z) < fZRange * 2) {
228 ids[i] = TMath::Gaus(abs(fInputHitsEvent->
GetZ(i) - z), 0, fZRange);
235 vector<stripehit> stripehits;
237 stripehit h = stripehit(fInputHitsEvent->GetHits(), i.first);
240 if (stripehits.size() > 0) {
241 for (
unsigned int k = 0; k < stripehits.size(); k++) {
242 stripehit& _h = stripehits[k];
243 if (_h.x1 == h.x1 && _h.x2 == h.x2 && _h.y1 == h.y1 && _h.y2 == h.y2) {
246 }
else if (k == stripehits.size() - 1) {
247 stripehits.push_back(h);
251 stripehits.push_back(h);
256 if (stripehits.size() > 50) {
257 RESTWarning <<
"(TRestDetectorHits3DReconstructionProcess) event id: " << fInputHitsEvent->GetID()
258 <<
", z: " << z <<
", bad frame, too much hits! (" << stripehits.size() <<
")"
264 cout <<
"stripehits: ";
265 for (
auto h : stripehits) {
266 cout << h.x1 <<
"," << h.y1 <<
"," << h.x2 <<
"," << h.y2 <<
"," << h.e <<
" ";
272 vector<regionhit> regionhits;
273 double layerambiguity = 0;
274 set<int> ambcalculatedstripes;
275 for (
unsigned int i = 0; i < stripehits.size(); i++) {
276 vector<int> crossedids;
277 double crossedsumenergy = 0;
278 for (
unsigned int j = 0; j < stripehits.size(); j++) {
279 if (i == j)
continue;
280 if (stripehit::IsRectCross(stripehits[i], stripehits[j])) {
281 crossedids.push_back(j);
282 crossedsumenergy += stripehits[j].e;
286 for (
int id : crossedids) {
287 regionhit rh = stripehits[i] & stripehits[id];
288 rh.e = stripehits[i].e * stripehits[id].e / crossedsumenergy;
290 if (regionhits.size() > 0) {
291 for (
unsigned int k = 0; k < regionhits.size(); k++) {
292 regionhit& _rh = regionhits[k];
293 if (_rh.x == rh.x && _rh.y == rh.y) {
296 }
else if (k == regionhits.size() - 1) {
297 regionhits.push_back(rh);
301 regionhits.push_back(rh);
306 if (ambcalculatedstripes.count(i) == 0 && crossedids.size() > 0) {
307 for (
int id : crossedids) {
308 ambcalculatedstripes.insert(
id);
311 auto line = stripehits[i];
312 auto newline = stripehits[crossedids[0]];
313 vector<int> crossedids_newline;
314 for (
unsigned int j = 0; j < stripehits.size(); j++) {
315 if (j == (
unsigned int)crossedids[0])
continue;
316 if (stripehit::IsRectCross(stripehits[j], newline)) {
317 crossedids_newline.push_back(j);
320 for (
int id : crossedids_newline) {
321 ambcalculatedstripes.insert(
id);
324 vector<pair<double, double>> V_PosE;
325 vector<pair<double, double>> H_PosE;
327 if (line.x1 == line.x2) {
330 for (
int id : crossedids_newline) {
331 H_PosE.push_back({stripehits[id].x1, stripehits[id].e});
333 }
else if (line.y1 == line.y2) {
336 for (
int id : crossedids_newline) {
337 V_PosE.push_back({stripehits[id].y1, stripehits[id].e});
341 if (newline.x1 == newline.x2) {
343 for (
int id : crossedids) {
344 H_PosE.push_back({stripehits[id].x1, stripehits[id].e});
346 }
else if (newline.y1 == newline.y2) {
348 for (
int id : crossedids) {
349 V_PosE.push_back({stripehits[id].y1, stripehits[id].e});
354 auto comp = [](
const pair<double, double>& p1,
const pair<double, double>& p2) ->
bool {
355 return p1.first < p2.first;
357 sort(V_PosE.begin(), V_PosE.end(), comp);
358 sort(H_PosE.begin(), H_PosE.end(), comp);
360 for (
unsigned int j = 1; j < V_PosE.size(); j++) {
361 if (V_PosE[j].first - V_PosE[j - 1].first > PITCH) {
362 V_PosE.insert(V_PosE.begin() + j, {(V_PosE[j].first + V_PosE[j - 1].first) / 2, 0});
366 V_PosE.insert(V_PosE.begin(), {(*V_PosE.begin()).first - 1, 0});
367 V_PosE.insert(V_PosE.end(), {(*(V_PosE.end() - 1)).first + 1, 0});
369 for (
unsigned int j = 1; j < H_PosE.size(); j++) {
370 if (H_PosE[j].first - H_PosE[j - 1].first > PITCH) {
371 H_PosE.insert(H_PosE.begin() + j, {(H_PosE[j].first + H_PosE[j - 1].first) / 2, 0});
375 H_PosE.insert(H_PosE.begin(), {(*H_PosE.begin()).first - 1, 0});
376 H_PosE.insert(H_PosE.end(), {(*(H_PosE.end() - 1)).first + 1, 0});
380 for (
auto pose : V_PosE) {
381 cout << pose.first <<
"," << pose.second <<
" ";
386 for (
auto pose : H_PosE) {
387 cout << pose.first <<
"," << pose.second <<
" ";
395 for (
unsigned int j = 1; j < V_PosE.size() - 1; j++) {
396 if (V_PosE[j].second > V_PosE[j - 1].second && V_PosE[j].second > V_PosE[j + 1].second) {
401 for (
unsigned int j = 1; j < H_PosE.size() - 1; j++) {
402 if (H_PosE[j].second > H_PosE[j - 1].second && H_PosE[j].second > H_PosE[j + 1].second) {
408 cout << i <<
" z: " << z <<
" " << Vsegs <<
" " << Hsegs <<
" -> "
409 << LogAmbiguity(Vsegs, Hsegs) << endl;
413 if (Vsegs == 0) Vsegs = 1;
414 if (Hsegs == 0) Hsegs = 1;
416 layerambiguity += LogAmbiguity(Vsegs, Hsegs);
420 totalambiguity += layerambiguity;
423 vector<regionhit> regionhits_merged;
424 set<int> mergeflag_all;
425 while (mergeflag_all.size() < regionhits.size()) {
429 for (
unsigned int i = 0; i < regionhits.size(); i++) {
430 if (mergeflag_all.count(i) == 0) {
431 if (regionhits[i].e > max) {
432 max = regionhits[i].e;
438 mergeflag.insert(maxposition);
441 for (
auto _i = mergeflag.begin(); _i != mergeflag.end();) {
448 bool inserted =
false;
449 for (
unsigned int j = 0; j < regionhits.size(); j++) {
450 if (regionhit::dist(regionhits[i], regionhits[j]) < PITCH &&
451 regionhits[i].e > regionhits[j].e && mergeflag_all.count(j) == 0) {
452 auto insertresult = mergeflag.insert(j);
453 inserted += insertresult.second;
458 _i = mergeflag.begin();
468 for (
auto i : mergeflag) {
469 mergeflag_all.insert(i);
470 h.add(regionhits[i]);
472 regionhits_merged.push_back(h);
474 if (regionhits_merged.size() > 0) Nlayers++;
476 if (fDraw && (layerambiguity > fDrawThres)) {
477 if (regionhits.size() == 0)
continue;
478 double maxenergyofregionhits =
479 (*min_element(regionhits.begin(), regionhits.end(),
480 [](
const regionhit& h1,
const regionhit& h2) ->
bool { return h1.e > h2.e; }))
484 htemp->SetStats(
false);
485 htemp->SetTitle(Form(
"event id %i, z = %.1f", fInputHitsEvent->GetID(), z));
487 vector<TLine*> lines;
488 vector<TMarker*> markers;
489 for (
auto h : stripehits) {
490 TLine* line =
new TLine(h.x1, h.y1, h.x2, h.y2);
492 lines.push_back(line);
494 for (
auto rh : regionhits) {
495 TMarker* marker =
new TMarker(rh.x, rh.y, kFullDotLarge);
496 marker->SetMarkerColor(kGreen);
497 marker->SetMarkerSize(rh.e / maxenergyofregionhits);
499 markers.push_back(marker);
501 for (
auto rhm : regionhits_merged) {
502 TMarker* marker =
new TMarker(rhm.x, rhm.y, kFullDotLarge);
503 marker->SetMarkerColor(kRed);
505 markers.push_back(marker);
510 for (
auto l : lines) {
513 for (
auto m : markers) {
519 for (
auto rh : regionhits_merged) {
520 fOutputHitsEvent->
AddHit(rh.x, rh.y, z, rh.e, 0, XYZ);
525 if (fDoEnergyScaling) {
526 double e1 = fInputHitsEvent->GetTotalEnergy();
527 double e2 = fOutputHitsEvent->GetTotalEnergy();
528 for (
auto h : *fOutputHitsEvent->GetHits()) {
529 h.e() = h.e() / e2 * e1;
535 if (fCompareProc !=
nullptr && fOutputHitsEvent->GetNumberOfHits() > 0 &&
538 auto hits1 = *fOutputHitsEvent->GetHits();
539 auto hits2 = *reference->GetHits();
541 double maxz1 = (*min_element(hits1.begin(), hits1.end(),
545 double maxz2 = (*min_element(hits2.begin(), hits2.end(),
549 double dz = maxz2 - maxz1;
551 double mindistmean = std::numeric_limits<Double_t>::max();
552 for (
double i = -fZRange * 2; i <= fZRange * 2; i += fZRange / 2) {
554 for (
auto hit : hits1) {
555 auto pos = TVector3(hit.x(), hit.y(), hit.z() + dz);
556 int n = hits2.GetClosestHit(pos);
557 auto dist2 = (pos - hits2.GetPosition(n)).Mag2();
560 double distmean = sqrt(distsum2) / fOutputHitsEvent->GetNumberOfHits();
561 if (distmean < mindistmean) {
562 mindistmean = distmean;
568 return fOutputHitsEvent;
571int TRestDetectorHits3DReconstructionProcess::Ambiguity(
const int& n,
const int& m) {
572 if (n <= 0 || m <= 0)
return 0;
573 if (n == 1 || m == 1)
return 1;
574 if (n > 5 && m > 5)
return 1e9;
575 if (n > 8 || m > 8)
return 1e9;
580 int ambiguity = pow(N, M);
581 for (
int i = N - 1; i > 0; i--) {
583 int combinations = Factorial(N) / Factorial(N - i) / Factorial(i);
584 ambiguity -= Ambiguity(i, M) * combinations;
589int TRestDetectorHits3DReconstructionProcess::Factorial(
const int& n) {
592 while (N > 1) result *= (--N);
608 fDoEnergyScaling = StringToBool(
GetParameter(
"scaleE",
"true"));
void Initialize() override
Making default settings.
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
void InitProcess() override
To be executed at the beginning of the run (outside event loop)
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void EndProcess() override
To be executed at the end of the run (outside event loop)
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.
TVector2 fCanvasSize
Canvas size.
void BeginPrintProcess()
[name, cut range]
TCanvas * fCanvas
< Canvas for some viewer event
void CreateCanvas()
Create the canvas.
bool fDynamicObs
It defines whether to use added observables only or all the observables appear in the code.
virtual RESTValue GetOutputEvent() const =0
Get pointer to output event. Must be implemented in the derived class.
std::map< std::string, int > fObservablesDefined
Stores the list of all the appeared process observables in the code.
void SetObservableValue(const std::string &name, const T &value)
Set observable value for AnalysisTree.
A base class for any REST event.
It saves a 3-coordinate position and an energy for each punctual deposition.
@ REST_Extreme
show everything
@ REST_Debug
+show the defined debug messages
Int_t GetChar(std::string hint="Press a KEY to continue ...")
Helps to pause the program, printing a message before pausing.
Double_t StringToDouble(std::string in)
Gets a double from a string.
TVector2 StringTo2DVector(std::string in)
Gets a 2D-vector from a string.