REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionTrueWolterOptics.cxx
1/******************** REST disclaimer ***********************************
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
113
114#include "TRestAxionTrueWolterOptics.h"
115
116using namespace std;
117#include <TAxis.h>
118#include <TGraph.h>
119#include <TH1F.h>
120
121#include <cmath>
122
123#include "TRestPhysics.h"
124using namespace REST_Physics;
126
131
136
152 : TRestAxionOptics(cfgFileName) {
153 RESTDebug << "Entering TRestAxionTrueWolterOptics constructor( cfgFileName, name )" << RESTendl;
154
155 Initialize();
156
158
160}
161
167
168 SetSectionName(this->ClassName());
169 SetLibraryVersion(LIBRARY_VERSION);
170
171 fR1 = GetR1();
172 fR2 = GetR2();
173 fR3 = GetR3();
174 fR4 = GetR4();
175 fR5 = GetR5();
176 fAlpha = GetAlpha();
178
179 fXSep.clear();
180 for (unsigned int n = 0; n < fAlpha.size(); n++)
181 if (fR2[n] == fR3[n]) {
182 fXSep.push_back(0);
183 } else {
184 fXSep.push_back(2 * (fR1[n] - fR3[n] - fMirrorLength * TMath::Sin(fAlpha[n])) /
185 TMath::Tan(fAlpha[n]));
186 }
187
188 if (fAlpha.size() == 0) return;
189
190 fCosAlpha.clear();
191 for (const auto& a : fAlpha) fCosAlpha.push_back(TMath::Cos(a));
192
193 fCosAlpha_3.clear();
194 for (const auto& a : fAlpha) fCosAlpha_3.push_back(TMath::Cos(3 * a));
195
196 fFrontVertex.clear();
197 for (unsigned int n = 0; n < fAlpha.size(); n++) fFrontVertex.push_back(fR3[n] / TMath::Tan(fAlpha[n]));
198
199 fBackVertex.clear();
200 for (unsigned int n = 0; n < fAlpha.size(); n++)
201 fBackVertex.push_back(fR3[n] / TMath::Tan(3. * fAlpha[n]));
202
204 if (fEntranceRingsMask) {
205 delete fEntranceRingsMask;
206 fEntranceRingsMask = nullptr;
207 }
209
210 std::vector<Double_t> inner, outer;
211 for (unsigned int n = 0; n < fR1.size() - 1; n++) {
212 inner.push_back(fR1[n] + fThickness[n]);
213 outer.push_back(fR1[n + 1]);
214 }
215 fEntranceRingsMask->SetRadii(inner, outer);
216
217 fEntranceMask->AddMask(fSpiderMask);
219
221 if (fMiddleRingsMask) {
222 delete fMiddleRingsMask;
223 fMiddleRingsMask = nullptr;
224 }
226
227 inner.clear();
228 outer.clear();
229 for (unsigned int n = 0; n < fR3.size() - 1; n++) {
230 inner.push_back(fR3[n] + fThickness[n]);
231 outer.push_back(fR3[n + 1]);
232 }
233 fMiddleRingsMask->SetRadii(inner, outer);
234
235 fMiddleMask->AddMask(fSpiderMask);
237
239 if (fExitRingsMask) {
240 delete fExitRingsMask;
241 fExitRingsMask = nullptr;
242 }
244
245 inner.clear();
246 outer.clear();
247 for (unsigned int n = 0; n < fR5.size() - 1; n++) {
248 inner.push_back(fR5[n] + fThickness[n]);
249 outer.push_back(fR5[n + 1]);
250 }
251 fExitRingsMask->SetRadii(inner, outer);
252
253 fExitMask->AddMask(fSpiderMask);
254 fExitMask->AddMask(fExitRingsMask);
255
256 if (fRandom != nullptr) {
257 delete fRandom;
258 fRandom = nullptr;
259 }
260 fRandom = new TRandom3(0);
261}
262
267Int_t TRestAxionTrueWolterOptics::FirstMirrorReflection(const TVector3& pos, const TVector3& dir) {
268 Int_t mirror = GetMirror();
269 RESTDebug << "--> Entering TRestAxionTrueWolterOptics::FirstMirrorReflection" << RESTendl;
270 RESTDebug << "Mirror: " << mirror << RESTendl;
271 if (mirror < 0) {
272 RESTError << "TRestAxionTrueWolterOptics::FirstMirrorReflection. Mirror index cannot be negative!"
273 << RESTendl;
274 return -1;
275 }
276
277 if (mirror >= 0 && (unsigned int)mirror >= fFrontVertex.size()) {
278 RESTError
279 << "TRestAxionTrueWolterOptics::FirstMirrorReflection. Mirror index above number of mirrors!"
280 << RESTendl;
281 return -1;
282 }
283
284 RESTDebug << "Vertex Z: " << fFrontVertex[mirror] << RESTendl;
285 RESTDebug << "Cos Alpha: " << fCosAlpha[mirror] << RESTendl;
286 TVector3 vertex(0, 0, fFrontVertex[mirror]);
287
290 REST_Physics::GetParabolicVectorIntersection(pos, dir, fAlpha[mirror], fR3[mirror]);
291
293 fFirstInteractionPosition.Z() > -(0.5 * fXSep[mirror])) {
294 RESTDebug << "TRestAxionTrueWolterOptics::FirstMirrorReflection. No interaction!" << RESTendl;
297 fFirstInteraction = false;
298 return 0;
299 }
300
301 TVector3 paraNormal =
303 RESTDebug << "Parabolic normal: (" << paraNormal.X() << ", " << paraNormal.Y() << ", " << paraNormal.Z()
304 << ")" << RESTendl;
305
307
308 RESTDebug << "<-- Exiting TRestAxionTrueWolterOptics::FirstMirrorReflection" << RESTendl;
309
310 fFirstInteraction = true;
311 return 1;
312}
313
318Int_t TRestAxionTrueWolterOptics::SecondMirrorReflection(const TVector3& pos, const TVector3& dir) {
319 Int_t mirror = GetMirror();
320 if (mirror < 0) {
321 RESTError << "TRestAxionTrueWolterOptics::FirstMirrorReflection. Mirror index cannot be negative!"
322 << RESTendl;
323 return 0;
324 }
325
326 if (mirror >= 0 && (unsigned int)mirror >= fFrontVertex.size()) {
327 RESTError
328 << "TRestAxionTrueWolterOptics::FirstMirrorReflection. Mirror index above number of mirrors!"
329 << RESTendl;
330 return 0;
331 }
332
333 TVector3 vertex(0, 0, fBackVertex[mirror]);
334 Double_t focal = fR3[mirror] / TMath::Tan(4 * fAlpha[mirror]);
335 Double_t beta = 3 * fAlpha[mirror];
336
339 REST_Physics::GetHyperbolicVectorIntersection(pos, dir, beta, fR3[mirror], focal);
340
342 fSecondInteractionPosition.Z() < (0.5 * fXSep[mirror])) {
343 RESTDebug << "TRestAxionTrueWolterOptics::SecondMirrorReflection. No interaction!" << RESTendl;
346 fSecondInteraction = false;
347 return 0;
348 }
349
350 TVector3 hyperNormal =
352
354
355 fSecondInteraction = true;
356 return 1;
357}
358
363 if (fSpiderMask) {
364 delete fSpiderMask;
365 fSpiderMask = nullptr;
366 }
367
368 fSpiderMask = (TRestSpiderMask*)this->InstantiateChildMetadata("TRestSpiderMask");
369
370 if (fSpiderMask == nullptr) {
371 RESTWarning << "TRestAxionTrueWolterOptics requires usually a TRestSpiderMask definition" << RESTendl;
372 } else {
373 }
374
376
377 // If we recover the metadata class from ROOT file we will need to call Initialize ourselves
378 this->Initialize();
379}
380
388 if (fR3.size() > 0) {
389 for (unsigned int n = 0; n < fR3.size(); n++) {
390 Double_t dR1 = fR1[n] - fR3[n] - fMirrorLength * TMath::Sin(fAlpha[n]);
391 Double_t dR5 = fR5[n] - fR3[n] + fMirrorLength * TMath::Sin(3 * fAlpha[n]);
392 if (n % 10 == 0)
393 std::cout << "## R1\tdelta R1\tR3\tR5\tdelta R5\talpha\tCosAlpha\tFrontVertex\tBackVertex"
394 << std::endl;
395
396 std::cout << fR1[n] << "\t" << dR1 << "\t" << fR3[n] << "\t" << fR5[n] << "\t" << dR5 << "\t"
397 << fAlpha[n] << "\t" << fCosAlpha[n] << "\t" << fFrontVertex[n] << "\t"
398 << fBackVertex[n] << std::endl;
399 }
400 }
401}
402
408}
409
414
420
421 fPad->cd();
422 std::vector<TGraph*> graphCollection;
423 for (unsigned int mirror = 0; mirror < fR3.size(); mirror++) {
424 TGraph* gr = new TGraph(); //"Mirror" + IntegerToString(mirror + 1));
425
426 Double_t lX = fMirrorLength * fCosAlpha[mirror];
427
428 gr->SetPoint(0, -lX, fR1[mirror]);
429 gr->SetPoint(1, 0, fR3[mirror]);
430 gr->SetPoint(2, lX, fR5[mirror]);
431
432 gr->GetXaxis()->SetLimits(-3.5 * lX, 3.5 * lX);
433 gr->GetHistogram()->SetMaximum(fR1.back() * 1.15);
434 gr->GetHistogram()->SetMinimum(fR1.front() * 0.8);
435
436 gr->GetXaxis()->SetTitle("Z [mm]");
437 gr->GetXaxis()->SetTitleSize(0.04);
438 gr->GetXaxis()->SetLabelSize(0.04);
439 gr->GetXaxis()->SetNdivisions(5);
440 gr->GetYaxis()->SetTitle("R [mm]");
441 gr->GetYaxis()->SetTitleOffset(1.4);
442 gr->GetYaxis()->SetTitleSize(0.04);
443 gr->GetYaxis()->SetLabelSize(0.04);
444 gr->SetLineWidth(6 * fThickness[mirror]);
445 gr->SetLineColor(20 + mirror % 20);
446 if (mirror == 0)
447 gr->Draw("AL");
448 else
449 gr->Draw("L");
450 }
451
452 return fPad;
453}
An abstract class to define common optics parameters and methods.
Bool_t fSecondInteraction
During the photon propagation it tells us if the photon interacted in the second mirror.
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.
TVector3 fMiddleDirection
The particle position at the optics plane middle.
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.
TRestCombinedMask * fExitMask
The exit optical mask that defines a pattern with regions identified with a number.
virtual void Initialize()
Initialization of TRestAxionOptics members.
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...
void InitFromConfigFile()
Initialization of TRestAxionOptics field members through a RML file.
TVector3 fSecondInteractionPosition
The particle position at the back mirror interaction point.
TRestCombinedMask * fEntranceMask
The entrance optical mask that defines a pattern with regions identified with a number.
TRandom3 * fRandom
Random number generator.
Int_t GetMirror()
It returns the mirror index to be used in the photon reflection.
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 that calculates the reflection path of X-rays through a Wolter 1 telescope.
std::vector< Double_t > fR5
Radius R5 in mm. See schematic figure.
Int_t SecondMirrorReflection(const TVector3 &pos, const TVector3 &dir) override
Implementation of first mirror interaction. It updates fSecondInteractionPosition and fExitDirection ...
Double_t GetEntrancePositionZ() override
It returns the entrance Z-position defined by the optical axis.
std::vector< Double_t > fThickness
Mirror thickness in mm. See schematic figure.
TRestSpiderMask * fSpiderMask
The spider structure to be used as an optical opaque mask (common to all planes)
std::vector< Double_t > fR2
Radius R2 in mm. See schematic figure.
std::vector< Double_t > fR1
Entrance radius R1 in mm. See schematic figure.
std::vector< Double_t > GetThickness()
It returns a vector with the values of mirror thickness.
void Initialize() override
Initialization of TRestAxionTrueWolterOptics members.
std::vector< Double_t > fAlpha
Mirror angle (alpha) in radians. See schematic figure.
std::vector< Double_t > GetR5()
It returns a vector with the values of R5.
std::vector< Double_t > GetR3()
It returns a vector with the values of R3.
std::vector< Double_t > fR3
Radius R3 in mm. See schematic figure.
~TRestAxionTrueWolterOptics()
Default destructor.
std::vector< Double_t > fCosAlpha_3
Mirror pre-calculated cosine angle (alpha). See schematic figure.
Int_t FirstMirrorReflection(const TVector3 &pos, const TVector3 &dir) override
Implementation of first mirror interaction. It updates fFirstInteractionPosition and fMiddleDirection...
void InitFromConfigFile() override
Initialization of TRestAxionTrueWolterOptics field members through a RML file.
std::vector< Double_t > fR4
Radius R4 in mm. See schematic figure.
void PrintParameters()
It prints out the Wolter (relevant) parameters extracted from the optics data file,...
void PrintMetadata() override
Prints on screen the information about the metadata members of TRestAxionTrueWolterOptics.
TRestRingsMask * fEntranceRingsMask
The rings structure to be used at entrance as an optical opaque mask.
Double_t GetExitPositionZ() override
It returns the exit Z-position defined by the optical axis.
std::vector< Double_t > GetR2()
It returns a vector with the values of R2.
void PrintSpider()
It prints out the spider mask common to all the optical planes.
std::vector< Double_t > GetR4()
It returns a vector with the values of R4.
TRestAxionTrueWolterOptics()
Default constructor.
std::vector< Double_t > fFrontVertex
The Z-position of the cone vertex defined by the front mirrors.
TRestRingsMask * fExitRingsMask
The rings structure to be used at entrance as an optical opaque mask.
std::vector< Double_t > fBackVertex
The Z-position of the cone vertex defined by the back mirrors.
std::vector< Double_t > GetAlpha()
It returns a vector with the values of alpha.
TRestRingsMask * fMiddleRingsMask
The rings structure to be used at entrance as an optical opaque mask.
TPad * DrawMirrors() override
A method to to draw an optics schematic including the mirrors geometry.
std::vector< Double_t > GetR1()
It returns a vector with the values of R1.
std::vector< Double_t > fCosAlpha
Mirror pre-calculated cosine angle (alpha). See schematic figure.
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.
TRestMetadata * InstantiateChildMetadata(int index, std::string pattern="")
This method will retrieve a new TRestMetadata instance of a child element of the present TRestMetadat...
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
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.
A class used to define a rings mask pattern.
void SetRadii(const std::vector< Double_t > &innerR, const std::vector< Double_t > &outterR)
It allows to redefine the inner and outter rings radii directly.
A class used to define and generate a spider structure mask.
void PrintMetadata() override
Prints on screen the complete information about the metadata members from this class.
@ REST_Info
+show most of the information for each steps
This namespace serves to define physics constants and other basic physical operations.
Definition: TRestPhysics.h:34
TVector3 GetHyperbolicVectorIntersection(const TVector3 &pos, const TVector3 &dir, const Double_t beta, const Double_t R3, const Double_t focal)
TVector3 MoveByDistance(const TVector3 &pos, const TVector3 &dir, Double_t d)
This method transports a position pos by a distance d in the direction defined by dir.
TVector3 GetParabolicVectorIntersection(const TVector3 &pos, const TVector3 &dir, const Double_t alpha, const Double_t R3)
TVector3 GetParabolicNormal(const TVector3 &pos, const Double_t alpha, const Double_t R3)
This method returns the normal vector on a parabolic surface pointing towards the inside of the parab...
TVector3 GetHyperbolicNormal(const TVector3 &pos, const Double_t beta, const Double_t R3, const Double_t focal)
This method returns the normal vector on a hyperbolic surface pointing towards the inside of the hype...
TVector3 GetVectorReflection(const TVector3 &dir, const TVector3 &n)
This method will return the reflected vector respect to a plane defined by its normal vector n....