REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestPatternMask.h
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
23#ifndef REST_TRestPatternMask
24#define REST_TRestPatternMask
25
26#include <TCanvas.h>
27#include <TRestMetadata.h>
28
31 private:
33 TVector2 fOffset = TVector2(0, 0); //<
34
36 Double_t fRotationAngle = 0; //<
37
39 std::string fPatternType = "None"; //<
40
42 TCanvas* fCanvas = nullptr;
43
44 protected:
46 Double_t fMaskRadius = 0; //<
47
49 Int_t fMaxRegions = 100; //<
50
52 void SetType(const std::string& type) { fPatternType = type; }
53
54 Int_t ApplyCommonMaskTransformation(Double_t& x, Double_t& y);
55
56 public:
57 Int_t GetMaxRegions() { return fMaxRegions; }
58
59 void SetMaxRegions(Int_t regions) { fMaxRegions = regions; }
60
61 Bool_t HitsPattern(Double_t x, Double_t y);
62
64 virtual Int_t GetRegion(Double_t& x, Double_t& y) {
65 if (ApplyCommonMaskTransformation(x, y) == 0) return 1;
66 return 0;
67 }
68
70 std::string GetType() { return fPatternType; }
71
73 Double_t GetRotationAngle() { return fRotationAngle; }
74
76 TVector2 GetOffset() { return fOffset; }
77
79 Double_t GetMaskRadius() { return fMaskRadius; }
80
82 void SetRotationAngle(const Double_t& angle) { fRotationAngle = angle; }
83
85 void SetOffset(const TVector2& offset) { fOffset = offset; }
86
88 void SetMaskRadius(const Double_t& radius) { fMaskRadius = radius; }
89
91
92 virtual void PrintMaskMembers() {}
93 virtual void PrintMask() { PrintCommonPatternMembers(); }
94
95 void PrintMetadata() override;
96
97 TCanvas* DrawMonteCarlo(Int_t nSamples = 10000);
98
100 TRestPatternMask(const char* cfgFileName, std::string name = "");
102
103 ClassDefOverride(TRestPatternMask, 1);
104};
105#endif
A base class for any REST metadata class.
Definition: TRestMetadata.h:70
An abstract class used to encapsulate different mask pattern class definitions.
Double_t GetMaskRadius()
It returns the mask radius.
TCanvas * DrawMonteCarlo(Int_t nSamples=10000)
It generates a Monte Carlo with positions and paints them the returned canvas. Each unique region is ...
void SetType(const std::string &type)
It defines the mask type. To be called by the inherited class constructor.
TVector2 fOffset
It is used to introduce an offset on the pattern.
TRestPatternMask()
Default constructor.
TVector2 GetOffset()
It returns the rotation angle.
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 SetRotationAngle(const Double_t &angle)
It defines the rotation angle.
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.
Double_t GetRotationAngle()
It returns the rotation angle.
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)
void SetOffset(const TVector2 &offset)
It defines the pattern offset.
std::string fPatternType
The pattern type (None/Stripped/Grid/Spider/Rings)
void SetMaskRadius(const Double_t &radius)
It defines the mask radius.
std::string GetType()
It returns the mask pattern type.
void PrintCommonPatternMembers()
Prints on screen the information about the metadata members without header.