Subversion Repositories f9daq

Rev

Rev 25 | Rev 70 | Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
25 f9daq 1
#include "include/guide.h"
2
 
3
#include <iostream>
4
 
5
// vector output shortcut
6
void printv(TVector3 v)
7
{
8
        printf("(x,y,z) = (%.4lf, %.4lf, %.4lf)\n", v.x(), v.y(), v.z());
9
}
10
// TVector3::Rotate does not seem accurate enough
11
TVector3 rotatey(TVector3 v, double theta)
12
{
13
        return TVector3(v.x() * TMath::Cos(theta) + v.z() * TMath::Sin(theta),
14
                                        v.y(),
15
                                        -v.x() * TMath::Sin(theta) + v.z() * TMath::Cos(theta));
16
}
17
// another shortcut not found in TMath
18
int sign(double in)
19
{
20
        if(in >= 0.0) return 1;
21
        else return -1;
22
}
23
//=================================================================================
24
 
25
//-----------------------------------------------------------------------------
26
void CRay::Set(TVector3 r0, TVector3 n0)
27
{
28
        r = r0; n = n0.Unit();
29
}
30
//-----------------------------------------------------------------------------
31
//void CRay::Set(double x0, double y0, double z0, double l0, double m0, double n0)
32
//{
33
        //r.SetXYZ(x0, y0, z0);
34
        //n.SetXYZ(l0, m0, n0); n = n.Unit();
35
//}
36
//-----------------------------------------------------------------------------
37
/*
38
CRay& CRay::operator = (const CRay& p)
39
{
40
        r.SetXYZ(p.GetR().x(), p.GetR().y(), p.GetR().z());
41
        //this->r.SetXYZ(p.x(), p.y(), p.z());
42
        n.SetXYZ(p.GetN().x(), p.GetN().y(), p.GetN().z());
43
        return *this;
44
} */
45
//-----------------------------------------------------------------------------
46
void CRay::Print()
47
{
48
        printf("---> CRay::Print() <---\n");
49
        printf("(x,y,z)=(%.2lf, %.2lf, %.2lf); (l,m,n)=(%.2lf, %.2lf, %.2lf)\n",
50
                r.x(), r.y(), r.z(), n.x(), n.y(), n.z());
51
}
52
//-----------------------------------------------------------------------------
53
void CRay::Draw()
54
{
55
double t = 50.0;
56
TPolyLine3D *line3d = new TPolyLine3D(2);
57
        //line3d->SetPoint(0, r.x() - t*n.x(), r.y() - t*n.y(), r.z() - t*n.z());
58
        line3d->SetPoint(0, r.x(), r.y(), r.z());
59
        line3d->SetPoint(1, r.x() + t*n.x(), r.y() + t*n.y(), r.z() + t*n.z());
60
        line3d->SetLineWidth(1);
61
        line3d->SetLineColor(color);
62
 
63
        line3d->Draw();
64
}
65
//-----------------------------------------------------------------------------
66
void CRay::Draw(double x_from, double x_to)
67
{
68
double A1, A2;
69
  TPolyLine3D *line3d = new TPolyLine3D(2);
70
 
71
        if(n.x() < MARGIN) {
72
                A1 = A2 = 0.0;
73
        } else {
74
                A1 = (x_from - r.x())/n.x();   
75
                A2 = (x_to - r.x())/n.x();
76
        }
77
 
78
        line3d->SetPoint(0, x_from, A1*n.y()+r.y(), A1*n.z()+r.z());
79
        line3d->SetPoint(1, x_to, A2*n.y()+r.y(), A2*n.z()+r.z());
80
        line3d->SetLineWidth(1);
81
        line3d->SetLineColor(color);
82
 
83
        line3d->Draw();
84
}
85
//-----------------------------------------------------------------------------
86
void CRay::DrawS(double x_from, double t)
87
{
88
        double A1;
89
        TPolyLine3D *line3d = new TPolyLine3D(2);
90
 
91
        if(n.x() < MARGIN)
92
                A1 = 0.0;
93
        else
94
                A1 = (x_from - r.x())/n.x();
95
 
96
        line3d->SetPoint(0, x_from, A1*n.y()+r.y(), A1*n.z()+r.z());
97
        line3d->SetPoint(1, r.x() + t*n.x(), r.y() + t*n.y(), r.z() + t*n.z());
98
        line3d->SetLineWidth(1);
99
        line3d->SetLineColor(color);
100
 
101
        line3d->Draw();
102
}
103
//=================================================================================
104
 
105
 
106
//=================================================================================
107
CPlane4::CPlane4() :
108
        n(TVector3(1.0, 0.0, 0.0)),
109
  A(0),
110
  B(0),
111
  C(0),
112
  D(0)
113
{ r[0] = TVector3(0.0,-1.0,-1.0);
114
        r[1] = TVector3(0.0,-1.0, 1.0);
115
        r[2] = TVector3(0.0, 1.0, 1.0);
116
        r[3] = TVector3(0.0, 1.0,-1.0);
117
        for(int i=0;i<4;i++) edge[i] = TVector3(0,0,0);
118
        for(int i=0;i<4;i++) angle_r[i] = 0;
119
};
120
//-----------------------------------------------------------------------------
121
CPlane4::CPlane4(TVector3 r1, TVector3 r2, TVector3 r3, TVector3 r4)
122
{
123
        //Set(r1, r2, r3, r4);
124
//}
125
//-----------------------------------------------------------------------------
126
// za izracun parametrov ravnine je en vektor prevec, vendar tega 
127
// rabim kot zadnji vogal poligona
128
//void CPlane4::Set(TVector3 r1, TVector3 r2, TVector3 r3, TVector3 r4)
129
//{
130
double x1,y1,z1, x2,y2,z2, x3,y3,z3;
131
 
132
        x1 = r1.x(); y1 = r1.y(); z1 = r1.z();
133
        x2 = r2.x(); y2 = r2.y(); z2 = r2.z();
134
        x3 = r3.x(); y3 = r3.y(); z3 = r3.z();
135
 
136
        A = y3*(z1 - z2) + y1*(z2 - z3) + y2*(z3 - z1);
137
        B = x3*(z2 - z1) + x1*(z3 - z2) + x2*(z1 - z3);
138
        C = x3*(y1 - y2) + x1*(y2 - y3) + x2*(y3 - y1);
139
        D = y3*(x1*z2 - x2*z1) + x3*(y2*z1 - y1*z2) + z3*(x2*y1 - x1*y2);
140
 
141
        r[0] = r1; r[1] = r2; r[2] = r3; r[3] = r4;
142
        n.SetXYZ(A, B, C);
143
        n = n.Unit();
144
 
145
        for(int i=0;i<4;i++)
146
                edge[i] = r[i-3 ? i+1 : 0] - r[i];
147
 
148
        for(int i=0;i<4;i++)
149
                angle_r[i] = TMath::ACos(/*TMath::Abs*/( ((-edge[i ? i-1 : 3]).Unit()) * (edge[i].Unit()) ));
150
};
151
 
152
void CPlane4::Set(TVector3 r1, TVector3 r2, TVector3 r3, TVector3 r4)
153
{
154
  double x1,y1,z1, x2,y2,z2, x3,y3,z3;
155
 
156
        x1 = r1.x(); y1 = r1.y(); z1 = r1.z();
157
        x2 = r2.x(); y2 = r2.y(); z2 = r2.z();
158
        x3 = r3.x(); y3 = r3.y(); z3 = r3.z();
159
 
160
        A = y3*(z1 - z2) + y1*(z2 - z3) + y2*(z3 - z1);
161
        B = x3*(z2 - z1) + x1*(z3 - z2) + x2*(z1 - z3);
162
        C = x3*(y1 - y2) + x1*(y2 - y3) + x2*(y3 - y1);
163
        D = y3*(x1*z2 - x2*z1) + x3*(y2*z1 - y1*z2) + z3*(x2*y1 - x1*y2);
164
 
165
        r[0] = r1; r[1] = r2; r[2] = r3; r[3] = r4;
166
        n.SetXYZ(A, B, C);
167
        n = n.Unit();
168
 
169
        for(int i=0;i<4;i++)
170
                edge[i] = r[i-3 ? i+1 : 0] - r[i];
171
 
172
        for(int i=0;i<4;i++)
173
                angle_r[i] = TMath::ACos(/*TMath::Abs*/( ((-edge[i ? i-1 : 3]).Unit()) * (edge[i].Unit()) ));
174
};
175
 
