REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestPatternMask.cxx
1 /*************************************************************************
2  * This file is part of the REST software framework. *
3  * *
4  * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5  * For more information see https://gifna.unizar.es/trex *
6  * *
7  * REST is free software: you can redistribute it and/or modify *
8  * it under the terms of the GNU General Public License as published by *
9  * the Free Software Foundation, either version 3 of the License, or *
10  * (at your option) any later version. *
11  * *
12  * REST is distributed in the hope that it will be useful, *
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15  * GNU General Public License for more details. *
16  * *
17  * You should have a copy of the GNU General Public License along with *
18  * REST in $REST_PATH/LICENSE. *
19  * If not, see https://www.gnu.org/licenses/. *
20  * For the list of contributors see $REST_PATH/CREDITS. *
21  *************************************************************************/
22 
70 
71 #include "TRestPatternMask.h"
72 
73 #include "TAxis.h"
74 #include "TColor.h"
75 #include "TGraph.h"
76 #include "TH1F.h"
77 #include "TRandom3.h"
78 
79 ClassImp(TRestPatternMask);
80 
85 
100 TRestPatternMask::TRestPatternMask(const char* cfgFileName, std::string name) : TRestMetadata(cfgFileName) {
102 }
103 
108 
116 Int_t TRestPatternMask::ApplyCommonMaskTransformation(Double_t& x, Double_t& y) {
117  if (fMaskRadius > 0 && x * x + y * y > fMaskRadius * fMaskRadius) return 0;
118 
119  TVector2 pos(x, y);
120 
121  pos = pos.Rotate(-fRotationAngle);
122  pos -= fOffset;
123 
124  x = pos.X();
125  y = pos.Y();
126 
127  return 1;
128 }
129 
134 bool TRestPatternMask::HitsPattern(Double_t x, Double_t y) {
135  if (GetRegion(x, y) > 0) return false;
136  return true;
137 }
138 
144 
146  RESTMetadata << "----" << RESTendl;
147 }
148 
153  RESTMetadata << " - Pattern type : " << fPatternType << RESTendl;
154  RESTMetadata << " - Mask radius : " << fMaskRadius << " mm" << RESTendl;
155  RESTMetadata << " - Max regions : " << fMaxRegions << RESTendl;
156  RESTMetadata << " - Offset : (" << fOffset.X() << ", " << fOffset.Y() << ") mm" << RESTendl;
157  RESTMetadata << " - Rotation angle : " << fRotationAngle * 180. / TMath::Pi() << " degrees" << RESTendl;
158 }
159 
165 TCanvas* TRestPatternMask::DrawMonteCarlo(Int_t nSamples) {
166  if (fCanvas != NULL) {
167  delete fCanvas;
168  fCanvas = NULL;
169  }
170  fCanvas = new TCanvas("canv", "This is the canvas title", 1200, 1200);
171  fCanvas->Draw();
172 
173  TPad* pad1 = new TPad("pad1", "This is pad1", 0.01, 0.02, 0.99, 0.97);
174  pad1->Draw();
175 
176  fCanvas->SetRightMargin(0.09);
177  fCanvas->SetLeftMargin(0.15);
178  fCanvas->SetBottomMargin(0.15);
179 
181  std::map<int, std::vector<TVector2> > points;
182 
183  TRandom3* rnd = new TRandom3(0);
184 
185  for (int n = 0; n < nSamples; n++) {
186  Double_t xO = 2.5 * (rnd->Rndm() - 0.5) * fMaskRadius;
187  Double_t yO = 2.5 * (rnd->Rndm() - 0.5) * fMaskRadius;
188  Double_t x = xO;
189  Double_t y = yO;
190 
191  Int_t id = GetRegion(x, y);
192 
193  if (points.count(id) == 0) {
194  std::vector<TVector2> a;
195  a.push_back(TVector2(xO, yO));
196  points[id] = a;
197  } else {
198  points[id].push_back(TVector2(xO, yO));
199  }
200  }
201 
202  std::vector<TGraph*> gridGraphs;
203  Int_t nGraphs = 0;
204  for (const auto& x : points) {
205  std::string grname = "gr" + IntegerToString(nGraphs);
206  TGraph* gr = new TGraph();
207  gr->SetName(grname.c_str());
208 
209  for (unsigned int n = 0; n < x.second.size(); n++) {
210  gr->SetPoint(gr->GetN(), x.second[n].X(), x.second[n].Y());
211  }
212 
213  if (nGraphs == 0) {
214  gr->SetLineColor(kBlack);
215  gr->SetMarkerColor(kBlack);
216  gr->SetMarkerSize(0.6);
217  gr->SetLineWidth(2);
218  } else {
219  gr->SetMarkerColor((nGraphs * 3) % 29 + 20);
220  gr->SetLineColor((nGraphs * 3) % 29 + 20);
221  gr->SetMarkerSize(0.6);
222  gr->SetLineWidth(2);
223  }
224 
225  gr->SetMarkerStyle(20);
226  gridGraphs.push_back(gr);
227  nGraphs++;
228  }
229 
230  gridGraphs[0]->GetXaxis()->SetLimits(-1.25 * fMaskRadius, 1.25 * fMaskRadius);
231  gridGraphs[0]->GetHistogram()->SetMaximum(1.25 * fMaskRadius);
232  gridGraphs[0]->GetHistogram()->SetMinimum(-1.25 * fMaskRadius);
233 
234  gridGraphs[0]->GetXaxis()->SetTitle("X [mm]");
235  gridGraphs[0]->GetXaxis()->SetTitleSize(0.04);
236  gridGraphs[0]->GetXaxis()->SetLabelSize(0.04);
237  gridGraphs[0]->GetYaxis()->SetTitle("Y[mm]");
238  gridGraphs[0]->GetYaxis()->SetTitleOffset(1.5);
239  gridGraphs[0]->GetYaxis()->SetTitleSize(0.04);
240  gridGraphs[0]->GetYaxis()->SetLabelSize(0.04);
241  pad1->SetLogy();
242  gridGraphs[0]->Draw("AP");
243  for (int n = 1; n < nGraphs; n++) gridGraphs[n]->Draw("P");
244 
245  return fCanvas;
246 }
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
virtual void Initialize()
Making default settings.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
An abstract class used to encapsulate different mask pattern class definitions.
TCanvas * DrawMonteCarlo(Int_t nSamples=10000)
It generates a Monte Carlo with positions and paints them the returned canvas. Each unique region is ...
TVector2 fOffset
It is used to introduce an offset on the pattern (not the mask, mask is always centered)
TRestPatternMask()
Default constructor.
TCanvas * fCanvas
A canvas for drawing.
Int_t fMaxRegions
The maximum number of regions allowed in each mask.
Bool_t HitsPattern(Double_t x, Double_t y)
Returns true if the pattern was hit. If (x,y) it is inside a region then, the pattern was not hit by ...
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestPatternMask.
Double_t fRotationAngle
An angle (in radians) used to introduce a rotation to the pattern.
Int_t ApplyCommonMaskTransformation(Double_t &x, Double_t &y)
It produces an effective mask rotation and translation for the point x,y.
virtual Int_t GetRegion(Double_t &x, Double_t &y)
To be implemented at the inherited class with the pattern and region identification logic.
~TRestPatternMask()
Default destructor.
Double_t fMaskRadius
The maximum mask radius in mm (if 0 it will be infinite)
std::string fPatternType
The pattern type (None/Stripped/Grid/Spider/Rings)
void PrintCommonPatternMembers()
Prints on screen the information about the metadata members without header.
@ REST_Info
+show most of the information for each steps
std::string IntegerToString(Int_t n, std::string format="%d")
Gets a string from an integer.