REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRadialStrippedMask.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
106
107#include "TRestRadialStrippedMask.h"
108
109#include "TRandom3.h"
110
112
117
132TRestRadialStrippedMask::TRestRadialStrippedMask(const char* cfgFileName, std::string name)
133 : TRestPatternMask(cfgFileName) {
134 Initialize();
135
137
139}
140
145
151 SetSectionName(this->ClassName());
152 SetType("RadialStripped");
153}
154
163Int_t TRestRadialStrippedMask::GetRegion(Double_t& x, Double_t& y) {
164 if (TRestPatternMask::GetRegion(x, y)) return 0;
165
166 Double_t d = TMath::Sqrt(x * x + y * y);
167
168 if (d < fInitialRadius) {
169 if (fInternalRegionRadius > 0 && d < fInternalRegionRadius) return 1;
170
171 return 0;
172 }
173
174 TVector2 point(x, y);
175 Double_t phi = point.Phi();
176
178 Int_t region = (Int_t)(phi / fStripsAngle);
179 region = 2 + region % fMaxRegions;
180
181 Double_t angle = 0;
183 while (angle < 2 * TMath::Pi()) {
184 if (point.Y() < fStripsThickness / 2. && point.Y() > -fStripsThickness / 2. && point.X() >= 0)
185 return 0;
186
187 point = point.Rotate(fStripsAngle);
188 angle += fStripsAngle;
189 }
190
191 return 1 + region % fModulus;
192}
193
199
201 RESTMetadata << "++++" << RESTendl;
202}
203
210 RESTMetadata << "----" << RESTendl;
212}
213
219 RESTMetadata << " - Strips angle : " << fStripsAngle * units("degrees") << " degrees" << RESTendl;
220 RESTMetadata << " - Strips thickness : " << fStripsThickness << " mm" << RESTendl;
221}
endl_t RESTendl
Termination flag object for TRestStringOutput.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string fConfigFileName
Full name of the rml file.
An abstract class used to encapsulate different mask pattern class definitions.
void SetType(const std::string &type)
It defines the mask type. To be called by the inherited class constructor.
Int_t fMaxRegions
The maximum number of regions allowed in each mask.
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestPatternMask.
virtual Int_t GetRegion(Double_t &x, Double_t &y)
To be implemented at the inherited class with the pattern and region identification logic.
void PrintCommonPatternMembers()
Prints on screen the information about the metadata members without header.
A class used to define a stripped mask pattern.
Double_t fStripsAngle
The periodity of the stripped structure in radians.
Double_t fInitialRadius
The spacers structure will be effective from this radius, in mm. Default is from 20 mm.
virtual Int_t GetRegion(Double_t &x, Double_t &y) override
It returns a number identifying the region where the particle with coordinates (x,...
void Initialize() override
Function to initialize input/output event members and define the section name.
Double_t fStripsThickness
The width of the stripped structure in mm.
Int_t fModulus
It defines the maximum number of cells/regions in each axis.
TRestRadialStrippedMask()
Default constructor.
void PrintMaskMembers() override
Prints on screen the information about the metadata members of TRestRingsMask, excluding common metad...
~TRestRadialStrippedMask()
Default destructor.
void PrintMask() override
Prints on screen the information about the metadata members of TRestRingsMask, including common patte...
void PrintMetadata() override
Prints on screen the complete information about the metadata members from this class.
Double_t fInternalRegionRadius
Radius of an internal circular region defined inside the fInitialRadius. If 0, there will be no regio...
@ REST_Info
+show most of the information for each steps