176
CPlane4::CPlane4(TVector3 *vr)
177
{
178
  double x1,y1,z1, x2,y2,z2, x3,y3,z3;
179
 
180
        x1 = vr[0].x(); y1 = vr[0].y(); z1 = vr[0].z();
181
        x2 = vr[1].x(); y2 = vr[1].y(); z2 = vr[1].z();
182
        x3 = vr[2].x(); y3 = vr[2].y(); z3 = vr[2].z();
183
 
184
        A = y3*(z1 - z2) + y1*(z2 - z3) + y2*(z3 - z1);
185
        B = x3*(z2 - z1) + x1*(z3 - z2) + x2*(z1 - z3);
186
        C = x3*(y1 - y2) + x1*(y2 - y3) + x2*(y3 - y1);
187
        D = y3*(x1*z2 - x2*z1) + x3*(y2*z1 - y1*z2) + z3*(x2*y1 - x1*y2);
188
 
189
        r[0] = vr[0]; r[1] = vr[1]; r[2] = vr[2]; r[3] = vr[3];
190
        n.SetXYZ(A, B, C);
191
        n = n.Unit();
192
 
193
        for(int i=0;i<4;i++)
194
                edge[i] = r[i-3 ? i+1 : 0] - r[i];
195
 
196
        for(int i=0;i<4;i++)
197
                angle_r[i] = TMath::ACos(/*TMath::Abs*/( ((-edge[i ? i-1 : 3]).Unit()) * (edge[i].Unit()) ));
198
};
199
//-----------------------------------------------------------------------------
200
// posce presecisce !neskoncne! ravnine s premico (class CRay)
201
// ce najde presecisce vrne 1
202
int CPlane4::GetIntersection(TVector3 *vec, CRay ray)
203
{
204
TVector3 N; //nenormirani vektor (A,B,C)
205
double num, den; //stevec, imenovalec
206
double t;
207
TVector3 tmp;
208
 
209
        N.SetXYZ(A,B,C);
210
 
211
        num = N*ray.GetR() + D;
212
        den = N*ray.GetN();
213
 
214
        if (dbg) printf("t = %6.3lf / %6.3lf =  %6.3lf\n", num, den, num/den);
215
 
216
        //if(den == 0)
217
        if(TMath::Abs(den) < MARGIN) {
218
                //if(num == 0)
219
                if(TMath::Abs(num) < MARGIN) {
220
                  if (dbg) printf("The ray is on the surface!\n");
221
                        return 0; //return 2; // premica lezi na ravnini
222
                }
223
                else {
224
                        if (dbg) printf("The ray is parallel to the surface!\n");
225
                        return 0; // ni presecisca
226
                }
227
        }
228
 
229
        t = num / den;
230
 
231
        tmp = ray.GetR();
232
        tmp -= t*ray.GetN();
233
        *vec = tmp;
234
        return 1;
235
}
236
//-----------------------------------------------------------------------------
237
// ali je vektor vec, ki lezi na ravnini skupaj z e1 in e2, med njima
238
// angle_r je kot med e1 in e2, vsi vektorji imajo skupno izhodisce
239
int CPlane4::IsInTri(TVector3 vec, TVector3 e1, TVector3 e2, double angle)
240
{
241
  double angle_ve1, angle_ve2;
242
 
243
        if(dbg) printf("--- CPlane4::IsInTri ---\n");
244
 
245
        angle_ve1 = TMath::ACos(/*TMath::Abs*/( (e1.Unit()) * (vec.Unit()) ));
246
        angle_ve2 = TMath::ACos(/*TMath::Abs*/( (e2.Unit()) * (vec.Unit()) ));
247
 
248
        if(dbg)
249
        {
250
                //printf("angle_ve1 = %lf\n", angle_ve1*DEGREE);
251
                //printf("angle_ve2 = %lf\n", angle_ve2*DEGREE);
252
                //printf("angle_sum = %lf\n", (angle_ve1 + angle_ve2)*DEGREE);
253
                printf("  angle_r   = %lf\n", angle*DEGREE);
254
        }
255
 
256
  bool difference = (MARGIN < TMath::Abs(angle - (angle_ve1 + angle_ve2)));
257
  if (dbg) printf("  MARGIN < Difference = %d\n", difference);
258
  return (int) !difference;
259
}
260
//-----------------------------------------------------------------------------
261
// ali je vektor vec, ki lezi na ravnini!, znotraj meja, ki jih definirajo
262
// strije vogali te ravnine r[i]
263
int CPlane4::IsVectorIn(TVector3 vec)
264
{
265
  int status;
266
 
267
        if(dbg) printf("--- CPlane4::IsVectorIn ---\n");
268
 
269
        for(int i=0;i<3;i++)
270
        {
271
                status = IsInTri(vec - r[i], edge[i], -edge[i ? i-1 : 3], angle_r[i]);
272
                if(dbg) printf("  [%d] vec is %s\n", i, status ? "inside" : "outside");
273
                if(!status) return 0;
274
        }
275
 
276
        return 1;
277
}
278
//-----------------------------------------------------------------------------
279
int CPlane4::TestIntersection(CRay in)
280
{
281
        TVector3 tmp;
282
 
283
        if( GetIntersection(&tmp, in) )
284
                if( IsVectorIn(tmp) )
285
                        return 1;
286
 
287
        return 0;
288
}
289
//-----------------------------------------------------------------------------
290
int CPlane4::TestIntersection(TVector3 *vec, CRay in)
291
{
292
        TVector3 tmp;
293
 
294
        if( GetIntersection(&tmp, in) )
295
                if( IsVectorIn(tmp) ) {
296
                        *vec = tmp;
297
                        return 1;
298
                }
299
 
300
        return 0;
301
}
302
//-----------------------------------------------------------------------------
303
void CPlane4::Print()
304
{
305
        printf("--- CPlane4::Print() ---\n");
306
        printf("  r=(%.2lf, %.2lf, %.2lf); n=(%.2lf, %.2lf, %.2lf); ",
307
                r[0].x(), r[0].y(), r[0].z(), n.x(), n.y(), n.z());
308
        printf(  "(A,B,C,D)=(%.2lf, %.2lf, %.2lf, %.2lf) \n", A, B, C, D);
309
        for(int i=0;i<4;i++) printf("  edge[%d] = (%lf, %lf, %lf)\n", i, edge[i].x(), edge[i].y(), edge[i].z());
310
        for(int i=0;i<4;i++) printf("  angle[%d] = %lf\n", i, angle_r[i]*DEGREE);
311
}
312
//-----------------------------------------------------------------------------
313
void CPlane4::Draw(int color, int width)
314
{
315
TPolyLine3D *line3d = new TPolyLine3D(5);
316
 
317
        for(int i=0;i<4;i++) line3d->SetPoint(i, r[i].x(), r[i].y(), r[i].z());
318
        line3d->SetPoint(4, r[0].x(), r[0].y(), r[0].z());
319
        line3d->SetLineWidth(width); line3d->SetLineColor(color);
320
 
321
        line3d->Draw();
322
}
323
//=================================================================================
324
 
325
 
326
//=================================================================================
327
CSurface::CSurface(int type0):
328
  type(type0)
