REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionOptics.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 http://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 http://www.gnu.org/licenses/. *
20  * For the list of contributors see $REST_PATH/CREDITS. *
21  *************************************************************************/
22 
23 #ifndef _TRestAxionOptics
24 #define _TRestAxionOptics
25 
26 #include <TRestAxionOpticsMirror.h>
27 #include <TRestCombinedMask.h>
28 #include <TRestMetadata.h>
29 
30 #include <iostream>
31 
32 #include "TRandom3.h"
33 
36  protected:
38  std::string fOpticsFile = "";
39 
41  Double_t fMirrorLength = 0; //<
42 
44  std::vector<std::vector<Double_t>> fOpticsData; //<
45 
48 
51 
54 
57 
59  TVector3 fOriginPosition;
60 
62  TVector3 fEntrancePosition;
63 
66 
68  TVector3 fMiddlePosition;
69 
72 
74  TVector3 fExitPosition;
75 
77  TVector3 fEntranceDirection;
78 
80  TVector3 fMiddleDirection;
81 
83  TVector3 fExitDirection;
84 
86  Double_t fFocal = -1;
87 
89  Int_t fCurrentMirror = -1;
90 
92  Bool_t fFirstInteraction = false;
93 
95  Bool_t fSecondInteraction = false;
96 
98  TRandom3* fRandom = nullptr;
99 
101  TRestAxionOptics(const char* cfgFileName, std::string name = "");
102 
103  private:
104  void ResetPositions();
105 
106  Double_t GetPhotonReflectivity(Double_t energy);
107 
108  protected:
110  TPad* fPad = nullptr;
111 
112  Int_t TransportToEntrance(const TVector3& pos, const TVector3& dir);
113  Int_t TransportToMiddle(const TVector3& pos, const TVector3& dir);
114  Int_t TransportToExit(const TVector3& pos, const TVector3& dir);
115 
117  Int_t GetMirror() { return fCurrentMirror; }
118 
120  virtual void SetMirror() = 0;
121 
122  public:
123  virtual void Initialize();
124 
126  virtual std::pair<Double_t, Double_t> GetRadialLimits() = 0;
127 
129  virtual Double_t GetEntrancePositionZ() = 0;
130 
132  virtual Double_t GetExitPositionZ() = 0;
133 
135  virtual Int_t FirstMirrorReflection(const TVector3& pos, const TVector3& dir) = 0;
136 
138  virtual Int_t SecondMirrorReflection(const TVector3& pos, const TVector3& dir) = 0;
139 
141  virtual TPad* DrawMirrors() = 0;
142 
144  Double_t GetEntranceAngle() { return TMath::ACos(fEntranceDirection.Dot(TVector3(0, 0, 1))); }
145 
146  virtual Double_t FindFocal(Double_t from, Double_t to, Double_t energy, Double_t precision = 1,
147  Bool_t recalculate = false, Int_t particles = 5000);
148 
149  Double_t CalculateSpotSize(Double_t energy, Double_t z, Int_t particles = 15000);
150 
151  TPad* CreatePad(Int_t nx = 1, Int_t ny = 1);
152 
153  TPad* DrawParticleTracks(Double_t deviation = 0, Int_t particles = 10);
154 
155  TPad* DrawScatterMaps(Double_t z, Double_t energy = 0, Double_t deviation = 0, Int_t particles = 1000,
156  Double_t focalHint = 7500);
157 
158  TPad* DrawDensityMaps(Double_t z, Double_t energy = 0, Double_t deviation = 0, Int_t particles = 1000,
159  Double_t focalHint = 7500);
160 
161  Double_t PropagatePhoton(const TVector3& pos, const TVector3& dir, Double_t energy);
162 
163  Double_t PropagateMonteCarloPhoton(Double_t energy, Double_t deviation);
164 
166  TVector3 GetEntrancePosition() { return fEntrancePosition; }
167 
169  TVector3 GetMiddlePosition() { return fMiddlePosition; }
170 
172  TVector3 GetExitPosition() { return fExitPosition; }
173 
176 
178  TVector3 GetMiddleDirection() { return fMiddleDirection; }
179 
181  TVector3 GetExitDirection() { return fExitDirection; }
182 
183  TVector3 GetLastGoodPosition();
184  TVector3 GetLastGoodDirection();
185 
188 
191 
192  Int_t GetNumberOfReflections();
193 
195  TRestCombinedMask* const& GetEntranceMask() const { return fEntranceMask; }
196 
198  TRestCombinedMask* const& GetMiddleMask() const { return fMiddleMask; }
199 
201  TRestCombinedMask* const& GetExitMask() const { return fExitMask; }
202 
205 
206  void PrintMetadata();
207 
208  void PrintMirror();
209  void PrintMasks();
210 
211  void PrintEntranceMask();
212  void PrintMiddleMask();
213  void PrintExitMask();
214 
216 
217  void InitFromConfigFile();
218 
220 
221  ClassDef(TRestAxionOptics, 1);
222 };
223 #endif
A metadata class accessing the Henke database to load reflectivity data.
An abstract class to define common optics parameters and methods.
TVector3 GetLastGoodPosition()
It returns the last valid particle position known in the particle tracking.
Bool_t fSecondInteraction
During the photon propagation it tells us if the photon interacted in the second mirror.
void PrintEntranceMask()
Prints on screen the mask used on the entrance optics plane.
virtual Double_t FindFocal(Double_t from, Double_t to, Double_t energy, Double_t precision=1, Bool_t recalculate=false, Int_t particles=5000)
It implements a generic method to identify the optimum focal point. It can be reimplemented at each s...
void PrintPhotonTrackingSummary()
Prints the positions taken by the photon after ray-tracing that should have been updated using the me...
TVector3 fMiddlePosition
The particle position at the optics plane middle.
Int_t fCurrentMirror
During the photon propagation it keeps track of the active mirror shell.
TPad * fPad
A pad pointer to be used by the drawing methods.
TRestCombinedMask * fMiddleMask
The middle optical mask that defines a pattern with regions identified with a number.
Bool_t IsSecondMirrorReflection()
It returns true if the photon got reflected in the second mirror.
std::vector< std::vector< Double_t > > fOpticsData
The optics data table extracted from fOpticsFile.
TVector3 fMiddleDirection
The particle position at the optics plane middle.
Double_t PropagatePhoton(const TVector3 &pos, const TVector3 &dir, Double_t energy)
Propagating photon.
virtual Int_t SecondMirrorReflection(const TVector3 &pos, const TVector3 &dir)=0
It updates the values fSecondInteractionPosition and fExitDirection. Returns 0 if is not in region.
Bool_t fFirstInteraction
During the photon propagation it tells us if the photon interacted in the first mirror.
TVector3 fExitDirection
The particle position at the optics plane exit.
Double_t PropagateMonteCarloPhoton(Double_t energy, Double_t deviation)
It will produce a MonteCarlo photon spatially distributed in XY as defined by the GetRadialLimits met...
TRestCombinedMask * fExitMask
The exit optical mask that defines a pattern with regions identified with a number.
Double_t GetEntranceAngle()
It returns the entrance angle to the optical axis (in radians).
virtual void Initialize()
Initialization of TRestAxionOptics members.
Double_t GetPhotonReflectivity(Double_t energy)
A prototype method to be implemented by specific optics to draw an schematic including the mirrors ge...
TRestAxionOpticsMirror * fMirrorProperties
The mirror properties.
Int_t TransportToEntrance(const TVector3 &pos, const TVector3 &dir)
It moves the incoming particle at the entrance of the optics plane and returns the region number wher...
Int_t TransportToMiddle(const TVector3 &pos, const TVector3 &dir)
It moves the incoming particle to the middle of the optics plane and returns the region number where ...
TRestAxionOptics()
Default constructor.
virtual Double_t GetEntrancePositionZ()=0
It returns the entrance Z-position defined by the optical axis.
Bool_t IsFirstMirrorReflection()
It returns true if the photon got reflected in the first mirror.
TVector3 fExitPosition
The particle position at the optics plane exit.
Double_t fMirrorLength
The mirror length. If all mirrors got the same length. Otherwise will be zero.
TPad * CreatePad(Int_t nx=1, Int_t ny=1)
A prototype method to be implemented by specific optics to draw an schematic including the mirrors ge...
Double_t fFocal
The calculated focal in mm. It is updated by the method TRestAxionOptics::FindFocal.
void PrintMasks()
Prints on screen the 3-optical masks used on the optics planes.
TVector3 fEntrancePosition
The particle position at the optics plane entrance.
TPad * DrawParticleTracks(Double_t deviation=0, Int_t particles=10)
A method to draw an optics schematic including the mirrors geometry, and few photon tracks....
TVector3 GetMiddleDirection()
Returns the middle position from the latest propagated photon.
std::string fOpticsFile
An optics file that contains all the specific optics parameters.
TVector3 GetExitDirection()
Returns the exit position from the latest propagated photon.
void InitFromConfigFile()
Initialization of TRestAxionOptics field members through a RML file.
TPad * DrawDensityMaps(Double_t z, Double_t energy=0, Double_t deviation=0, Int_t particles=1000, Double_t focalHint=7500)
It implements a generic method to identify the optimum focal point. It can be reimplemented at each s...
void PrintMiddleMask()
Prints on screen the mask used on the middle optics plane.
TVector3 GetLastGoodDirection()
It returns the last valid particle direction known in the particle tracking.
TVector3 fSecondInteractionPosition
The particle position at the back mirror interaction point.
TRestCombinedMask *const & GetMiddleMask() const
Returns a pointer to access directly the middle mask information.
TRestAxionOpticsMirror *const & GetMirrorProperties() const
Returns a pointer to access directly the exit mask information.
virtual std::pair< Double_t, Double_t > GetRadialLimits()=0
It returns the lower/higher radius range where photons are allowed.
TVector3 GetEntranceDirection()
Returns the entrance position from the latest propagated photon.
TRestCombinedMask * fEntranceMask
The entrance optical mask that defines a pattern with regions identified with a number.
TRestCombinedMask *const & GetExitMask() const
Returns a pointer to access directly the exit mask information.
virtual void SetMirror()=0
It must be implemented at the inherited optics, making use of fEntrancePosition.
virtual Double_t GetExitPositionZ()=0
It returns the exit Z-position defined by the optical axis.
TVector3 GetMiddlePosition()
Returns the middle position from the latest propagated photon.
void ResetPositions()
It reinitializes particle positions and directions at the different optical regions.
void PrintMirror()
Prints on screen the 3-optical masks used on the optics planes.
TVector3 GetExitPosition()
Returns the exit position from the latest propagated photon.
TRandom3 * fRandom
Random number generator.
void PrintExitMask()
Prints on screen the mask used on the exit optics plane.
Int_t TransportToExit(const TVector3 &pos, const TVector3 &dir)
It moves the incoming particle to the exit of the optics plane and returns the region number where th...
TVector3 GetEntrancePosition()
Returns the entrance position from the latest propagated photon.
virtual Int_t FirstMirrorReflection(const TVector3 &pos, const TVector3 &dir)=0
It updates the values fFirstInteractionPosition and fMiddleDirection. Returns 0 if is not in region.
TRestCombinedMask *const & GetEntranceMask() const
Returns a pointer to access directly the entrance mask information.
virtual TPad * DrawMirrors()=0
It draws the mirrors using a TGraph. To be implemented at the inherited class.
TPad * DrawScatterMaps(Double_t z, Double_t energy=0, Double_t deviation=0, Int_t particles=1000, Double_t focalHint=7500)
It implements a generic method to identify the optimum focal point. It can be reimplemented at each s...
Int_t GetMirror()
It returns the mirror index to be used in the photon reflection.
~TRestAxionOptics()
Default destructor.
TVector3 fOriginPosition
The particle position at the origin.
Double_t CalculateSpotSize(Double_t energy, Double_t z, Int_t particles=15000)
It measures the spot size through Monte Carlo at a given plane given by z. If z=0 this method will ch...
Int_t GetNumberOfReflections()
It returns the total number of reflections. Considering maximum 1-reflection per mirror.
TVector3 fFirstInteractionPosition
The particle position at the front mirror interaction point.
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionOptics.
TVector3 fEntranceDirection
The particle position at the optics plane entrance.
A class used to define and generate a combined structure mask.
A base class for any REST metadata class.
Definition: TRestMetadata.h:74