REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionOpticsMirror.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 _TRestAxionOpticsMirror
24#define _TRestAxionOpticsMirror
25
26#include <TCanvas.h>
27#include <TRestMetadata.h>
28
29#include <iostream>
30
33 private:
35 std::string fMirrorType = "Single"; //<
36
38 std::string fLayerTop = "C"; //<
39
41 std::string fLayerThicknessTop = "30"; //<
42
44 std::string fSigmaTop = "0"; //<
45
47 std::string fLayerBottom = ""; //<
48
50 std::string fLayerThicknessBottom = ""; //<
51
53 std::string fSigmaBottom = ""; //<
54
56 std::string fSubstrate = "SiO2"; //<
57
59 std::map<std::string, std::string> fHenkeKeys;
60
62 std::vector<std::vector<Float_t> > fReflectivityTable;
63
65 std::vector<std::vector<Float_t> > fTransmissionTable;
66
68 TCanvas* fCanvas = nullptr;
69
70 std::string DownloadHenkeFile();
71
72 Int_t ExportTables();
73
74 std::string GetReflectivityFilename();
75 std::string GetTransmissionFilename();
76
77 public:
78 void Initialize();
79
80 void SetMirrorType(const std::string& type) { fMirrorType = type; }
81 void SetLayerMaterial(const std::string& layer) { fLayerTop = layer; }
82 void SetLayerThickness(const std::string& thickness) { fLayerThicknessTop = thickness; }
83 void SetLayerRoughness(const std::string& roughness) { fSigmaTop = roughness; }
84 void SetTopLayerMaterial(const std::string& layer) { fLayerTop = layer; }
85 void SetTopLayerThickness(const std::string& thickness) { fLayerThicknessTop = thickness; }
86 void SetTopLayerRoughness(const std::string& roughness) { fSigmaTop = roughness; }
87 void SetBottomLayerMaterial(const std::string& layer) { fLayerBottom = layer; }
88 void SetBottomLayerThickness(const std::string& thickness) { fLayerThicknessBottom = thickness; }
89 void SetBottomLayerRoughness(const std::string& roughness) { fSigmaBottom = roughness; }
90 void SetSubstrateMaterial(const std::string& substrate) { fSubstrate = substrate; }
91
92 void LoadTables();
93
94 Double_t GetReflectivity(const Double_t angle, const Double_t energy);
95
96 Double_t GetTransmission(const Double_t angle, const Double_t energy);
97
98 void PrintMetadata();
99
100 TCanvas* DrawOpticsProperties(std::string options = "", Double_t lowRange = 1.e-9,
101 Double_t lowRange2 = 1.e-3);
102
104 TRestAxionOpticsMirror(const char* cfgFileName, std::string name = "");
106
107 ClassDef(TRestAxionOpticsMirror, 1);
108};
109#endif
A metadata class accessing the Henke database to load reflectivity data.
std::string fSigmaTop
Layer surface roughness in nm.
std::string fLayerThicknessBottom
The bottom layer thickness in nm. Only used for "Bilayer" mirror type.
std::string fSubstrate
The substrate material.
std::string GetTransmissionFilename()
It returns the corresponding transmission filename for the mirror properties defined in the data memb...
std::string fSigmaBottom
Bottom layer surface roughness in nm. Only used for "Bilayer" mirror type.
Double_t GetReflectivity(const Double_t angle, const Double_t energy)
It returns the interpolated reflectivity for a given angle (in degrees) and a given energy (in keV).
TCanvas * DrawOpticsProperties(std::string options="", Double_t lowRange=1.e-9, Double_t lowRange2=1.e-3)
A method that creates a canvas where the mirror optics properties are drawn. It generates two plots,...
std::string fLayerThicknessTop
The layer thickness in nm.
void LoadTables()
Loads the reflectivity table.
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionOpticsMirror.
std::vector< std::vector< Float_t > > fReflectivityTable
The reflectivity loaded as a table with angle versus energy.
std::map< std::string, std::string > fHenkeKeys
A set of key-value pairs sent to the Henke website for data request.
Int_t ExportTables()
It is a private method to export the tables to a binary file once the tables have been downloaded fro...
TRestAxionOpticsMirror()
Default constructor.
void Initialize()
Initialization of TRestAxionOpticsMirror members.
std::string fLayerBottom
The mirror bottom layer material (chemical formula). Only used for "Bilayer" mirror type.
TCanvas * fCanvas
A canvas to insert the optic properties drawing.
~TRestAxionOpticsMirror()
Default destructor.
Double_t GetTransmission(const Double_t angle, const Double_t energy)
It returns the interpolated transmission for a given angle (in degrees) and a given energy (in keV).
std::string fLayerTop
The mirror layer material (chemical formula).
std::string GetReflectivityFilename()
It returns the corresponding reflectivity filename for the mirror properties defined in the data memb...
std::string fMirrorType
The mirror type (Thick, Single, Bilayer, Multilayer). Only Single/Bilayer is supported now.
std::string DownloadHenkeFile()
It downloads the reflectivity file for the present mirror properties defined at the metadata members.
std::vector< std::vector< Float_t > > fTransmissionTable
The transmission loaded as a table with angle versus energy.
A base class for any REST metadata class.
Definition: TRestMetadata.h:70