329
{
330
TVector3 vr[4];
331
TDatime now;
332
 
333
        vr[0].SetXYZ(0.0,-1.0,-1.0);
334
        vr[1].SetXYZ(0.0,-1.0, 1.0);
335
        vr[2].SetXYZ(0.0, 1.0, 1.0);
336
        vr[3].SetXYZ(0.0, 1.0,-1.0);
337
        //CPlane4::Set(vr);
338
        SetIndex(1.0, 1.5);
339
 
340
        reflection = c_reflectivity;   
341
        rand.SetSeed(now.Get());
342
 
343
        SetFresnel();
344
}
345
//-----------------------------------------------------------------------------
346
CSurface::CSurface(int type0, TVector3 r1, TVector3 r2, TVector3 r3, TVector3 r4, double n10, double n20, double reflectivity)
347
{
348
TDatime now;
349
 
350
        type = type0; CPlane4::Set(r1, r2, r3, r4);
351
        SetIndex(n10, n20);
352
 
353
        reflection = reflectivity;
354
        rand.SetSeed(now.Get());
355
 
356
        SetFresnel();
357
}
358
//-----------------------------------------------------------------------------
359
CSurface::CSurface(int type0, TVector3 *vr, double n10, double n20, double reflectivity)
360
{
361
TDatime now;
362
 
363
        type = type0; CPlane4::Set(vr);
364
        SetIndex(n10, n20);
365
 
366
        reflection = reflectivity;
367
        rand.SetSeed(now.Get());
368
 
369
        SetFresnel();
370
}
371
//-----------------------------------------------------------------------------
372
void CSurface::SetIndex(double n10, double n20)
373
{
374
        n1 = n10; n2 = n20; n1_n2 = n1/n2;
375
 
376
        if(n1 > n2)
377
                cosTtotal = TMath::Sqrt( 1 - TMath::Power(n2/n1, 2) );
378
        else
379
                cosTtotal = 0.0;
380
}
381
//-----------------------------------------------------------------------------
382
// sprejme zarek, vrne uklonjen/odbit zarek in presecisce
383
// vrne 0 ce ni presecisca; 1 ce se je lomil
384
// 2 ce se je odbil; -2 ce se je absorbiral
385
int CSurface::PropagateRay(CRay in, CRay *out, TVector3 *intersection)
386
{
387
  if (dbg) printf("--- CSurface::PropagateRay ---\n");
388
  double cosTi, cosTt, p_ref;
389
  TVector3 intersect, transmit;
390
 
391
        if( !(GetIntersection(&intersect, in) == 1) )
392
                return 0;
393
 
394
        *intersection = intersect;
395
        if( !IsVectorIn(intersect) )
396
                return 0;
397
 
398
        // --------------- Fresnel ----------------------------------------------------
399
        // R_f = a_te * R_te  +  a_tm * R_tm
400
        double r_te=0;
401
        double r_tm=0;
402
        double R_te=0;
403
        double R_tm=0;
404
        double R_f = 0.0;
405
        double a_te = 0.0;
406
        double a_tm = 0.0;
407
        TVector3 v_te; // polarization perpendicular to the plane of incidence
408
        TVector3 v_tm; // inbound polarization parallel with the plane of incidence
409
        TVector3 v_tm_t;// transmited polarization parallel with the plane of incidence
410
        TVector3 pol_t = in.GetP(); // transmited polarization
411
        int sign_n; // sign of normal direction vs. inbound ray
412
        double cosTN; // debug
413
 
414
        if(fresnel) {  
415
                v_te = n.Cross(in.GetN()); v_te = v_te.Unit();
416
                v_tm = -v_te.Cross(in.GetN()); v_tm = v_tm.Unit();
417
                if(dbg) {
418
                        printf("  v_te = "); printv(v_te);
419
                        printf("  v_tm = "); printv(v_tm);
420
                }
421
 
422
                double cosAf = v_te * in.GetP();
423
                if(dbg) printf("  cosAf = %lf (Af = %lf)\n", cosAf, TMath::ACos(cosAf)*DEGREE);
424
 
425
                a_te = cosAf;
426
                a_tm = TMath::Sqrt(1 - cosAf*cosAf);
427
                if(dbg) printf("  a_te = %lf, a_tm = %lf\n", a_te, a_tm);
428
        }
429
        // ----------------------------------------------------------------------------
430
 
431
        if(type == 3) type = SURF_REFRA; //SURF_TOTAL -> SURF_REFRA
432
        switch(type){
433
                // ----------------------------------------------------------------------------
434
                // --------------- refraction from n1 to n2 -----------------------------------
435
                // ----------------------------------------------------------------------------
436
                case SURF_REFRA:
437
                        cosTi = in.GetN() * n;
438
                        if(dbg) printf("  cosTi = %lf (Ti = %lf)\n", cosTi, TMath::ACos(cosTi)*DEGREE);
439
                        sign_n = -sign(cosTi);
440
                        if(dbg) printf("  sign_n = %d\n", sign_n);
441
                        cosTi = TMath::Abs(cosTi);
442
 
443
                        // Check if there can be total reflection: n1 > n2
444
                        if(N1_N2(-sign_n) < 1.0)
445
                                cosTtotal = TMath::Sqrt( 1 - TMath::Power(N1_N2(-sign_n), 2) );
446
                        else
447
                                cosTtotal = 0.0;
448
 
449
                        if(dbg) printf("  cosTtotal = %lf (Ttotal = %lf)\n", cosTtotal, TMath::ACos(cosTtotal)*DEGREE);
54 f9daq 450
                        // reflection dependance on polarization missing
451
                        // reflection hardcoded to 0.96
452
                        p_ref = rand.Uniform(0.0, 1.0);
453
                        if (dbg) printf("   reflection probability = %f\n", p_ref);    
25 f9daq 454
 
54 f9daq 455
                        // Total reflection
456
                        /*
25 f9daq 457
                        if( (cosTi <= cosTtotal) && (p_ref < reflection) ) { // totalni odboj z verjetnostjo "reflection"
458
                                if(dbg) printf("  TOTAL\n");
459
                                transmit = in.GetN() + sign_n*2*cosTi*n;
460
 
461
                                if(dbg) {
462
                                        cosTN = TMath::Abs(transmit.Unit() * n);
463
                                        printf("  cosTN = %lf (TN = %lf) (Abs(TN) = %lf)\n", cosTN, TMath::ACos(cosTN)*DEGREE, TMath::ACos(TMath::Abs(cosTN))*DEGREE);
464
                                }      
465
                                out->Set(intersect, transmit);
466
 
467
                                pol_t = -in.GetP() + sign_n*2*cosTi*n;
468
                                out->SetPolarization(pol_t);
54 f9daq 469
                                return REFLECTION;
470
                        } else { */
471
                          // reflection or refraction according to Fresnel equations
25 f9daq 472
                                if(dbg) printf("  REFRACTION\n");
473
                                if(dbg) printf("  N1_N2(sign_n) = %lf\n", N1_N2(sign_n));      
474
                                cosTt = TMath::Sqrt(1 - TMath::Power(N1_N2(sign_n), 2)*(1 - TMath::Power(cosTi, 2)));                          
475
                                if(dbg) printf("  cosTt = %lf (Tt = %lf) \n", cosTt, TMath::ACos(cosTt)*DEGREE);
476
 
477
                                transmit = N1_N2(sign_n)*in.GetN() + sign_n*(N1_N2(sign_n)*cosTi - cosTt)*n;
478
                                if(dbg) {printf("  transmit.Unit() = "); printv(transmit.Unit());}
479
 
480
                                if(dbg) {
481
                                        cosTN = transmit.Unit() * n;
482
                                        printf("  cosTN = %lf (TN = %lf) (Abs(TN) = %lf)\n", cosTN, TMath::ACos(cosTN)*DEGREE, TMath::ACos(TMath::Abs(cosTN))*DEGREE); 
483
                                }
484
                                //if(cosTi<=cosTtotal) cosTt = TMath::Sqrt(1 - TMath::Power(N1_N2(sign_n), 2)*(1 - TMath::Power(cosTi, 2)));
54 f9daq 485
                                //if(fresnel) {
25 f9daq 486
                                        r_te = (n1*cosTi - n2*cosTt)/(n1*cosTi + n2*cosTt); // transverse
487
                                        r_tm = (n2*cosTi - n1*cosTt)/(n1*cosTt + n2*cosTi); // paralel
488
 
489
                                        if(dbg) printf("  r_te = %lf, r_tm = %lf\n", r_te, r_tm);
490
 
491
                                        v_tm_t = -v_te.Cross(transmit);
492
                                        v_tm_t = v_tm_t.Unit();
493
                                        pol_t = a_te * (1.0 -  TMath::Abs(r_te)) * v_te  +  a_tm * (1.0 -  TMath::Abs(r_tm)) * v_tm_t;
494
 
495
                                        if(dbg) {
496
                                                printf("  v_tm_t = "); printv(v_tm_t);
497
                                                printf("  pol_t = "); printv(pol_t);
498
                                        }
499
 
500
                                        R_te = TMath::Power(r_te, 2);
501
                                        R_tm = TMath::Power(r_tm, 2);
502
                                        R_f = a_te*a_te*R_te + a_tm*a_tm*R_tm;
503
 
54 f9daq 504
                                        if (dbg) printf("  R_te = %lf, R_tm = %lf, R_f = %lf\n", R_te, R_tm, R_f);
505
                                //}
506
 
25 f9daq 507
                                if(p_ref >= R_f) { // se lomi
54 f9daq 508
                                  if (dbg) printf("   SURFACE REFRACTED. Return.\n");
25 f9daq 509
                                        out->Set(intersect, transmit);
510
                                        out->SetPolarization(pol_t);
54 f9daq 511
                                        return REFRACTION;
25 f9daq 512
                                } else { // se odbije
54 f9daq 513
                                  if (dbg) printf("   SURFACE REFLECTED. p_ref=%f, R_f=%f\n", p_ref, R_f);
25 f9daq 514
                                        transmit = in.GetN() + sign_n*2*cosTi*n;
515
                                        out->Set(intersect, transmit);
516
                                        pol_t = -in.GetP() + sign_n*2*cosTi*n;
517
                                        out->SetPolarization(pol_t);
54 f9daq 518
                                        return REFLECTION;
519
                                }      
520
 
521
                        //}
25 f9daq 522
                        break;
54 f9daq 523
 
25 f9daq 524
                // ----------------------------------------------------------------------------
525
                // --------------- reflection at "reflection" probability ---------------------
526
                // ----------------------------------------------------------------------------
527
                case SURF_REFLE:
528
                        p_ref = rand.Uniform(0.0, 1.0);
529
                        if(p_ref < reflection) { // se odbije
530
                                cosTi = in.GetN() * n;
531
                                transmit = in.GetN() - 2*cosTi*n;
532
                                out->Set(intersect, transmit);
54 f9daq 533
                                return REFLECTION; //sdhfvjhsdbfjhsdbcvjhsb
25 f9daq 534
                        } else { // se ne odbije
535
                                transmit = in.GetN();
536
                                out->Set(intersect, transmit);         
54 f9daq 537
                                return ABSORBED;
25 f9daq 538
                        }                                              
539
                        break;
540
 
541
                // total reflection from n1 to n2 with R probbability
542
                case SURF_IMPER:
543
                        p_ref = rand.Uniform(0.0, 1.0);        
544
                        if(p_ref < reflection) { // se odbije
545
                                cosTi = in.GetN() * n;                 
546
                                if(TMath::Abs(cosTi) < cosTtotal) { // totalni odboj
547
                                        transmit = in.GetN() - 2*cosTi*n;
548
                                        out->Set(intersect, transmit);
549
                                } else { // ni tot. odboja
550
                                        transmit = in.GetN();
551
                                        out->Set(intersect, transmit);
54 f9daq 552
                                        return ABSORBED;
25 f9daq 553
                                }
554
                        } else { // se ne odbije
555
                                transmit = in.GetN();
556
                                out->Set(intersect, transmit);                 
54 f9daq 557
                                return ABSORBED;
25 f9daq 558
                        }                                              
559
                        break;
560
 
561
                default:
562
                        *out = in;
563
                        break;
564
        }
565
 
54 f9daq 566
        return REFRACTION;
25 f9daq 567
}
568
//=================================================================================
569
 
