source: trunk/Mars/msimreflector/MFresnelLens.cc@ 19924

Last change on this file since 19924 was 19753, checked in by tbretz, 5 years ago
Allow setting SurfaceRoughness (note that this is DANGEROUS to use) and fixed some compiler warnings.
File size: 53.2 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of CheObs, the Modular Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appears in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 6/2019 <mailto:tbretz@physik.rwth-aachen.de>
19!
20! Copyright: CheObs Software Development, 2000-2019
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MFresnelLens
28//
29// For some details on definitions please refer to
30// https://application.wiley-vch.de/berlin/journals/op/07-04/OP0704_S52_S55.pdf
31//
32// The HAWC's Eye lens is an Orafol SC943
33// https://www.orafol.com/en/europe/products/optic-solutions/productlines#pl1
34//
35// A good description on ray-tracing can be found here
36// https://graphics.stanford.edu/courses/cs148-10-summer/docs/2006--degreve--reflection_refraction.pdf
37//
38//////////////////////////////////////////////////////////////////////////////
39#include "MFresnelLens.h"
40
41#include <fstream>
42#include <errno.h>
43
44#include <TRandom.h>
45
46#include "MQuaternion.h"
47#include "MReflection.h"
48
49#include "MMath.h"
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54ClassImp(MFresnelLens);
55
56using namespace std;
57
58// ==========================================================================
59
60enum exception_t
61{
62 kValidRay = 0,
63
64 kStrayUpgoing,
65 kOutsideRadius,
66 kNoSurfaceFound,
67 kStrayDowngoing,
68 kAbsorbed,
69
70 kFoundSurfaceUnavailable,
71
72 kInvalidOrigin,
73 kTransitionError,
74
75 kEnter = 1000,
76 kLeave = 2000,
77};
78
79enum surface_t
80{
81 kPhotonHasLeft = 0,
82
83 kEntrySurface,
84 kSlopeSurface,
85 kDraftSurface,
86 kExitSurface,
87
88 kMaterial = 5,
89
90 kNoSurface = 9
91};
92
93
94class raytrace_exception : public runtime_error
95{
96protected:
97 int fError;
98 int fOrigin;
99 int fSurface;
100
101public:
102 raytrace_exception(const int &_id, const int &_origin, const int &_surface, const string& what_arg) :
103 runtime_error(what_arg), fError(_id), fOrigin(_origin), fSurface(_surface)
104 {
105 }
106
107 raytrace_exception(const int &_id, const int &_origin, const int &_surface, const char* what_arg) :
108 runtime_error(what_arg), fError(_id), fOrigin(_origin), fSurface(_surface)
109 {
110 }
111
112 int id() const { return fError + fSurface*10 + fOrigin*100; }
113 int error() const { return fError; }
114 int origin() const { return fOrigin; }
115 int surface() const { return fSurface; }
116};
117
118class raytrace_error : public raytrace_exception
119{
120public:
121 raytrace_error(const int &_id, const int &_origin, const int &_surface, const string& what_arg) :
122 raytrace_exception(_id, _origin, _surface, what_arg) { }
123 raytrace_error(const int &_id, const int &_origin, const int &_surface, const char* what_arg) :
124 raytrace_exception(_id, _origin, _surface, what_arg) { }
125};
126class raytrace_info : public raytrace_exception
127{
128public:
129 raytrace_info(const int &_id, const int &_origin, const int &_surface, const string& what_arg) :
130 raytrace_exception(_id, _origin, _surface, what_arg) { }
131 raytrace_info(const int &_id, const int &_origin, const int &_surface, const char* what_arg) :
132 raytrace_exception(_id, _origin, _surface, what_arg) { }
133};
134class raytrace_user : public raytrace_exception
135{
136public:
137 raytrace_user(const int &_id, const int &_origin, const int &_surface, const string& what_arg) :
138 raytrace_exception(_id, _origin, _surface, what_arg) { }
139 raytrace_user(const int &_id, const int &_origin, const int &_surface, const char* what_arg) :
140 raytrace_exception(_id, _origin, _surface, what_arg) { }
141};
142
143// ==========================================================================
144
145
146// --------------------------------------------------------------------------
147//
148// Default constructor
149//
150MFresnelLens::MFresnelLens(const char *name, const char *title) :
151 fPSF(0), fSlopeAbsorption(false), fDraftAbsorption(false),
152 fBottomReflection(true), fDisableMultiEntry(false), fFresnelReflection(true),
153 fMinHits(0), fMaxHits(0)
154{
155 fName = name ? name : "MFresnelLens";
156 fTitle = title ? title : "Parameter container storing a collection of several mirrors (reflector)";
157
158 // Default: Orafol SC943
159
160 DefineLens();
161}
162
163// ==========================================================================
164
165// --------------------------------------------------------------------------
166//
167// Default ORAFOL SC943
168//
169// Focal Length: F = 50.21 cm
170// Diameter: D = 54.92 cm
171// Groove width: w = 0.01 cm
172// Lens thickness: h = 0.25 cm
173//
174// Default wavelength: 546 nm
175//
176void MFresnelLens::DefineLens(double F, double D, double w, double h, double lambda)
177{
178 fR = D/2; // [cm] Lens radius
179 fW = w; // [cm] Width of a single groove
180 fH = h; // [cm] Thickness of lens
181 fF = F; // [cm] focal length (see also MGeomCamFAMOUS!)
182
183 fLambda = lambda;
184
185 fN = MFresnelLens::RefractiveIndex(fLambda); // Lens
186
187 // Velocity of light within the lens material [cm/ns]
188 // FIXME: Note that for the correct conversion in Transmission()
189 // also the speed in the surrounding medium has to be taken correctly
190 // into account (here it is assumed to be air with N=1
191 fVc = fN/(TMath::C()*100/1e9); // cm/ns
192
193 InitGeometry(fR, fW, fN, fF, fH);
194}
195
196// --------------------------------------------------------------------------
197//
198// Precalculate values such as the intersection points inside the grooves,
199// the angle of the slope and draft surface and the corresponding tangents.
200//
201void MFresnelLens::InitGeometry(double maxr, double width, double N0, double F, double d)
202{
203 const uint32_t num = TMath::CeilNint(maxr/width);
204
205 fGrooves.resize(num);
206
207 for (uint32_t i=0; i<num; i++)
208 {
209 const double r0 = i*width;
210 const double rc = i*width + width/2;
211 const double r1 = i*width + width;
212
213 // Slope angle of the reflecting surface alpha
214 // Angle of the draft surface psi
215 const double alpha = -MFresnelLens::SlopeAngle(rc, F, N0, d); // w.r.t. x [30]
216 const double psi = MFresnelLens::DraftAngle(r1); // w.r.t. z [ 5]
217
218 const double tan_alpha = tan(alpha);
219 const double tan_psi = tan(psi);
220
221 fGrooves[i].slope.z = r0*tan_alpha;
222 fGrooves[i].draft.z = -r1/tan_psi;
223
224 fGrooves[i].slope.theta = TMath::Pi()/2-alpha; // w.r.t. +z [ 60]
225 fGrooves[i].draft.theta = -psi; // w.r.t. +z [- 5]
226
227 fGrooves[i].slope.tan_theta = tan(fGrooves[i].slope.theta);
228 fGrooves[i].draft.tan_theta = tan(fGrooves[i].draft.theta);
229
230 fGrooves[i].slope.tan_theta2 = fGrooves[i].slope.tan_theta*fGrooves[i].slope.tan_theta;
231 fGrooves[i].draft.tan_theta2 = fGrooves[i].draft.tan_theta*fGrooves[i].draft.tan_theta;
232
233 fGrooves[i].slope.theta_norm = TMath::Pi()/2-fGrooves[i].slope.theta; // [ 30]
234 fGrooves[i].draft.theta_norm = TMath::Pi()/2-fGrooves[i].draft.theta; // [ 95]
235
236 const double dr = width/(tan_alpha*tan_psi+1);
237
238 fGrooves[i].r = r0 + dr;
239
240 const double z = -dr*tan_alpha;
241
242 fGrooves[i].slope.h = z;
243 fGrooves[i].draft.h = z;
244
245 if (z<-fH)
246 *fLog << warn << "Groove " << i << " deeper (" << z << ") than thickness of lens material (" << fH << ")." << endl;
247 }
248
249 fMaxR = (num+1)*width;
250}
251
252// --------------------------------------------------------------------------
253//
254// Reads the transmission curve from a file
255// (tranmission in percent versus wavelength in nanometers)
256//
257// The transmission curve is used to calculate the absorption lengths.
258// Therefore the thickness for which the tranission curve is valid is
259// required (in cm).
260//
261// The conversion can correct for fresnel reflection at the entry and exit
262// surface assuming that the outside material during the measurement was air
263// (n=1.0003) and the material in PMMA. Correction is applied when
264// correction is set to true <default>.
265//
266// If no valid data was read, 0 is returned. -1 is returned if any tranmission
267// value read from the file is >1. If the fresnel correction leads to a value >1,
268// the value is set to 1. The number of valid data points is returned.
269//
270Int_t MFresnelLens::ReadTransmission(const TString &file, float thickness, bool correction)
271{
272 TGraph transmission(file);
273
274 /*
275 double gx_min, gx_max, gy_min, gy_max;
276 absorption.ComputeRange(gx_min, gy_min, gx_max, gy_max);
277 if (lambda<gx_min || lambda>gx_max)
278 {
279 cout << "Invalid wavelength" << endl;
280 return;
281 }*/
282
283 if (transmission.GetN()==0)
284 return 0;
285
286 for (int i=0; i<transmission.GetN(); i++)
287 {
288 // Correct transmission for Fresnel reflection on the surface
289 const double lambda = transmission.GetX()[i];;
290
291 double trans = transmission.GetY()[i];
292 if (trans>1)
293 {
294 *fLog << err << "Transmission larger than 1." << endl;
295 return -1;
296 }
297
298 if (correction)
299 {
300 // Something like this is requried if correction
301 // for optical boundaries is necessary
302 const double n0 = MFresnelLens::RefractiveIndex(lambda);
303
304 // FIXME: Make N_air a variable
305 const double r0 = (n0-1.0003)/(n0+1.0003);
306 const double r2 = r0*r0;
307
308 trans *= (1+r2)*(1+r2);
309
310 if (trans>1)
311 {
312 *fLog << warn << "Transmission at " << lambda << "nm (" << trans << ") after Fresnel correction larger than 1." << endl;
313 trans = 1;
314 }
315 }
316
317 // convert to absorption length (FIMXE: Sanity check)
318 transmission.GetY()[i] = -thickness/log(trans>0.999 ? 0.999 : trans);
319 }
320
321 fAbsorptionLength = MSpline3(transmission);
322
323 return fAbsorptionLength.GetNp();
324}
325
326Int_t MFresnelLens::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
327{
328 Bool_t rc = kFALSE;
329
330 if (IsEnvDefined(env, prefix, "SurfaceRoughness", print))
331 {
332 rc = kTRUE;
333 if (!GetEnvValue(env, prefix, "SurfaceRoughness", fPSF))
334 return kERROR;
335 }
336
337 const int correction = GetEnvValue(env, prefix, "Transmission.FresnelCorrection", -1);
338 const float thickness = GetEnvValue(env, prefix, "Transmission.Thickness", -1.0); // [cm]
339 const TString fname = GetEnvValue(env, prefix, "Transmission.FileName", "");
340
341 const bool correction_valid = correction>=0;
342 const bool thickness_valid = thickness>0;
343 const bool fname_valid = !fname.IsNull();
344
345 if (!correction_valid && !thickness_valid && !fname_valid)
346 return rc;
347
348 if (correction_valid && thickness_valid && fname_valid)
349 return ReadTransmission(fname, thickness, correction) >= 0 || rc;
350
351 *fLog << err << "Reading transmission file required FileName, Thickness and FresnelCorrection." << endl;
352 return kERROR;
353}
354
355// ==========================================================================
356
357// --------------------------------------------------------------------------
358//
359// Refractive Index of PMMA, according to
360// https://refractiveindex.info/?shelf=organic&book=poly(methyl_methacrylate)&page=Szczurowski
361//
362// n^2-1=\frac{0.99654 l^2}{l^2-0.00787}+\frac{0.18964 l^2}{l^2-0.02191}+\frac{0.00411 l^2}{l^2-3.85727}
363//
364// Returns the refractive index n as a function of wavelength (in nanometers)
365//
366double MFresnelLens::RefractiveIndex(double lambda)
367{
368 const double l2 = lambda*lambda;
369
370 const double c0 = 0.99654/(1-0.00787e6/l2);
371 const double c1 = 0.18964/(1-0.02191e6/l2);
372 const double c2 = 0.00411/(1-3.85727e6/l2);
373
374 return sqrt(1+c0+c1+c2);
375}
376
377// --------------------------------------------------------------------------
378//
379// A Fresnel lens with parabolic surface calculated with the sagittag
380// function (k=-1) and a correction for the thickness of the lens
381// on the curvature. See also PhD thesis, Tim Niggemann ch. 7.1.1.
382//
383// see also W.J.Smith, Modern Optical Engineering, 2.8 The "Thin Lens"
384// 1/f = (n-1)/radius Eq. 2.36 with thickness t = 0
385// bfl = f Eq. 2.37 and R2 = inf (c2 = 0)
386//
387// Parameters are:
388// The distance from the center r
389// The focal length to be achieved F
390// The refractive index of the outer medium (usually air) n0
391// The refractive index of the lens material (e.g. PMMA) n1
392// The thichness of the lens d
393//
394// r, F and d have to be in the same units.
395//
396// Return the slope angle alpha [rad]. The Slope angle is defined with
397// respect to the plane of the lens. (0 at the center, decreasing
398// with increasing radial distance)
399//
400double MFresnelLens::SlopeAngleParabolic(double r, double F, double n0, double n1, double d)
401{
402 // In the datasheet, it looks as if F is calculated
403 // towards the center of the lens. It seems things are more
404 // consistent if the thickness correction in caluating the
405 // slope angle is omitted and the focal distance is measured
406 // from the entrance of the lens => FIXME: To be checked
407 const double rn = n1/n0;
408 const double c = (rn - 1) * (F + d/rn); // FIXME: try and error with a large d
409 return -atan(r/c);
410
411 // F = 50.21
412 // d= 10 d=20
413 // -: 47 43.7
414 // 0: 53.5 57.0
415 // +: 60.3 70.3
416}
417
418// --------------------------------------------------------------------------
419//
420// A Fresnel lens with an optimized parabolic surface calculated with
421// the sagittag function (k=-1) and fitted coefficients according
422// to Master thesis, Eichler.
423//
424// Note that for this setup other parameters must be fixed
425//
426// Parameters are:
427// The distance from the center r
428//
429// r is in cm.
430//
431// Return the slope angle alpha [rad]. The Slope angle is defined with
432// respect to the plane of the lens. (0 at the center, decreasing
433// with increasing radial distance)
434//
435double MFresnelLens::SlopeAngleAspherical(double r)
436{
437 // Master, Eichler [r/cm]
438 return -atan( r/26.47
439 +2*1.18e-4 * 1e1*r
440 +4*1.34e-9 * 1e3*r*r*r
441 +6*9.52e-15 * 1e5*r*r*r*r*r
442 -8*2.04e-19 * 1e7*r*r*r*r*r*r*r);
443}
444
445// --------------------------------------------------------------------------
446//
447// Ideal angle of the Fresnel surfaces at a distance r from the center
448// to achieve a focal distance F for a positive Fresnel lens made
449// from a material with a refractive index n.
450// A positive Fresnel lens is one which focuses light from infinity
451// (the side with the grooves) to a point (the flat side of the lens).
452//
453// The calculation follows
454// https://shodhganga.inflibnet.ac.in/bitstream/10603/131007/13/09_chapter%202.pdf
455// Here, a thin lens is assumed
456//
457// sin(omega) = r / sqrt(r^2+F^2)
458// tan(alpha) = sin(omega) / [ 1 - sqrt(n^2-sin(omega)^2) ]
459//
460// Return alpha [rad] as a function of the radial distance r, the
461// focal length F and the refractive index n. r and F have to have
462// the same units. The Slope angle is defined with respect to the plane
463// of the lens. (0 at the center, decreasing with increasing radial
464// distance)
465//
466double MFresnelLens::SlopeAngleOptimized(double r, double F, double n)
467{
468 // Use F+d/2
469 double so = r / sqrt(r*r + F*F);
470 return atan(so / (1-sqrt(n*n - so*so))); // alpha<0, Range [0deg; -50deg]
471}
472
473// --------------------------------------------------------------------------
474//
475// Currently calles SlopeAngleParabolic(r, F, 1, n, d)
476//
477double MFresnelLens::SlopeAngle(double r, double F, double n, double d)
478{
479 return SlopeAngleParabolic(r, F, 1.0003, n, d);
480}
481
482
483//
484// Draft angle of the Orafol SC943 According to the thesis of Eichler
485// and NiggemannTim Niggemann:
486//
487// The surface of the lens follows the shape of a parabolic lens to compensate spherical aberration
488// Draft angle: psi(r) = 3deg + r * 0.0473deg/mm
489//
490// The draft angle is returned in radians and is defined w.r.t. to the
491// normal of the lens surface. (almost 90deg at the center,
492// decreasing with increasing radial distance)
493//
494double MFresnelLens::DraftAngle(double r)
495{
496 return (3 + r*0.473)*TMath::DegToRad(); // Range [0deg; 15deg]
497}
498
499// ==========================================================================
500
501// --------------------------------------------------------------------------
502//
503// Return the total Area of all mirrors. Note, that it is recalculated
504// with any call.
505//
506Double_t MFresnelLens::GetA() const
507{
508 return fMaxR*fMaxR*TMath::Pi();
509}
510
511// --------------------------------------------------------------------------
512//
513// Check with a rough estimate whether a photon can hit the reflector.
514//
515Bool_t MFresnelLens::CanHit(const MQuaternion &p) const
516{
517 // p is given in the reflectory coordinate frame. This is meant as a
518 // fast check without lengthy calculations to omit all photons which
519 // cannot hit the reflector at all
520 return p.R2()<fMaxR*fMaxR;
521}
522
523// ==========================================================================
524
525// FIXME: The rays could be 'reflected' inside the material
526// (even though its going out) or vice versa
527static double RandomTheta(double psf)
528{
529 return psf>0 ? MMath::RndmPSF(psf)/2 : 0;
530}
531
532// FIXME: The rays could be 'reflected' inside the material
533// (even though its going out) or vice versa
534static double RandomPhi(double r, double psf)
535{
536 return psf>0 ? MMath::RndmPSF(psf)/2 : 0;
537}
538
539
540// --------------------------------------------------------------------------
541//
542// Calculate the intersection point beweteen a line defined by the position p
543// and the direction u and a cone defined by the object cone.
544//
545// Z: position of peak of cone
546// theta: opening angle of cone
547//
548// Distance r of cone surface at given z from z-axis
549// r_cone(z) = (Z-z)*tan(theta)
550//
551// Equalition of line
552// (x) (p.x) (u.x/u.z)
553// (y) = (p.y) + dz * (u.y/u.z)
554// (z) (p.z) ( 1 )
555//
556// Normalization
557// U.x := u.x/u.z
558// U.y := u.y/u.z
559//
560// Distance of line at given z from z-axis
561// r_line(z) = sqrt(x^2 + y^2) = sqrt( (p.x+dz*u.x)^2 + (p.y+dz*u.y)^2) with dz = z-p.z
562//
563// Equation to be solved
564// r_cone(z) = r_line(z)
565//
566// Solved with wxmaxima:
567//
568// [0] solve((px+(z-pz)*Ux)^2+(py+(z-pz)*Uy)^2= ((Z-z)*t)^2, z);
569//
570// z= (sqrt(((Uy^2+Ux^2)*pz^2+(-2*Uy*py-2*Ux*px-2*Z*Uy^2-2*Z*Ux^2)*pz+py^2+2*Z*Uy*py+px^2+2*Z*Ux*px+Z^2*Uy^2+Z^2*Ux^2)*t^2-Ux^2*py^2+2*Ux*Uy*px*py-Uy^2*px^2)+Z*t^2+(-Uy^2-Ux^2)*pz+Uy*py+Ux*px)/(t^2-Uy^2-Ux^2),
571// z=-(sqrt(((Uy^2+Ux^2)*pz^2+(-2*Uy*py-2*Ux*px-2*Z*Uy^2-2*Z*Ux^2)*pz+py^2+2*Z*Uy*py+px^2+2*Z*Ux*px+Z^2*Uy^2+Z^2*Ux^2)*t^2-Ux^2*py^2+2*Ux*Uy*px*py-Uy^2*px^2)-Z*t^2+( Uy^2+Ux^2)*pz-Uy*py-Ux*px)/(t^2-Uy^2-Ux^2)
572//
573double MFresnelLens::CalcIntersection(const MQuaternion &p, const MQuaternion &u, const Cone &cone) const
574{
575 const double &Z = cone.z;
576
577 const double Ux = u.X()/u.Z();
578 const double Uy = u.Y()/u.Z();
579
580 const double px = p.X();
581 const double py = p.Y();
582 const double pz = p.Z();
583
584 //const double &t = cone.tan_theta;
585 const double &t2 = cone.tan_theta2;
586
587 const double Ur2 = Ux*Ux + Uy*Uy;
588 const double pr2 = px*px + py*py;
589 const double Up2 = Ux*px + Uy*py;
590
591 const double cr2 = Ux*py - Uy*px;
592
593 const double a = t2 - Ur2;
594 const double b = Ur2*pz - Up2 - Z*t2;
595
596 const double h = Z-pz;
597 const double h2 = h*h;
598
599 // [ -b +-sqrt(b^2 - 4 ac) ] / [ 2a ]
600
601 const double radix = (Ur2*h2 + 2*Up2*h + pr2)*t2 - cr2*cr2;
602 if (radix<0)
603 return 0;
604
605 const double sqrt_radix = sqrt(radix);
606
607 const double dz[2] =
608 {
609 (+sqrt_radix - b)/a,
610 (-sqrt_radix - b)/a
611 };
612
613 // Return the closest solution inside the allowed range
614 // which is in the direction of movement
615
616 const double &H = cone.h;
617
618 const bool is_inside0 = dz[0]>=H && dz[0]<0;
619 const bool is_inside1 = dz[1]>=H && dz[1]<0;
620
621 // FIXME: Simplify!
622 if (!is_inside0 && !is_inside1)
623 return 0;
624
625 // Only dz[0] is in the right z-range
626 if (is_inside0 && !is_inside1)
627 {
628 // Check if dz[0] is in the right direction
629 if ((u.Z()>=0 && dz[0]>=p.Z()) ||
630 (u.Z()< 0 && dz[0]< p.Z()))
631 return dz[0];
632
633 return 0;
634 }
635
636 // Only dz[1] is in the right z-range
637 if (!is_inside0 && is_inside1)
638 {
639 // Check if dz[1] is in the right direction
640 if ((u.Z()>=0 && dz[1]>=p.Z()) ||
641 (u.Z()< 0 && dz[1]< p.Z()))
642 return dz[1];
643
644 return 0;
645 }
646
647 /*
648 if (is_inside0^is_inside1)
649 {
650 if (u.Z()>=0)
651 return dz[0]>p.Z() ? dz[0] : dz[1];
652 else
653 return dz[0]<p.Z() ? dz[0] : dz[1];
654 }*/
655
656
657 // dz[0] and dz[1] are in the right range
658 // return the surface which is hit first
659
660 // moving upwards
661 if (u.Z()>=0)
662 {
663 // Both solution could be correct
664 if (dz[0]>=p.Z() && dz[1]>=p.Z())
665 return std::min(dz[0], dz[1]);
666
667 // only one solution can be correct
668 return dz[0]>=p.Z() ? dz[0] : dz[1];
669 }
670 else
671 {
672 // Both solution could be correct
673 if (dz[0]<p.Z() && dz[1]<p.Z())
674 return std::max(dz[0], dz[1]);
675
676 // only one solution can be correct
677 return dz[0]<p.Z() ? dz[0] : dz[1];
678 }
679}
680
681// --------------------------------------------------------------------------
682//
683// Find the peak (draft+slope) which will be hit by the photon which
684// is defined by position p and direction u. ix gives the index of the groove
685// to originate the search from.
686//
687// Returns the index of the groove to which the surface belongs, -1 if no
688// matching surface was found.
689//
690int MFresnelLens::FindPeak(size_t ix, const MQuaternion &p, const MQuaternion &u) const
691{
692 // ---------------------------
693 // check for first groove first
694 if (ix==0)
695 {
696 const auto test = p.fVectorPart + (fGrooves[0].slope.h-p.Z())/u.Z()*u.fVectorPart;
697 if (test.XYvector().Mod()<fGrooves[0].r)
698 return 0;
699 }
700
701 // r = sqrt( (px + t*ux) + (py + t*uy)^2 )
702 // dr/dt = (2*uy*(dz*uy+py)+2*ux*(dz*ux+px))/(2*sqrt((dz*uy+py)^2+(dz*ux+px)^2))
703 // dr/dt = (uy*py + ux*px)/sqrt(py^2+px^2)
704 const bool outgoing = u.X()*p.X() + u.Y()*p.Y() > 0; // r is (at least locally) increasing
705
706 // ---------------------------
707 const double Ux = u.X()/u.Z();
708 const double Uy = u.Y()/u.Z();
709
710 const double px = p.X();
711 const double py = p.Y();
712 const double pz = p.Z();
713
714 const double Ur2 = Ux*Ux + Uy*Uy;
715 const double cr2 = Ux*py - Uy*px;
716 const double pr2 = px*px + py*py;
717 const double Up2 = Ux*px + Uy*py;
718
719 //for (int i=1; i<fGrooves.size(); i++)
720
721 // To speed up the search, search first along the radial moving direction of
722 // the photon. If that was not successfull, try in the opposite direction.
723 // FIXME: This could still fail in some very rare cases, for some extremely flat trajectories
724 for (int j=0; j<2; j++)
725 {
726 const bool first = j==0;
727
728 const int step = outgoing ^ !first ? 1 : -1;
729 const int end = outgoing ^ !first ? fGrooves.size() : 1;
730 const int beg = std::max<size_t>(j==0 ? ix : ix+step, 1);
731
732 for (int i=beg; i!=end; i+=step)
733 {
734 const Groove &groove1 = fGrooves[i-1];
735 const Groove &groove2 = fGrooves[i];
736
737 const double &z1 = groove1.draft.h;
738 const double &z2 = groove2.slope.h;
739
740 const double &r1 = groove1.r;
741 const double &r2 = groove2.r;
742
743 Cone cone;
744 cone.tan_theta = -(r2-r1)/(z2-z1);
745 cone.tan_theta2 = cone.tan_theta*cone.tan_theta;
746 cone.z = z1 + r1/cone.tan_theta;
747
748 const double &Z = cone.z;
749 const double &t2 = cone.tan_theta2;
750
751 const double a = t2 - Ur2;
752 const double b = Ur2*pz - Up2 - Z*t2;
753
754 const double h = Z-pz;
755 const double h2 = h*h;
756
757 // [ -b +-sqrt(b^2 - 4 ac) ] / [ 2a ]
758
759 const double radix = (Ur2*h2 + 2*Up2*h + pr2)*t2 - cr2*cr2;
760 if (radix<0)
761 continue;
762
763 const double sqrt_radix = sqrt(radix);
764
765 const double dz[2] =
766 {
767 (+sqrt_radix - b)/a,
768 (-sqrt_radix - b)/a
769 };
770
771 if (dz[0]>=z2 && dz[0]<=z1)
772 return i;
773
774 if (dz[1]>=z2 && dz[1]<=z1)
775 return i;
776 }
777 }
778
779 return -1;
780}
781
782// --------------------------------------------------------------------------
783//
784// If no transmission was given returns true. Otherwaise calculates the
785// absorption length for a flight time dt in the material and a photon
786// with wavelength lambda. The flight time is converted to a geometrical
787// using the speed of light in the medium.
788//
789// Returns true if the poton passed, false if it was absorbed.
790//
791bool MFresnelLens::Transmission(double dt, double lambda) const
792{
793 if (fAbsorptionLength.GetNp()==0)
794 return true;
795
796 // FIXME: Speed up!
797 const double alpha = fAbsorptionLength.Eval(lambda);
798
799 // We only have the travel time, thus we have to convert back to distance
800 // Note that the transmission coefficients are w.r.t. to geometrical
801 // distance not light-travel distance. Thus the distance has to be corrected
802 // for the corresponding refractive index of the material.
803 const double cm = dt/fVc;
804
805 const double trans = exp(-cm/alpha);
806 return gRandom->Uniform()<trans;
807}
808
809/*
810// surface=0 : incoming ray
811// surface=1 : slope
812// surface=2 : draft
813// surface=3 : bottom
814int MFresnelLens::EnterGroove(int surface, double n0, double lambda, MQuaternion &pos, MQuaternion &dir) const
815{
816 const double rx = pos.R();
817
818 if (surface==3)
819 {
820 //cout << "Bottom as origin invalid" << endl;
821 throw -1;
822
823 }
824 if (rx>=fR)
825 {
826 //cout << "Left the lens radius (enter)" << endl;
827 throw -2;
828 }
829 //if (dir.Z()>0)
830 //{
831 // cout << "Upgoing, outside of the material" << endl;
832 // PropagateZ(pos, dir, dir.Z()>0 ? 3 : -3);
833 // return -1;
834 //}
835
836
837 // Calculate the ordinal number of the groove correpsonding to rx
838 const int ix = TMath::FloorNint(rx/fW);
839
840 // Photons was just injected (test both surfaces) or came from the other surface
841 if (surface==0 || surface==2)
842 {
843 // Get the possible intersection point with the slope angle
844 const double z1 = CalcIntersection(pos, dir, fGrooves[ix].slope);
845
846 // We hit the slope angle
847 if (z1!=0)
848 {
849 // Move photon to new hit position
850 pos.PropagateZ(dir, z1);
851
852 if (fSlopeAbsorption)
853 throw -100;
854
855 // Get the normal vector of the surface which was hit
856 const VectorNorm norm(fGrooves[ix].slope.theta_norm+RandomTheta(fPSF),
857 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
858
859 // Get the optical transition of the direction vector
860 const int ret = MOptics::ApplyTransition(dir, norm, 1, n0);
861
862 // Transition was Reflection - try again
863 if (ret==1 || ret==2)
864 return EnterGroove(1, n0, lambda, pos, dir)+1;
865
866 // Transition was Refraction - enter
867 if (ret>=3)
868 return LeavePeak(1, n0, lambda, pos, dir, pos.T())+1;
869
870 // Error occured (see ApplyTransition for details)
871 //cout << "ERR[TIR1]" << endl;
872 throw -3;
873 }
874 }
875
876 // Photons was just injected (test both surfaces) or came from the other surface
877 if (surface==0 || surface==1)
878 {
879 const double z2 = CalcIntersection(pos, dir, fGrooves[ix].draft);
880
881 // We hit the draft angle
882 if (z2!=0)
883 {
884 // Move photon to new hit position
885 pos.PropagateZ(dir, z2);
886
887 if (fDraftAbsorption)
888 throw -101;
889
890 // Get the normal vector of the surface which was hit
891 const VectorNorm norm(fGrooves[ix].draft.theta_norm+RandomTheta(fPSF),
892 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
893
894 // Get the optical transition of the direction vector
895 const int ret = MOptics::ApplyTransition(dir, norm, 1, n0);
896
897 // Transition was Reflection - try again
898 if (ret==1 || ret==2)
899 return EnterGroove(2, n0, lambda, pos, dir)+1;
900
901 // Transition was Refraction - enter
902 if (ret>=3)
903 return -LeavePeak(2, n0, lambda, pos, dir, pos.T())+1;
904
905 // Error occured (see ApplyTransition for details)
906 //cout << "ERR[TIR2]" << endl;
907 throw -4;
908 }
909 }
910
911 if (dir.Z()>0)
912 {
913 //cout << "Upgoing, outside of the material" << endl;
914 //pos.PropagateZ(dir, dir.Z()>0 ? 3 : -3);
915 throw -5;
916 }
917
918 // The ray has left the peak at the bottom(?)
919 //cout << "ERR[N/A]" << endl;
920 throw -6;
921}
922*/
923
924
925// surface=0 : incoming ray
926// surface=1 : slope
927// surface=2 : draft
928// surface=3 : bottom
929int MFresnelLens::EnterGroove(int surface, double n0, MQuaternion &pos, MQuaternion &dir) const
930{
931 const double rx = pos.R();
932
933 if (surface==kExitSurface)
934 throw raytrace_error(kEnter+kInvalidOrigin, surface, -1,
935 "EnterGroove - Bottom as origin invalid");
936
937 if (rx>=fR) // This is an error as the direction vector is now invalid
938 throw raytrace_error(kEnter+kOutsideRadius, surface, -1,
939 "EnterGroove - Surface hit outside allowed radius");
940
941 /*
942 if (dir.Z()>0)
943 return -1;
944 }*/
945
946
947 // FIXME: There is a very tiny chance that a ray hits the same surface twice for
948 // very horizontal rays. Checking this needs to make sure that the same
949 // solution is not just found again.
950
951 // Calculate the ordinal number of the groove correpsonding to rx
952 const int ix = TMath::FloorNint(rx/fW);
953
954 // Photons was just injected (test both surfaces) or came from the other surface
955 if (surface==kEntrySurface || surface==kDraftSurface)
956 {
957 // Get the possible intersection point with the slope angle
958 const double z1 = CalcIntersection(pos, dir, fGrooves[ix].slope);
959
960 // We hit the slope angle
961 if (z1!=0)
962 {
963 // Move photon to new hit position
964 pos.PropagateZ(dir, z1);
965 if (fSlopeAbsorption)
966 throw raytrace_user(kEnter+kAbsorbed, surface, kSlopeSurface,
967 "EnterGroove - Photon absorbed by slope surface");
968
969 // Get the normal vector of the surface which was hit
970 const VectorNorm norm(fGrooves[ix].slope.theta_norm+RandomTheta(fPSF),
971 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
972
973 // Get the optical transition of the direction vector
974 const int ret = MOptics::ApplyTransition(dir, norm, 1, n0, fFresnelReflection);
975
976 // Transition was Reflection - try again
977 if (ret==1 || ret==2)
978 return kSlopeSurface;//EnterGroove(1, n0, lambda, pos, dir)+1;
979
980 // Transition was Refraction - enter
981 if (ret>=3)
982 return -kSlopeSurface;//LeavePeak(1, n0, lambda, pos, dir, pos.T())+1;
983
984 // Error occured (see ApplyTransition for details)
985 throw raytrace_error(kEnter+kTransitionError, surface, kSlopeSurface,
986 "EnterGroove - MOptics::ApplyTransition failed for slope surface");
987 }
988 }
989
990 // Photons was just injected (test both surfaces) or came from the other surface
991 if (surface==kEntrySurface || surface==kSlopeSurface)
992 {
993 const double z2 = CalcIntersection(pos, dir, fGrooves[ix].draft);
994
995 // We hit the draft angle
996 if (z2!=0)
997 {
998 // Move photon to new hit position
999 pos.PropagateZ(dir, z2);
1000 if (fDraftAbsorption)
1001 throw raytrace_user(kEnter+kAbsorbed, surface, kDraftSurface,
1002 "EnterGroove - Photon absorbed by draft surface");
1003
1004 // Get the normal vector of the surface which was hit
1005 const VectorNorm norm(fGrooves[ix].draft.theta_norm+RandomTheta(fPSF),
1006 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
1007
1008 // Get the optical transition of the direction vector
1009 const int ret = MOptics::ApplyTransition(dir, norm, 1, n0, fFresnelReflection);
1010
1011 // Transition was Reflection - try again
1012 if (ret==1 || ret==2)
1013 return kDraftSurface;//EnterGroove(2, n0, lambda, pos, dir)+1;
1014
1015 // Transition was Refraction - enter
1016 if (ret>=3)
1017 return -kDraftSurface;//LeavePeak(2, n0, lambda, pos, dir, pos.T())+1;
1018
1019 // Error occured (see ApplyTransition for details)
1020 throw raytrace_error(kEnter+kTransitionError, surface, kDraftSurface,
1021 "EnterGroove - MOptics::ApplyTransition failed for draft surface");
1022 }
1023 }
1024
1025 if (dir.Z()>0)
1026 {
1027 // We have missed both surfaces and we are upgoing...
1028 // ... ray can be discarded
1029 throw raytrace_info(kEnter+kStrayUpgoing, surface, kNoSurface,
1030 "EnterGroove - Particle is upgoing and has hit no surface");
1031 }
1032
1033 // The ray has left the peak at the bottom(?)
1034 throw raytrace_error(kEnter+kStrayDowngoing, surface, kNoSurface,
1035 "EnterGroove - Particle is downgoing and has hit no surface");
1036}
1037
1038/*
1039// Leave the peak from inside the material, either thought the draft surface or the
1040// slope surface or the bottom connecting the valley of both
1041int MFresnelLens::LeavePeak(int surface, double n0, double lambda, MQuaternion &pos, MQuaternion &dir, double T0) const
1042{
1043 const double rx = pos.R();
1044
1045 if (rx>=fR)
1046 {
1047 //cout << "Left the lens radius (leave)" << endl;
1048 throw -10;
1049 }
1050
1051 if (dir.Z()>0 && surface!=3) // && surface!=4)
1052 {
1053 //cout << "Upgoing, inside of the material" << endl;
1054 //pos.PropagateZ(dir, dir.Z()>0 ? 3 : -3);
1055 throw -11;
1056 }
1057
1058 if (surface!=1 && surface!=2 && surface!=3) // && surface!=4)
1059 {
1060 //cout << "Surface of origin invalid" << endl;
1061 throw -12;
1062 }
1063
1064
1065 // Calculate the ordinal number of the groove correpsonding to rx
1066 const int ix = TMath::FloorNint(rx/fW);
1067
1068 // FIXME: The Z-coordinate (cone.h) is actually a line through two points!!!
1069
1070 Cone slope = fGrooves[ix].slope;
1071 Cone draft = fGrooves[ix].draft;
1072
1073 const bool is_draft = rx>fGrooves[ix].r;
1074 if (is_draft)
1075 {
1076 // We are in the volume under the draft angle... taking the slope from ix+1
1077 if (ix<fGrooves.size()-1) // FIXME: Does that make sense?
1078 slope = fGrooves[ix+1].slope;
1079 }
1080 else
1081 {
1082 // We are in the volume under the slope angle... taking the draft from ix-1
1083 if (ix>0) // FIXME: Check whether this is correct
1084 draft = fGrooves[ix-1].draft;
1085 }
1086
1087 if (is_draft+1!=surface && (surface==1 || surface==2))
1088 cout << "SURFACE: " << is_draft+1 << " " << surface << endl;
1089
1090 if (surface==3)
1091 {
1092 //cout << "Upgoing, coming from the bottom of the lens" << endl;
1093 // Find out which triangle (peak) the photon is going to enter
1094 // then proceed...
1095 throw -13;
1096 }
1097
1098
1099 // We are inside the material and downgoing, so if we come from a slope surface,
1100 // we can only hit a draft surface after and vice versa
1101 if (is_draft || surface==3)
1102 {
1103 const double z1 = CalcIntersection(pos, dir, slope);
1104
1105 // We hit the slope angle and are currently in the volume under the draft surface
1106 if (z1!=0)
1107 {
1108 // Move photon to new hit position
1109 pos.PropagateZ(dir, z1);
1110
1111 if (fSlopeAbsorption)
1112 throw -200;
1113
1114 // Get the normal vector of the surface which was hit
1115 const VectorNorm norm(slope.theta_norm+RandomTheta(fPSF),
1116 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
1117
1118 // Get the optical transition of the direction vector
1119 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1);
1120
1121 // Transition was Reflection - try again
1122 if (ret==1 || ret==2)
1123 return LeavePeak(1, n0, lambda, pos, dir, T0)+1;
1124
1125 // Transition was Refraction - leave
1126 if (ret>=3)
1127 {
1128 if (!Transmission(pos.T()-T0, lambda))
1129 throw -14;
1130
1131 return EnterGroove(1, n0, lambda, pos, dir)+1;
1132 }
1133
1134 // Error occured (see ApplyTransition for details)
1135 //cout << "ERR[TIR3]" << endl;
1136 throw -15;
1137 }
1138 }
1139
1140 if (!is_draft || surface==3)
1141 {
1142 const double z2 = CalcIntersection(pos, dir, draft);
1143
1144 // We hit the draft angle from the inside and are currently in the volume under the slope angle
1145 if (z2!=0)
1146 {
1147 // Move photon to new hit position
1148 pos.PropagateZ(dir, z2);
1149
1150 if (fDraftAbsorption)
1151 throw -201;
1152
1153 // Get the normal vector of the surface which was hit
1154 const VectorNorm norm(draft.theta_norm+RandomTheta(fPSF),
1155 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
1156
1157 // Get the optical transition of the direction vector
1158 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1);
1159
1160 // Transition was Reflection - try again
1161 if (ret==1 || ret==2)
1162 return LeavePeak(2, n0, lambda, pos, dir, T0)+1;
1163
1164 // Transition was Refraction - leave
1165 if (ret>=3)
1166 {
1167 if (!Transmission(pos.T()-T0, lambda))
1168 throw -16;
1169
1170 return EnterGroove(2, n0, lambda, pos, dir)+1;
1171 }
1172
1173 // Error occured (see ApplyTransition for details)
1174 //cout << "ERR[TIR4]" << endl;
1175 throw -17;
1176 }
1177 }
1178
1179 if (surface==3)// || surface==4)
1180 {
1181 //cout << ix << " Lost bottom reflected ray " << surface << endl;
1182 throw -18;
1183 }
1184
1185 // The ray has left the peak at the bottom
1186
1187 // FIXME: There is a tiny chance to escape to the side
1188 // As there is a slope in the bottom surface of the peak
1189
1190 // Move photon to new hit position
1191 pos.PropagateZ(dir, -fH);
1192
1193 if (pos.R()>fR)
1194 {
1195 //cout << "Left the lens radius (bottom)" << endl;
1196 throw -19;
1197 }
1198
1199 // Get the normal vector of the surface which was hit
1200 const VectorNorm norm(RandomTheta(fPSF), gRandom->Uniform(0, TMath::TwoPi()));
1201
1202 // Get the optical transition of the direction vector
1203 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1);
1204
1205 // Transition was Reflection
1206 // (Photon scattered back from the bottom of the lens)
1207 if (ret==1 || ret==2)
1208 return LeavePeak(3, n0, lambda, pos, dir, T0)+1;
1209
1210 // Transition was Refraction
1211 // (Photon left at the bottom of the lens)
1212 if (ret>=3)
1213 {
1214 if (!Transmission(pos.T()-T0, lambda))
1215 throw -20;
1216
1217 return 0;
1218 }
1219
1220 // Error occured (see ApplyTransition for details)
1221 //cout << "ERR[TIR5]" << endl;
1222 throw -21;
1223}*/
1224
1225// Leave the peak from inside the material, either thought the draft surface or the
1226// slope surface or the bottom connecting the valley of both
1227int MFresnelLens::LeavePeak(int surface, double n0, MQuaternion &pos, MQuaternion &dir, double T0) const
1228{
1229 const double rx = pos.R();
1230
1231 if (rx>=fR) // This is an error as the direction vector is now invalid
1232 throw raytrace_error(kLeave+kOutsideRadius, surface, kNoSurface,
1233 "LeavePeak - Surface hit outside allowed radius");
1234
1235 // FIXME: Can we track them further?
1236 if (fDisableMultiEntry && dir.Z()>0 && surface!=3/* && surface!=4*/)
1237 throw raytrace_info(kLeave+kStrayUpgoing, surface, kNoSurface,
1238 "LeavePeak - Particle is upgoing inside the material and does not come from the bottom");
1239
1240 if (surface!=kSlopeSurface && surface!=kDraftSurface && surface!=kExitSurface/* && surface!=4*/)
1241 throw raytrace_error(kLeave+kInvalidOrigin, surface, kNoSurface,
1242 "LeavePeak - Invalid surface of origin");
1243
1244
1245 // Calculate the ordinal number of the groove correpsonding to rx
1246 const uint32_t ix = TMath::FloorNint(rx/fW);
1247
1248 // FIXME: The Z-coordinate (cone.h) is actually a line through two points!!!
1249
1250 Cone slope = fGrooves[ix].slope;
1251 Cone draft = fGrooves[ix].draft;
1252
1253 //if (is_draft+1!=surface && (surface==1 || surface==2))
1254 // cout << "SURFACE: " << is_draft+1 << " " << surface << endl;
1255
1256 const bool is_draft = rx>fGrooves[ix].r;
1257 if (is_draft)
1258 {
1259 // We are in the volume under the draft angle... taking the slope from ix+1
1260 if (ix<fGrooves.size()-1) // FIXME: Does that make sense?
1261 slope = fGrooves[ix+1].slope;
1262 }
1263 else
1264 {
1265 // We are in the volume under the slope angle... taking the draft from ix-1
1266 if (ix>0) // FIXME: Check whether this is correct
1267 draft = fGrooves[ix-1].draft;
1268 }
1269
1270 if (surface==kExitSurface)
1271 {
1272 if (!fBottomReflection)
1273 throw raytrace_user(kLeave+kAbsorbed, surface, kExitSurface,
1274 "LeavePeak - Particle absorbed on the bottom");
1275
1276 const int in = FindPeak(ix, pos, dir);
1277
1278 // This might happen if the ray is very flat and leaving
1279 // the lens before hitting the border boundary of the grooves
1280 if (in<0)
1281 throw raytrace_error(kLeave+kNoSurfaceFound, kExitSurface, kNoSurface,
1282 "LeavePeak - No hit surface found for particle reflected at the bottom");
1283
1284 slope = fGrooves[in].slope;
1285 draft = fGrooves[in==0 ? 0 : in-1].draft;
1286 }
1287
1288 // FIXME: There is a chance that we can hit the same surface twice (for very horizontal rays
1289 // but this requires a proper selection of the hit point
1290
1291 // We are inside the material and downgoing, so if we come from a slope surface,
1292 // we can only hit a draft surface after and vice versa
1293 if (is_draft || surface==kExitSurface)
1294 {
1295 const double z1 = CalcIntersection(pos, dir, slope);
1296
1297 // We hit the slope angle and are currently in the volume under the draft surface
1298 if (z1!=0)
1299 {
1300 // Move photon to new hit position
1301 pos.PropagateZ(dir, z1);
1302
1303 if (fSlopeAbsorption)
1304 throw raytrace_user(kLeave+kAbsorbed, surface, kSlopeSurface,
1305 "LeavePeak - Photon absorbed by slope surface");
1306
1307 // Get the normal vector of the surface which was hit
1308 const VectorNorm norm(slope.theta_norm+RandomTheta(fPSF),
1309 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
1310
1311 // Get the optical transition of the direction vector
1312 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1, fFresnelReflection);
1313
1314 // Transition was Reflection - try again
1315 if (ret==1 || ret==2)
1316 return -kSlopeSurface;//LeavePeak(1, n0, lambda, pos, dir, T0)+1;
1317
1318 // Transition was Refraction - leave
1319 if (ret>=3) // Transmission
1320 return kSlopeSurface;//EnterGroove(1, n0, lambda, pos, dir)+1;
1321
1322 // Error occured (see ApplyTransition for details)
1323 throw raytrace_error(kLeave+kTransitionError, surface, kSlopeSurface,
1324 "LeavePeak - MOptics::ApplyTransition failed for slope surface");
1325 }
1326 }
1327
1328 if (!is_draft || surface==kExitSurface)
1329 {
1330 const double z2 = CalcIntersection(pos, dir, draft);
1331
1332 // We hit the draft angle from the inside and are currently in the volume under the slope angle
1333 if (z2!=0)
1334 {
1335 // Move photon to new hit position
1336 pos.PropagateZ(dir, z2);
1337
1338 if (fDraftAbsorption)
1339 throw raytrace_user(kLeave+kAbsorbed, surface, kDraftSurface,
1340 "LeavePeak - Photon absorbed by draft surface");
1341
1342 // Get the normal vector of the surface which was hit
1343 const VectorNorm norm(draft.theta_norm+RandomTheta(fPSF),
1344 pos.XYvector().Phi()+RandomPhi(pos.R(), fPSF));
1345
1346 // Get the optical transition of the direction vector
1347 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1, fFresnelReflection);
1348
1349 // Transition was Reflection - try again
1350 if (ret==1 || ret==2)
1351 return -kDraftSurface;//LeavePeak(2, n0, lambda, pos, dir, T0)+1;
1352
1353 // Transition was Refraction - leave
1354 if (ret>=3) // Transmission
1355 return kDraftSurface;//EnterGroove(2, n0, lambda, pos, dir)+1;
1356
1357 // Error occured (see ApplyTransition for details)
1358 //cout << "ERR[TIR4]" << endl;
1359 throw raytrace_error(kLeave+kTransitionError, surface, kDraftSurface,
1360 "LeavePeak - MOptics::ApplyTransition failed for draft surface");
1361 }
1362 }
1363
1364 if (surface==kExitSurface/* || surface==4*/)
1365 throw raytrace_error(kLeave+kFoundSurfaceUnavailable, kExitSurface, is_draft?kSlopeSurface:kDraftSurface,
1366 "LeavePeak - Ray reflected on the bottom did not hit the found surface");
1367
1368 // The ray has left the peak at the bottom
1369
1370 // FIXME: There is a tiny chance to escape to the side
1371 // As there is a slope in the bottom surface of the peak
1372
1373 // FIXME: Theoretically, a ray can hit the same surface twice
1374
1375 // Move photon to new hit position
1376 pos.PropagateZ(dir, -fH);
1377
1378 if (pos.R()>fR)
1379 throw raytrace_info(kLeave+kOutsideRadius, surface, kExitSurface,
1380 "LeavePeak - Hit point at the bottom surface is beyond allowed radius");
1381
1382 // Get the normal vector of the surface which was hit
1383 const VectorNorm norm(RandomTheta(fPSF), gRandom->Uniform(0, TMath::TwoPi()));
1384
1385 // Get the optical transition of the direction vector
1386 const int ret = MOptics::ApplyTransition(dir, norm, n0, 1, fFresnelReflection);
1387
1388 // Transition was Reflection
1389 // (Photon scattered back from the bottom of the lens)
1390 if (ret==1 || ret==2)
1391 return -kExitSurface;//LeavePeak(3, n0, lambda, pos, dir, T0)+1;
1392
1393 // Transition was Refraction
1394 // (Photon left at the bottom of the lens)
1395 if (ret>=3) // Transmission
1396 return kPhotonHasLeft;
1397
1398 // Error occured (see ApplyTransition for details)
1399 throw raytrace_error(kLeave+kTransitionError, surface, kExitSurface, "LeavePeak - MOptics::ApplyTransition failed for bottom surface");
1400}
1401
1402
1403// Differences:
1404// Returns a 'reflected' vector at z=0
1405// Does not propagate to z=0 at the beginning
1406Int_t MFresnelLens::ExecuteOptics(MQuaternion &p, MQuaternion &u, const Short_t &wavelength) const
1407{
1408 // Corsika Coordinates are in cm!
1409
1410 const double lambda = wavelength==0 ? fLambda : wavelength;
1411 if (fAbsorptionLength.GetNp()!=0 &&
1412 (lambda<fAbsorptionLength.GetXmin() || lambda>fAbsorptionLength.GetXmax()))
1413 {
1414 *fLog << err << "Wavelength " << lambda << "nm out of absorption range [" << fAbsorptionLength.GetXmin() << "nm;" << fAbsorptionLength.GetXmax() << "nm]" << endl;
1415 return -1;
1416 }
1417
1418 const double n0 = MFresnelLens::RefractiveIndex(lambda);
1419
1420 try
1421 {
1422 int last_surface = kEntrySurface;//EnterGroove(kEntrySurface, n0, p, u);
1423
1424 // last_surface that was hit (photon originates from)
1425 // 0 entrance (Z=0) or exit (Z=-fH) surface
1426 // 1 slope
1427 // 2 draft
1428 // 3 bottom
1429 // positive: photon is outside of material --> Try to enter
1430 // nagative: photon is inside of material --> Try to leave
1431
1432 double T0 = 0;//last_surface<0 ? p.T() : 0;
1433
1434 // The general assumption is: no surface can be hit twice in a row
1435
1436 int cnt = -1;
1437 while (last_surface!=0)
1438 {
1439 cnt ++;
1440
1441 // photon is outside of material --> try to enter
1442 if (last_surface>0)
1443 {
1444 last_surface = EnterGroove(last_surface, n0, p, u);
1445
1446 // successfully entered --> remember time of entrance to calculate transimission
1447 if (last_surface<0)
1448 T0 = p.T();
1449
1450 continue;
1451 }
1452
1453 // photon is inside of material --> try to leave
1454 if (last_surface<0)
1455 {
1456 last_surface = LeavePeak(-last_surface, n0, p, u, T0);
1457
1458 // successfully left --> apply transmission
1459 if (last_surface>=0)
1460 {
1461 if (!Transmission(p.T()-T0, lambda))
1462 throw raytrace_error(kAbsorbed, last_surface, kMaterial,
1463 "TraceRay - Ray absorbed in material");
1464 }
1465
1466 continue;
1467 }
1468 }
1469
1470 // To make this consistent with a mirror system,
1471 // we now change our coordinate system
1472 // Rays from the lens to the camera are up-going (positive sign)
1473 u.fVectorPart.SetZ(-u.Z());
1474
1475 // In the datasheet, it looks as if F is calculated
1476 // towards the center of the lens. It seems things are more
1477 // consistent if the thickness correction in caluating the
1478 // slope angle is omitted and the focal distance is measured
1479 // from the entrance of the lens => FIXME: To be checked
1480 // (Propagating to F means not propagating a distance of F-H from the exit)
1481 //p.fVectorPart.SetZ(fH-fH/2/fN);//fH/2); Found by try-and-error
1482
1483 // We are already at -H, adding F and setting Z=0 means going to -(F+H)
1484 p.fVectorPart.SetZ(0);//fH/2); Found by try-and-error
1485
1486 return uint32_t(cnt)>=fMinHits && (fMaxHits==0 || uint32_t(cnt)<=fMaxHits) ? cnt : -1;;
1487 }
1488 catch (const raytrace_exception &e)
1489 {
1490 return -e.id();
1491 }
1492
1493/*
1494 try
1495 {
1496 const int cnt = EnterGroove(0, n0, lambda, p, u);
1497
1498 // To make this consistent with a mirror system,
1499 // we now change our coordinate system
1500 // Rays from the lens to the camera are up-going (positive sign)
1501 u.fVectorPart.SetZ(-u.Z());
1502
1503 // In the datasheet, it looks as if F is calculated
1504 // towards the center of the lens
1505 // (Propagating to F means not propagating a distance of F-H/2)
1506 p.fVectorPart.SetZ(0);
1507
1508 return cnt>=fMinHits && (fMaxHits==0 || cnt<=fMaxHits) ? cnt : -1;
1509
1510 }
1511 catch (const int &rc)
1512 {
1513 return rc;
1514 }
1515*/
1516}
1517
1518// Differences:
1519// Does propagate to z=0 at the beginning
1520Int_t MFresnelLens::TraceRay(vector<MQuaternion> &vec, MQuaternion &p, MQuaternion &u, const Short_t &wavelength, bool verbose) const
1521{
1522 // Corsika Coordinates are in cm!
1523
1524 const double lambda = wavelength==0 ? fLambda : wavelength;
1525 if (fAbsorptionLength.GetNp()!=0 &&
1526 (lambda<fAbsorptionLength.GetXmin() || lambda>fAbsorptionLength.GetXmax()))
1527 {
1528 *fLog << err << "Wavelength " << lambda << "nm out of absorption range [" << fAbsorptionLength.GetXmin() << "nm;" << fAbsorptionLength.GetXmax() << "nm]" << endl;
1529 return -1;
1530 }
1531
1532 const double n0 = MFresnelLens::RefractiveIndex(lambda);
1533
1534 // Photon must be at the lens surface
1535 p.PropagateZ(u, 0);
1536 vec.push_back(p);
1537
1538 try
1539 {
1540 int last_surface = kEntrySurface;//EnterGroove(kEntrySurface, n0, p, u);
1541
1542 // last_surface that was hit (photon originates from)
1543 // 0 entrance (Z=0) or exit (Z=-fH) surface
1544 // 1 slope
1545 // 2 draft
1546 // 3 bottom
1547 // positive: photon is outside of material --> Try to enter
1548 // nagative: photon is inside of material --> Try to leave
1549
1550 double T0 = 0;
1551
1552 // The general assumption is: no surface can be hit twice in a row
1553
1554 int cnt = -1;
1555 while (last_surface!=0)
1556 {
1557 cnt ++;
1558 vec.push_back(p);
1559
1560 // photon is outside of material --> try to enter
1561 if (last_surface>0)
1562 {
1563 last_surface = EnterGroove( last_surface, n0, p, u);
1564 //cout << "enter = " << last_surface << endl;
1565
1566 // successfully entered --> remember time of entrance to calculate transimission
1567 if (last_surface<0)
1568 T0 = p.T();
1569
1570 continue;
1571 }
1572
1573 // photon is inside of material --> try to leave
1574 if (last_surface<0)
1575 {
1576 last_surface = LeavePeak(-last_surface, n0, p, u, T0);
1577 //cout << "leave = " << last_surface << endl;
1578
1579 // successfully left --> apply transmission
1580 if (last_surface>=0)
1581 {
1582 if (!Transmission(p.T()-T0, lambda))
1583 throw raytrace_error(kAbsorbed, last_surface, kMaterial,
1584 "TraceRay - Ray absorbed in material");
1585 }
1586
1587 continue;
1588 }
1589 }
1590
1591 vec.push_back(p);
1592 return cnt;
1593 }
1594 catch (const raytrace_exception &e)
1595 {
1596 if (verbose)
1597 *fLog << all << e.id() << ": " << e.what() << endl;
1598
1599 // Hit point at bottom surface beyond allowed range
1600 // FIXME: Only if surface is kExitSurface
1601 if (e.id()==2342)
1602 vec.push_back(p);
1603
1604 return -e.id();
1605 }
1606}
Note: See TracBrowser for help on using the repository browser.