570
 
571
//=================================================================================
572
Guide::Guide(TVector3 center0, DetectorParameters &parameters)
573
{
574
double t;
575
 
576
        TDatime now; rand.SetSeed(now.Get());
577
 
578
        center = center0;
579
        double b = parameters.getB();
580
        double a = parameters.getA();
581
        _d = parameters.getD();
582
        n1 = parameters.getN1();
583
        n2 = parameters.getN2();
584
        // if PlateOn, then n0 = n3 (optical grease), else = n1 (air)
585
        double n0 = (parameters.getPlateOn() ? parameters.getN3(): n1);
586
        n3 = parameters.getN3();
587
        _r = c_reflectivity;
588
        int fresnel = parameters.getFresnel();
589
 
590
        t = b/2.0;
591
        vodnik_edge[0].SetXYZ(0.0, t,-t); vodnik_edge[1].SetXYZ(0.0, t, t);
592
        vodnik_edge[2].SetXYZ(0.0,-t, t); vodnik_edge[3].SetXYZ(0.0,-t,-t);
593
        t = a/2.0;
594
        vodnik_edge[4].SetXYZ(_d, t,-t); vodnik_edge[5].SetXYZ(_d, t, t);
595
        vodnik_edge[6].SetXYZ(_d,-t, t); vodnik_edge[7].SetXYZ(_d,-t,-t);
596
 
597
        for(int i = 0; i<8; i++) vodnik_edge[i] += center;
598
 
599
        s_side[0] = new CSurface(SURF_REFRA, vodnik_edge, n0, n2, _r);
600
        s_side[0]->FlipN();
601
 
602
        s_side[1] = new CSurface(SURF_REFRA, vodnik_edge[3], vodnik_edge[2], vodnik_edge[6], vodnik_edge[7], n2, n1, _r);
603
        s_side[2] = new CSurface(SURF_REFRA, vodnik_edge[2], vodnik_edge[1], vodnik_edge[5], vodnik_edge[6], n2, n1, _r);
604
        s_side[3] = new CSurface(SURF_REFRA, vodnik_edge[1], vodnik_edge[0], vodnik_edge[4], vodnik_edge[5], n2, n1, _r);
605
        s_side[4] = new CSurface(SURF_REFRA, vodnik_edge[0], vodnik_edge[3], vodnik_edge[7], vodnik_edge[4], n2, n1, _r);
606
 
607
        s_side[5] = new CSurface(SURF_REFRA, &vodnik_edge[4], n2, n3, _r);
608
        s_side[5]->FlipN();
609
 
610
        if(fresnel) for(int i=0; i<6; i++) s_side[i]->SetFresnel(1);
611
 
612
        hfate = (TH1F*)gROOT->FindObject("hfate"); if(hfate) delete hfate;
613
        hfate = new TH1F("hfate", "Ray fate", 8, -3.5, 4.5);
614
        (hfate->GetXaxis())->SetBinLabel(1, "Back Ref");
54 f9daq 615
        (hfate->GetXaxis())->SetBinLabel(2, "No Ref");
616
        (hfate->GetXaxis())->SetBinLabel(3, "Refrac");
25 f9daq 617
        (hfate->GetXaxis())->SetBinLabel(4, "LG Miss");
618
        (hfate->GetXaxis())->SetBinLabel(5, "Exit");
619
        (hfate->GetXaxis())->SetBinLabel(6, "Enter");
620
        (hfate->GetXaxis())->SetBinLabel(7, "Rays");
621
        (hfate->GetXaxis())->SetBinLabel(8, "Absorb");
622
 
623
        hnodb_all = (TH1F*)gROOT->FindObject("hnodb_all"); if(hnodb_all) delete hnodb_all;
54 f9daq 624
        hnodb_all = new TH1F("hnodb_all", "N reflected", MAX_REFLECTIONS, -0.5, MAX_REFLECTIONS-0.5);
25 f9daq 625
 
626
        hnodb_exit = (TH1F*)gROOT->FindObject("hnodb_exit"); if(hnodb_exit) delete hnodb_exit;
54 f9daq 627
        hnodb_exit = new TH1F("hnodb_exit", "N reflected and exit", MAX_REFLECTIONS, -0.5, MAX_REFLECTIONS-0.5);
25 f9daq 628
 
629
        int nBins = nch + 1;
630
        hin = (TH2F*)gROOT->FindObject("hin"); if(hin) delete hin;
631
        hin = new TH2F("hin", "Guide entrance window", nBins, -b/2.0, +b/2.0, nBins, -b/2.0, +b/2.0);
632
 
633
        hout = (TH2F*)gROOT->FindObject("hout"); if(hout) delete hout;
634
        hout = new TH2F("hout", "Guide exit window", nBins, -a/2.0, +a/2.0, nBins, -a/2.0, +a/2.0);
635
 
636
        absorption = 0;
637
        A = 0;
638
}
639
//-----------------------------------------------------------------------------
640
// Sledi zarku skozi vodnik. Vrne:                                             
641
//  0, ce zgresi vstopno ploskev                                               
642
//  1, ce zadane izstopno ploskev                                              
643
// -1, ce se v vodniku ne odbije totalno 
644
//  2, enter the light guide, bin 2 of hfate = refraction                                     
645
// -2, ce se ne odbije zaradi koncnega R stranic                               
646
// -3, ce se odbije nazaj in gre nazaj ven skozi sprednjo ploskev              
647
// +4, ce se absorbira v materialu                                             
648
Fate Guide::PropagateRay(CRay in, CRay *out, int *n_points, TVector3 *points)
649
{
650
  if (dbg) printf("--- GUIDE::PropagateRay ---\n");
651
  CRay ray0;
652
  CRay ray1;
653
  TVector3 vec0, vec1;
654
  int inters_i = 0;
655
 
656
        ray0 = in;
657
        int n_odb = 0;
658
        int last_hit = 0;
659
        int propagation = 0;
54 f9daq 660
        int result = s_side[0]->PropagateRay(ray0, &ray1, &vec1);
661
        if( !(result) ) {
25 f9daq 662
                // ce -NI- presecisca z vstopno
663
                if (dbg) printf("  GUIDE: missed the light guide\n");
664
                fate = missed;
54 f9daq 665
                //hfate->Fill(0);
666
        } else if(result == REFLECTION) {
667
          if (dbg) printf(" REFLECTED on the entry surface!\n");
668
          fate = backreflected;
669
          //hfate->Fill(-3);
25 f9daq 670
        } else {
671
          if (dbg) printf("  GUIDE: ray entered\n");
672
                points[0] = ray1.GetR();
673
                hfate->Fill(2); // enter
674
                hin->Fill(vec1.y(), vec1.z());
675
                if (dbg) printf("  GUIDE: n_odb = %d\n", n_odb);
54 f9daq 676
 
677
                while (n_odb++ < MAX_REFLECTIONS) {
25 f9daq 678
                  if (dbg) printf("  GUIDE: Boundary test: %d\n",n_odb);
679
                        ray0 = ray1;
680
                        vec0 = vec1;
681
                        propagation = 11;
682
                        for(inters_i=0; inters_i<6; inters_i++) {
683
                          if (dbg) printf("  GUIDE: Test intersection with surface %d \n", inters_i);
684
                                if( inters_i != last_hit) {
685
                                  int testBoundary = s_side[inters_i]->TestIntersection(&vec1, ray1);
686
                                        if( testBoundary ) {
687
                                          if (dbg) printf("  GUIDE: ray intersects with LG surface %d\n",inters_i);
688
                                                break;
689
                                  }
690
                                }  
691
                        }                      
692
                        points[n_odb] = vec1;
693
                        if(inters_i == 0) {
694
                          fate = backreflected;
54 f9daq 695
                          //hfate->Fill(backreflected); 
25 f9daq 696
                          break;
697
                        } // backreflection
54 f9daq 698
 
699
                        // the passage is possible, test propagation                    
25 f9daq 700
                        propagation = s_side[inters_i]->PropagateRay(ray0, &ray1, &vec1);
54 f9daq 701
 
25 f9daq 702
                        if (dbg) printf("  GUIDE: surface = %d, propagation = %d\n", inters_i, propagation);
54 f9daq 703
 
704
                        if(propagation == REFRACTION) {
705
                          fate = refracted;
706
                          n_odb++;
707
                          points[n_odb] = vec1;
708
                          ray0 = ray1;
709
                          break;
710
                        } // no total reflection when should be
711
                        if(propagation == ABSORBED) {
712
                          fate = noreflection;
713
                          break;
714
                        } //refraction due to finite reflectivity
715
 
25 f9daq 716
                        if(inters_i == 5) { // successfull exit
717
                        // check on which side the vector is?
718
                          TVector3 ray = ray1.GetN();
719
                          TVector3 exitNormal = s_side[5]->GetN();
720
                          //printf("theta(ray) = %lf, theta(normal5) = %lf ", ray.Theta()*DEGREE, exitNormal.Theta()*DEGREE);
721
                          //printf("phi(ray) = %lf, phi(normal5) = %lf\n", ray.Phi()*DEGREE, exitNormal.Phi()*DEGREE);
722
                          if (dbg) printf("ray*n_5 = %lf\n", ray*exitNormal);
723
                          if (ray*exitNormal > 0) {
724
                            if (dbg) printf("  GUIDE: ray is backreflected from exit window.\n");
725
                            fate = backreflected;
726
                            n_odb++;
727
                            points[n_odb] = vec1;
728
                            ray0 = ray1;
729
                            break;
730
                          }
54 f9daq 731
                                fate =  hitExit;
25 f9daq 732
                                hout->Fill(vec1.y(), vec1.z());
733
                                hnodb_exit->Fill(n_odb-1);
734
                                n_odb++;
735
                                points[n_odb] = vec1;
736
                                ray0 = ray1;
737
                                break;
738
                        }
739
                        last_hit = inters_i;
740
                }      
741
        }
742
 
743
        //--- material absorption ---
744
        if(absorption) {
745
                double travel = 0.0;
746
                printf("n_odb = %d\n", n_odb); //dbg   
747
                for(int point = 0; point < n_odb-1; point++) {
748
                        travel += (points[point] - points[point+1]).Mag();
749
                        printf("travel = %lf\n", travel); //dbg   
750
                }
751
                double T_abs = TMath::Exp(-travel/A);
752
                printf("T_abs = %lf\n", T_abs); //dbg   
753
                double p_abs = rand.Uniform(0.0, 1.0); 
754
                printf("p_abs = %lf\n", p_abs); //dbg   
755
 
756
                if(p_abs > T_abs) fate = absorbed; // absorption
757
        }
758
        //--- material absorption ---
759
 
760
        hfate->Fill(fate);
54 f9daq 761
        hfate->Fill(rays);
25 f9daq 762
        hnodb_all->Fill(n_odb-2);
763
        *n_points = n_odb+1;
764
        *out = ray0;
765
        return fate;
766
}
767
//-----------------------------------------------------------------------------
768
void Guide::GetVFate(int *out)
769
{
770
        for(int i=0;i<7;i++) out[i] = (int)hfate->GetBinContent(i+1);
771
}
772
//-----------------------------------------------------------------------------
773
void Guide::Draw(int color, int width)
774
{
775
        for(int i = 0; i<6; i++) s_side[i]->Draw(color, width);
776
}
777
//-----------------------------------------------------------------------------
778
void Guide::DrawSkel(int color, int width)
779
{
780
        TPolyLine3D *line3d = new TPolyLine3D(2);
781
        line3d->SetLineWidth(width); line3d->SetLineColor(color);
782
 
783
        for(int i=0; i<4; i++) {
784
                line3d->SetPoint(0, vodnik_edge[i+0].x(), vodnik_edge[i+0].y(), vodnik_edge[i+0].z());
785
                line3d->SetPoint(1, vodnik_edge[i+4].x(), vodnik_edge[i+4].y(), vodnik_edge[i+4].z());
786
                line3d->DrawClone();
787
        }
788
}
789
//=================================================================================
790
 
791
//=================================================================================
792
int CPlaneR::TestIntersection(TVector3 *vec, CRay ray)
793
{
794
        double num, den; //stevec, imenovalec
795
        double t;
796
        TVector3 tmp;
797
 
798
        if(dbg) printf("---> CPlaneR::TestIntersection <---\n");
799
        if(dbg) {printf("c = "); printv(center); printf(" | n = "); printv(n); printf("\n");}
800
 
801
        double D = - n*center; 
802
        num = n*ray.GetR() + D;
803
        den = n*ray.GetN();
804
 
805
        if(dbg) printf("D = %.4lf | num = %.4lf | den = %.4lf\n", D, num, den);
806
 
807
        if(TMath::Abs(den) < MARGIN) {
808
                if(TMath::Abs(num) < MARGIN)
809
                        return 0;
810
                else
811
                        return 0;
812
        }
813
 
814
        t = num / den;
815
 
816
        if(dbg) printf("t = %.4lf | ", t);
817
 
818
        tmp = ray.GetR();
819
        tmp -= t*ray.GetN();
820
        *vec = tmp;
821
 
822
        if(dbg) {printv(tmp); printf(" | Rv = %.4lf <> R = %.4lf\n", ((tmp - center).Mag()), _r);}
823
 
824
 
825
        if( ((tmp - center).Mag()) < _r )
826
                return 1;
827
        else
828
                return 0;
829
}
830
//-----------------------------------------------------------------------------
831
void CPlaneR::Draw(int color, int width)
832
{
833
        const int NN = 32;     
834
        double phi, x, y;
835
 
836
        TPolyLine3D *arc;
837
        arc = new TPolyLine3D(NN+1);
838
        arc->SetLineWidth(width);
839
        arc->SetLineColor(color);
840
 
841
        for(int i=0; i<=NN; i++) {
842
                phi = i*2.0*TMath::Pi()/NN;
843
                x = _r*TMath::Cos(phi);
844
                y = _r*TMath::Sin(phi);
845
                arc->SetPoint(i, center.x(),  x,  y);
846
        }
847
        arc->Draw();
848
}
849
//=================================================================================
850
 
851
 
852
//=================================================================================
853
CDetector::CDetector(TVector3 center0, DetectorParameters& parameters) :
854
  center(center0),
855
  glass_on(parameters.getGlassOn()),
856
  glass_d(parameters.getGlassD()),
857
  //x_gap(parameters.getGap().X()),
858
  //y_gap(parameters.getGap().Y()),
859
  //z_gap(parameters.getGap().Z()),
860
        //glass(new CSurface),
861
        //glass_circle(new CPlaneR),
862
        //active(new CPlane4),
863
        col_in(2),
864
        col_lg(8),
865
        col_out(4),
866
        col_rgla(6),
867
        col_LG(1),
868
        col_glass(4),
869
        col_active(7),
870
        guide_on(parameters.getGuideOn()),
871
        //window_R( parameters.getB() ),
872
        //window_d(0),
873
  guide(new Guide(center0, parameters)),
874
  plate(new Plate(parameters)),
875
  _plateWidth(parameters.getPlateWidth()),
54 f9daq 876
  _plateOn(parameters.getPlateOn()),
877
  offsetY(parameters.getOffsetY()),
878
  offsetZ(parameters.getOffsetZ())
25 f9daq 879
  {
880
//  };
881
 
882
//-----------------------------------------------------------------------------
883
//void CDetector::Init()
884
//{
885
  double d = parameters.getD();
886
        double x_offset;
887
        if(guide_on) x_offset = center.x();
888
        else x_offset = center.x() - d;
889
 
890
        //guide = new CVodnik(center, SiPM, M, d, type_in, type_side, type_out, n1, n2, n3, reflectivity, fresnel, absorption, A);
891
 
892
        double b = parameters.getB();
893
        double n1 = parameters.getN1();
894
        double n2 = parameters.getN2();
895
        //double n3 = parameters.getN3();
896
        double reflectivity = c_reflectivity; // for faster simulation, not using Fresnel eqs.
897
        double x_gap = parameters.getGap().X();
898
        double y_gap = parameters.getGap().Y();
899
        double z_gap = parameters.getGap().Z();
900
 
901
        TVector3 plane_v[4];
902
        int nBins = nch + 1;
903
        double p_size = b/2.0;
904
        plane_v[0].SetXYZ(x_offset+d+glass_d, y_gap + p_size, z_gap - p_size);
905
        plane_v[1].SetXYZ(x_offset+d+glass_d, y_gap + p_size, z_gap + p_size);
906
        plane_v[2].SetXYZ(x_offset+d+glass_d, y_gap - p_size, z_gap + p_size);
907
        plane_v[3].SetXYZ(x_offset+d+glass_d, y_gap - p_size, z_gap - p_size);
908
        glass = new CSurface(SURF_REFRA, plane_v, n2, n1, reflectivity); glass->FlipN();
909
 
910
        glass_circle = new CPlaneR(TVector3(x_offset+d+glass_d, y_gap, z_gap), TVector3(-1.0, 0.0, 0.0), b);
911
 
912
        hglass = (TH2F*)gROOT->FindObject("hglass"); if(hglass) delete hglass;
913
        hglass = new TH2F("hglass", "Hits glass",
914
                                          nBins, y_gap - p_size, y_gap + p_size,
915
                  nBins, z_gap - p_size, z_gap + p_size);
916
 
917
 
918
        p_size = parameters.getActive()/2.0;
919
        //cout<<"SiPM active length "<<detectorActive<<endl;
920
        //p_size = 1.0/2.0;
921
        plane_v[0].SetXYZ(x_offset+d+x_gap, y_gap + p_size, z_gap - p_size);
922
        plane_v[1].SetXYZ(x_offset+d+x_gap, y_gap + p_size, z_gap + p_size);
923
        plane_v[2].SetXYZ(x_offset+d+x_gap, y_gap - p_size, z_gap + p_size);
924
        plane_v[3].SetXYZ(x_offset+d+x_gap, y_gap - p_size, z_gap - p_size);
925
        active = new CPlane4(plane_v);
926
 
927
        hactive = (TH2F*)gROOT->FindObject("hactive"); if(hactive) delete hactive;
54 f9daq 928
        //hactive = new TH2F("hactive", "Active area hits", nBins, y_gap - p_size, y_gap + p_size, nBins, z_gap - p_size, z_gap + p_size);
929
        hactive = new TH2F("hactive", "Active area hits", nBins, y_gap - p_size + offsetY, y_gap + p_size + offsetY, nBins, z_gap - p_size + offsetZ, z_gap + p_size + offsetZ);
25 f9daq 930
 
931
        p_size = b/2.0;
932
        //p_size = 2.5;
933
        //p_size = M*0.6;
934
        hlaser = (TH2F*)gROOT->FindObject("hlaser"); if(hlaser) delete hlaser;
54 f9daq 935
        hlaser = new TH2F("hlaser", ";x [mm]; y [mm]", nBins, -p_size+offsetY, p_size+offsetY, nBins, -p_size+offsetZ, p_size+offsetZ);
25 f9daq 936
 
937
 
54 f9daq 938
        p_size = 1.4*b/2.0;
25 f9daq 939
        plane_v[0].SetXYZ(x_offset+d+x_gap, y_gap + p_size, z_gap - p_size);
940
        plane_v[1].SetXYZ(x_offset+d+x_gap, y_gap + p_size, z_gap + p_size);
941
        plane_v[2].SetXYZ(x_offset+d+x_gap, y_gap - p_size, z_gap + p_size);
942
        plane_v[3].SetXYZ(x_offset+d+x_gap, y_gap - p_size, z_gap - p_size);
943
        detector = new CPlane4(plane_v);
944
 
945
        hdetector = (TH2F*)gROOT->FindObject("hdetector"); if(hdetector) delete hdetector;
54 f9daq 946
        //hdetector = new TH2F("hdetector", "Hits detector plane", nBins, y_gap - p_size, y_gap + p_size, nBins, z_gap - p_size, z_gap + p_size);
947
        hdetector = new TH2F("hdetector", ";x [mm]; y [mm]", nBins, y_gap-p_size + offsetY, y_gap + p_size + offsetY, nBins, z_gap - p_size + offsetZ, z_gap + p_size + offsetZ);
25 f9daq 948
 
949
        /*
950
        window_circle = new CPlaneR(TVector3(x_offset+d+window_d, y_gap, z_gap), TVector3(-1.0, 0.0, 0.0), window_R);  
951
 
952
        p_size = M*a;
953
        plane_v[0].SetXYZ(x_offset+d+window_d, y_gap + p_size, z_gap - p_size);
954
        plane_v[1].SetXYZ(x_offset+d+window_d, y_gap + p_size, z_gap + p_size);
955
        plane_v[2].SetXYZ(x_offset+d+window_d, y_gap - p_size, z_gap + p_size);
956
        plane_v[3].SetXYZ(x_offset+d+window_d, y_gap - p_size, z_gap - p_size);
957
        window = new CSurface(SURF_REFRA, plane_v, n1, n2, reflectivity); window->FlipN();
958
 
959
        hwindow = (TH2F*)gROOT->FindObject("hwindow"); if(hwindow) delete hwindow;
960
        hwindow = new TH2F("hwindow", "Hits Window", nch, y_gap - window_R, y_gap + window_R, nch, z_gap - window_R, z_gap + window_R);
961
        */
962
        p_size = b/2.0;
963
        histoPlate = (TH2F*)gROOT->FindObject("histoPlate"); if(histoPlate) delete histoPlate;
964
        histoPlate = new TH2F("histoPlate", "Hits on glass plate", nBins, -p_size, +p_size, nBins, -p_size, +p_size);
965
}
54 f9daq 966
 
25 f9daq 967
//-----------------------------------------------------------------------------
968
// vrne 1 ce je zadel aktvino povrsino
969
// vrne <1 ce jo zgresi
970
int CDetector::Propagate(CRay in, CRay *out, int draw)
54 f9daq 971
 
25 f9daq 972
// Sledi zarku skozi vodnik. Vrne:                                             
973
//  0, ce zgresi vstopno ploskev MISSED                                              
974
//  1, ce zadane izstopno ploskev HIT                                             
975
// -1, ce se v vodniku ne odbije totalno REFRACTED
976
//  2, enter the light guide, bin 2 of hfate EXIT                                     
977
// -2, ce se ne odbije zaradi koncnega R stranic - no total reflection REFRACTED                             
978
// -3, ce se odbije nazaj in gre nazaj ven skozi sprednjo ploskev BACK_REFLECTED             
979
// +4, ce se absorbira v materialu ABSORBED
980
{
981
  if (dbg) printf("--- Detector::Propagate ---\n");
982
        //CRay *ray0 = new CRay; ray0->Set(in.GetR(), in.GetN()); ray0->SetColor(col_in);
983
        CRay *rayin = new CRay(in);
984
        rayin->SetColor(col_in);
985
        CRay *rayout = new CRay;
986
 
987
        const int max_n_points = guide->GetMAXODB() + 2;
988
        TVector3 pointsPlate[max_n_points];
989
        //TVector3 intersection;
54 f9daq 990
        Fate fatePlate;
25 f9daq 991
        int nPointsPlate;
992
        TPolyLine3D *line3d = new TPolyLine3D(2);
993
        line3d->SetLineWidth(1);
994
        line3d->SetLineColor(4);
995
 
996
        // Draw the plate and propagate the ray through
997
        // check if the ray should be reflected??
998
 
999
        if(_plateOn) {
1000
 
1001
          fatePlate = plate->propagateRay(*rayin, rayout, &nPointsPlate, pointsPlate);
1002
          if(draw) rayin->DrawS(center.x()- _plateWidth, -10.0);
1003
 
1004
          if(draw) {
54 f9daq 1005
            if(fatePlate == missed) {
25 f9daq 1006
              rayout->SetColor(col_in);
1007
              rayout->DrawS(center.x() - _plateWidth, -10.0);
1008
              }
1009
              else {
1010
                int p_i;
1011
                for(p_i = 0; p_i < nPointsPlate-1; p_i++) {
1012
                                          line3d->SetPoint(0, pointsPlate[p_i].x(), pointsPlate[p_i].y(), pointsPlate[p_i].z());
1013
                                          line3d->SetPoint(1, pointsPlate[p_i+1].x(), pointsPlate[p_i+1].y(), pointsPlate[p_i+1].z());
1014
                                          line3d->DrawClone();
1015
                                  }
54 f9daq 1016
                                  if(fatePlate != noreflection) {
25 f9daq 1017
                                          //rayout->SetColor(7);
1018
                                          rayout->DrawS(pointsPlate[p_i].x(), -0.1);
1019
                                  }
1020
                          }
1021
            }
1022
 
54 f9daq 1023
          if(! (fatePlate == hitExit or fatePlate == refracted) ) {
25 f9daq 1024
 
54 f9daq 1025
                  if (dbg) printf("CDetector::propagate Simulated ray missed the entry surface!\n");
25 f9daq 1026
                        return fatePlate;
1027
                }
54 f9daq 1028
                // missing: if refracted at plate sides
1029
                //if (fatePlate == refracted) return fatePlate;
25 f9daq 1030
                histoPlate->Fill(pointsPlate[0].y(), pointsPlate[0].z()); // entry point
1031
         }
1032
         else {
1033
          rayout = rayin;
1034
          if(draw) rayout->DrawS(center.x(), -10.0);
1035
          }
1036
 
1037
        // If the ray is not reflected in the plate
1038
        // Draw the light guide and propagate the ray through
1039
 
1040
        //const int max_n_points = guide->GetMAXODB() + 2;
1041
        TVector3 points[max_n_points];
1042
        TVector3 presecisce;
54 f9daq 1043
 
25 f9daq 1044
        int n_points;
1045
        int fate_glass;
1046
        CRay *ray0 = new CRay(*rayout);
1047
        // delete rayout; -> creates dangling reference when tries to delete ray0!
54 f9daq 1048
        delete rayin;
25 f9daq 1049
        CRay *ray1 = new CRay;
1050
 
1051
        fate = guide->PropagateRay(*ray0, ray1, &n_points, points);
54 f9daq 1052
        if (dbg) {
1053
          if (fate == backreflected) printf("DETECTOR::backreflected\n");
1054
          }
25 f9daq 1055
 
1056
        line3d->SetLineColor(col_lg);  
1057
        int p_i;
1058
        if(guide_on) {
1059
                if(draw) {
1060
                        if(fate == missed) {
1061
                          if (dbg) printf("Detector: fate=missed\n");
1062
                          TVector3 r = ray1->GetR();
1063
                          TVector3 n = ray1->GetN();
1064
                          ray1->Set(r,n);
1065
                                ray1->DrawS(center.x(), 10.0);
1066
                        } else {
1067
                                for(p_i = 0; p_i < n_points-1; p_i++) {
1068
                                        line3d->SetPoint(0, points[p_i].x(), points[p_i].y(), points[p_i].z());
1069
                                        line3d->SetPoint(1, points[p_i+1].x(), points[p_i+1].y(), points[p_i+1].z());
1070
                                        line3d->DrawClone();
1071
                                }
1072
                                if(fate != noreflection) {
1073
                                  if (dbg) printf("Detector: fate != noreflection, fate = %d\n", (int)fate);
1074
                                        if(glass_on) {/*if(fate == 1)*/ ray1->Draw(points[p_i].x(), center.x() + guide->getD() + glass_d);}
1075
                                        else {
1076
                                        ray1->SetColor(col_out);
1077
                                        ray1->DrawS(points[p_i].x(), 10.0);
1078
                                        }
1079
                                }
1080
                        }
1081
                }
1082
 
1083
 
54 f9daq 1084
                if(! (fate == hitExit or fate == refracted) ) {
1085
                  if (dbg) printf("Detector: fate != hit, refracted\n");
25 f9daq 1086
                        *out = *ray1;
1087
                        return fate;
1088
                }
1089
        } else {
1090
          if (dbg) printf("Detector: fate = hit or refracted");
1091
                ray1 = ray0;
1092
                if(draw) {
1093
                        if(glass_on) ray1->Draw(center.x(), center.x() /*+ window_d*/);
1094
                        else ray1->DrawS(center.x(), 10.0);
1095
                }
1096
        }
1097
 
1098
/*     
1099
        TVector3 pres_wind;
1100
        fate = window_circle->TestIntersection(&pres_wind, *ray1);
1101
        if(fate == 1) {
1102
                hwindow->Fill(pres_wind.y(), pres_wind.z());
1103
 
1104
                if(!guide_on) {
1105
                        window->PropagateRay(*ray0, ray1, &presecisce);
1106
                        if(draw) ray1->Draw(center.x() + window_d, center.x() + glass_d);
1107
                        *ray0 = *ray1;
1108
                }
1109
        */     
1110
                fate = missed; // zgresil aktivno povrsino
1111
                if(glass_on) {                 
1112
                        *ray0 = *ray1; ray1->SetColor(col_rgla);
1113
                        fate_glass = glass->PropagateRay(*ray0, ray1, &presecisce);
1114
                        if(fate_glass == 1) {
1115
                                hglass->Fill(presecisce.y(), presecisce.z());
1116
                                if(draw) ray1->DrawS(presecisce.x(), 10.0);
1117
                                if(active->TestIntersection(&presecisce, *ray1)) {
54 f9daq 1118
                                        fate = hitExit;
1119
                                        hactive->Fill(offsetY + presecisce.y(), offsetZ + presecisce.z());
25 f9daq 1120
                                        hlaser->Fill((in.GetR()).y(), (in.GetR()).z());
1121
                                }
1122
                                if(detector->TestIntersection(&presecisce, *ray1))
1123
                                        hdetector->Fill(presecisce.y(), presecisce.z());
1124
                        } else if(fate_glass == 2) {
1125
                                if(draw) ray1->DrawS(presecisce.x(), 10.0);
1126
                        }
1127
                } else {
54 f9daq 1128
                  // Main test: ray and SiPM surface
25 f9daq 1129
                        if(active->TestIntersection(&presecisce, *ray1)) {
54 f9daq 1130
                                fate = hitExit;
1131
                                hactive->Fill(offsetY + presecisce.y(), offsetZ + presecisce.z());
1132
                                hlaser->Fill((in.GetR()).y() + offsetY, (in.GetR()).z() + offsetZ);
25 f9daq 1133
                        }
54 f9daq 1134
                        // If it is on the same plane as SiPM
25 f9daq 1135
                        if(detector->TestIntersection(&presecisce, *ray1))
54 f9daq 1136
                                hdetector->Fill(offsetY + presecisce.y(), offsetZ + presecisce.z());
25 f9daq 1137
                }
1138
        //} else {
1139
                //if(draw) ray1->Draw(presecisce.x(), center.x()+d+window_d);
1140
        //}
1141
 
1142
        *out = *ray1;
54 f9daq 1143
        delete ray0;
1144
        delete ray1;
1145
        delete rayout;
25 f9daq 1146
        return fate;
1147
}
1148
//-----------------------------------------------------------------------------
1149
void CDetector::Draw(int width)
1150
{
1151
        if(guide_on) {
1152
                if( TMath::Abs(guide->getN1()-guide->getN2()) < MARGIN ) {
1153
                  if(_plateOn) plate->drawSkel(col_LG, width);
1154
                        guide->DrawSkel(col_LG, width);
1155
                        }
1156
                else {
1157
                  if(_plateOn) plate->draw(4, width);
1158
                        guide->Draw(col_LG, width);
1159
                        }
1160
        }
1161
 
1162
        if(glass_on) glass_circle->Draw(col_glass, width);
1163
        //window_circle->Draw(col_glass, width);
1164
        active->Draw(col_active, width);
1165
}
1166
//=================================================================================
1167
 
1168
Plate::Plate(DetectorParameters& parameters)
1169
{
1170
  TVector3 center = CENTER;
1171
  const double b = parameters.getB();
1172
  const double n1 = parameters.getN1();
1173
        const double n2 = parameters.getN2();
1174
        const double n3 = parameters.getN3();
1175
        const double reflectivity = c_reflectivity;
1176
        const double t = b/2.;
1177
        const double plateWidth = parameters.getPlateWidth();
1178
        center.SetX( CENTER.X() - plateWidth );
1179
        plate_edge[0].SetXYZ(0.0, t,-t); plate_edge[1].SetXYZ(0.0, t, t);
1180
        plate_edge[2].SetXYZ(0.0,-t, t); plate_edge[3].SetXYZ(0.0,-t,-t);
1181
        plate_edge[4].SetXYZ(plateWidth, t,-t); plate_edge[5].SetXYZ(plateWidth, t, t);
1182
        plate_edge[6].SetXYZ(plateWidth,-t, t); plate_edge[7].SetXYZ(plateWidth,-t,-t);
1183
 
1184
        for(int i = 0; i<8; i++) plate_edge[i] += center;
1185
 
1186
        sides[0] = new CSurface(SURF_REFRA, plate_edge, n1, n2, reflectivity);
1187
        sides[0]->FlipN();
1188
 
1189
        sides[1] = new CSurface(SURF_REFRA, plate_edge[3], plate_edge[2], plate_edge[6], plate_edge[7], n2, n2, reflectivity);
1190
        sides[2] = new CSurface(SURF_REFRA, plate_edge[2], plate_edge[1], plate_edge[5], plate_edge[6], n2, n2, reflectivity);
1191
        sides[3] = new CSurface(SURF_REFRA, plate_edge[1], plate_edge[0], plate_edge[4], plate_edge[5], n2, n2, reflectivity);
1192
        sides[4] = new CSurface(SURF_REFRA, plate_edge[0], plate_edge[3], plate_edge[7], plate_edge[4], n2, n2, reflectivity);
1193
 
1194
        sides[5] = new CSurface(SURF_REFRA, &plate_edge[4], n2, n3, reflectivity);
1195
        sides[5]->FlipN();
1196
 
1197
        for(int i=0; i<6; i++) sides[i]->SetFresnel(1);
1198
}
1199
 
1200
void Plate::draw(int color, int width)
1201
{
1202
        for(int i = 0; i<6; i++) sides[i]->Draw(color, width);
1203
}
1204
 
1205
void Plate::drawSkel(int color, int width)
1206
{
1207
        TPolyLine3D line3d(2);
1208
        line3d.SetLineWidth(width);
1209
        line3d.SetLineColor(color);
1210
 
1211
        for(int i=0; i<4; i++) {
1212
                line3d.SetPoint(0, plate_edge[i+0].x(), plate_edge[i+0].y(), plate_edge[i+0].z());
1213
                line3d.SetPoint(1, plate_edge[i+4].x(), plate_edge[i+4].y(), plate_edge[i+4].z());
1214
                line3d.DrawClone();
1215
        }
1216
}
1217
 
54 f9daq 1218
Fate Plate::propagateRay(CRay in, CRay *out, int *n_points, TVector3 *points)
25 f9daq 1219
{
1220
  CRay ray0;
1221
  CRay ray1;
1222
  TVector3 vec0, vec1;
54 f9daq 1223
  Fate fate = enter;
25 f9daq 1224
  int inters_i = 0;
1225
 
1226
        ray0 = in;
1227
        int n_odb = 0;
1228
        int last_hit = 0;
1229
        int propagation = 0;
1230
 
54 f9daq 1231
        int result = sides[0]->PropagateRay(ray0, &ray1, &vec1);
1232
        if( !result ) {
25 f9daq 1233
                // ce -NI- presecisca z vstopno
54 f9daq 1234
                fate = missed;
1235
        } else if(result == REFLECTION) {
1236
          if (dbg) printf("PLATE: reflected\n");
1237
          fate = backreflected;
25 f9daq 1238
        } else {
1239
                points[0] = ray1.GetR();
1240
                //hfate->Fill(2);
1241
                //hin->Fill(vec1.y(), vec1.z());
54 f9daq 1242
                while (n_odb++ < MAX_REFLECTIONS) {
25 f9daq 1243
                        ray0 = ray1;
1244
                        vec0 = vec1;
1245
                        propagation = 11;
1246
                        for(inters_i=0; inters_i<6; inters_i++) {
1247
                                if( inters_i != last_hit) {
1248
                                        if( sides[inters_i]->TestIntersection(&vec1, ray1) ) break;
1249
                                        }
1250
                                }                                                      
1251
                        points[n_odb] = vec1;
1252
                        if(inters_i == 0) {
54 f9daq 1253
                          fate = backreflected;
25 f9daq 1254
                          break;} // backreflection
1255
 
1256
                        propagation = sides[inters_i]->PropagateRay(ray0, &ray1, &vec1);
1257
                        if(inters_i == 5) { // successfull exit
54 f9daq 1258
                                fate = hitExit;
25 f9daq 1259
                                //hout->Fill(vec1.y(), vec1.z()); 
1260
                                //hnodb_exit->Fill(n_odb-1); 
1261
                                n_odb++;
1262
                                points[n_odb] = vec1;
1263
                                ray0 = ray1;
1264
                                break;
1265
                        }
1266
                        if(propagation == 1) {
54 f9daq 1267
                          fate = refracted; //at side?
25 f9daq 1268
                          n_odb++;
1269
                          points[n_odb] = vec1;
1270
                          ray0 = ray1;
1271
                          break;} // no total reflection when should be
1272
 
54 f9daq 1273
                        if(propagation == -2) {
1274
                          fate = noreflection;
1275
                          break;
1276
                          } // absorption due to finite reflectivity
1277
 
25 f9daq 1278
                        last_hit = inters_i;
1279
                }      
1280
        }
1281
 
1282
        *n_points = n_odb+1;
1283
        *out = ray0;
1284
        return fate;
1285
};
1286
//=============================================================================================================================== <<<<<<<<
1287
 
1288