Ignore:
Timestamp:
09/02/19 18:33:42 (5 years ago)
Author:
tbretz
Message:
A first simplified implementation of a fresnel lens.
Location:
trunk/Mars/msimreflector
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/Mars/msimreflector/MLens.cc

    r19539 r19592  
    3333#include <errno.h>
    3434
    35 #include <stdlib.h> // atof (Ubuntu 8.10)
    36 
    37 #include <TClass.h>
    38 #include <TSystem.h>
    39 #include <TEllipse.h>
    40 
    4135#include "MQuaternion.h"
    42 #include "MMirror.h"
    4336
    4437#include "MLog.h"
    4538#include "MLogManip.h"
    4639
    47 #include "MH.h"
    48 
    4940ClassImp(MLens);
    5041
     
    6051    fTitle = title ? title : "Parameter container storing a collection of several mirrors (reflector)";
    6152
    62     fMirrors.SetOwner();
    63 }
    64 
    65 // --------------------------------------------------------------------------
    66 //
    67 // Set the SigmaPSF of all mirrors currently stored.
    68 //
    69 void MLens::SetSigmaPSF(Double_t psf)
    70 {
    71     fMirrors.R__FOR_EACH(MMirror, SetSigmaPSF)(psf);
    72 }
    73 
    74 // --------------------------------------------------------------------------
    75 //
    76 // Calculate the maximum radius of th ereflector. This is not meant as
    77 // a precise number but as a rough estimate e.g. to bin a histogram.
    78 //
    79 void MLens::InitMaxR()
    80 {
    81     fMaxR = 0;
    82 
    83     TIter Next(&fMirrors);
    84     MMirror *m = 0;
    85     while ((m=static_cast<MMirror*>(Next())))
    86     {
    87         // Take into account the maximum incident angle 8eg 10deg) and
    88         // the theta-angle of the mirror and the z-distance.
    89         const Double_t r = m->GetDist()+1.5*m->GetMaxR();
    90         if (r > fMaxR)
    91             fMaxR = r;
    92     }
     53    fMaxR = 25;
    9354}
    9455
     
    10061Double_t MLens::GetA() const
    10162{
    102     Double_t A = 0;
    103 
    104     TIter Next(&fMirrors);
    105     MMirror *m = 0;
    106     while ((m=static_cast<MMirror*>(Next())))
    107         A += m->GetA();
    108 
    109     return A;
    110 }
    111 
    112 // --------------------------------------------------------------------------
    113 //
    114 // Get the pointer to the first mirror. This is a very dangerous way of
    115 // access, but the fastest possible. because it is the most often called
    116 // function in ExecuteReflector we have to have a very fast access.
    117 //
    118 const MMirror **MLens::GetFirstPtr() const
    119 {
    120     return (const MMirror**)fMirrors.GetObjectRef(0);
    121 }
    122 
    123 // --------------------------------------------------------------------------
    124 //
    125 // Get number of mirrors. There should be no holes in the array!
    126 //
    127 const UInt_t MLens::GetNumMirrors() const
    128 {
    129     return fMirrors.GetEntriesFast();
     63    return fMaxR*fMaxR*TMath::Pi();
    13064}
    13165
     
    14276}
    14377
    144 // Jeder Spiegel sollte eine Liste aller andern Spiegel in der
    145 // reihenfolge Ihrer Entfernung enthalten. Wir starten mit der Suche
    146 // immer beim zuletzt getroffenen Spiegel!
    147 //
    148 // --------------------------------------------------------------------------
    149 //
    150 // Loops over all mirrors of the reflector. After doing a rough check
    151 // whether the mirror can be hit at all the reflection is executed
    152 // calling the ExecuteMirror function of the mirrors.
    153 //
    154 // If a mirror was hit its index is retuened, -1 otherwise.
    155 //
    156 // FIXME: Do to lopping over all mirrors for all photons this is the
    157 // most time consuming function in teh reflector simulation. By a more
    158 // intelligent way of finding the right mirror then just testing all
    159 // this could be accelerated a lot.
    160 //
    161 Int_t MLens::ExecuteOptics(MQuaternion &p, MQuaternion &u) const
    162 {
    163     //static const TObjArray *arr = &((MMirror*)fMirrors[0])->fNeighbors;
    164 
    165     // This way of access is somuch faster than the program is
    166     // a few percent slower if accessed by UncheckedAt
    167     const MMirror **s = GetFirstPtr();
    168     const MMirror **e = s+GetNumMirrors();
    169     //const MMirror **s = (const MMirror**)fMirrors.GetObjectRef(0);
    170     //const MMirror **e = s+fMirrors.GetEntriesFast();
    171     //const MMirror **s = (const MMirror**)arr->GetObjectRef(0);
    172     //const MMirror **e = s+arr->GetEntriesFast();
    173 
    174     // Loop over all mirrors
    175     for (const MMirror **m=s; m<e; m++)
     78struct Groove
     79{
     80    double fTheta;
     81    double fPeakZ;
     82
     83    double fTanTheta;
     84    double fTanTheta2[3];
     85
     86    Groove() { }
     87    Groove(double theta, double z) : fTheta(theta), fPeakZ(z)
    17688    {
    177         const MMirror &mirror = **m;
    178 
    179         // FIXME: Can we speed up using lookup tables or
    180         //        indexed tables?
    181 
    182         // MirrorShape: Check if this mirror can be hit at all
    183         // This is to avoid time consuming calculation if there is no
    184         // chance of a coincidence.
    185         // FIXME: Inmprove search algorithm (2D Binary search?)
    186         if (!mirror.CanHit(p))
    187             continue;
    188 
    189         // Make a local copy of position and direction which can be
    190         // changed by ExecuteMirror.
    191         MQuaternion q(p);
    192         MQuaternion v(u);
    193 
    194         // Check if this mirror is hit, and if it is hit return
    195         // the reflected position and direction vector.
    196         // If the mirror is missed we go on with the next mirror.
    197         if (!mirror.ExecuteMirror(q, v))
    198             continue;
    199 
    200         // We hit a mirror. Restore the local copy of position and
    201         // direction back into p und u.
    202         p = q;
    203         u = v;
    204 
    205         //arr = &mirror->fNeighbors;
    206 
    207         return m-s;
     89        fTanTheta = tan(theta);
     90
     91        fTanTheta2[0] = fTanTheta*fTanTheta;
     92        fTanTheta2[1] = fTanTheta*fTanTheta * fPeakZ;
     93        fTanTheta2[2] = fTanTheta*fTanTheta * fPeakZ*fPeakZ;
    20894    }
    209 
    210     return -1;
    211 }
    212 
    213 // --------------------------------------------------------------------------
    214 //
    215 // I/O helper for ReadFile to avoid calling "delete arr" before every return
    216 //
    217 MMirror *MLens::EvalTokens(TObjArray &arr, Double_t defpsf) const
    218 {
    219     if (arr.GetEntries()<9)
     95};
     96
     97float CalcIntersection(const MQuaternion &p, const MQuaternion &u, const Groove &g, const double &fThickness)
     98{
     99   // p is the position  in the plane  of the lens (px, py, 0,  t)
     100   // u is the direction of the photon             (ux, uy, dz, dt)
     101
     102   // Assuming a cone with peak at z and opening angle theta
     103
     104   // Radius at distance z+dz from the tip and at distance r from the axis of the cone
     105   // r = (z-dz) * tan(theta)
     106
     107   // Defining
     108   // zc := z+dz
     109
     110   // Surface of the cone
     111   // (X)         ( tan(theta) * sin(phi) )
     112   // (Y)  = zc * ( tan(theta) * cos(phi) )
     113   // (Z)         (          1            )
     114
     115   // Line
     116   // (X)        (u.x)   (p.x)
     117   // (Y) = dz * (u.y) + (p.y)
     118   // (Z)        (u.z)   ( 0 )
     119
     120   // Equality
     121   //
     122   // X^2+Y^2 = r^2
     123   //
     124   // (dz*u.x+p.x)^2 + (dz*u.y+p.y)^2 = (z-dz)^2 * tan(theta)^2
     125   //
     126   // dz^2*ux^2 + px^2 + 2*dz*ux*px + dz^2*uy^2 + py^2 + 2*dz*uy*py = (z^2*tan(theta)^2 + dz^2*tan(theta)^2 - 2*z*dz*tan(theta)^2
     127   //
     128   // dz^2*(ux^2 + uy^2 - tan(theta)^2)  +  2*dz*(ux*px + uy*py + z*tan(theta)^2)  +  (px^2+py^2 - z^2*tan(theta)^2)  =  0
     129
     130   // t   := tan(theta)
     131   // a   := ux^2    + uy^2    -     t^2
     132   // b/2 := ux*px   + uy*py   - z  *t^2 := B
     133   // c   :=    px^2 +    py^2 - z^2*t^2
     134
     135   // dz = [ -b +- sqrt(b^2 - 4*a*c) ] / [2*a]
     136   // dz = [ -B +- sqrt(B^2 -   a*c) ] /    a
     137
     138    const double a = u.R2()                    - g.fTanTheta2[0];
     139    const double B = u.XYvector()*p.XYvector() + g.fTanTheta2[1];
     140    const double c = p.R2()                    - g.fTanTheta2[2];
     141
     142    const double B2 = B*B;
     143    const double ac = a*c;
     144
     145    if (B2<ac)
     146        return 0;
     147
     148    const double sqrtBac = sqrt(B2 - a*c);
     149
     150    if (a==0)
     151        return 0;
     152
     153    const double dz[2] =
    220154    {
    221         *fLog << err << "ERROR - Not enough arguments..." << endl;
    222         return 0;
    223     }
    224 
    225     const TVector3 pos(atof(arr[0]->GetName()),
    226                        atof(arr[1]->GetName()),
    227                        atof(arr[2]->GetName()));
    228 
    229     const TVector3 norm(atof(arr[3]->GetName()),
    230                         atof(arr[4]->GetName()),
    231                         atof(arr[5]->GetName()));
    232 
    233     UInt_t n = 6;
    234 
    235     TString six = arr[n]->GetName();
    236 
    237     Char_t shape = 'S';
    238     if (!six.IsFloat())
    239     {
    240         shape = six[0];
    241         n++;
    242     }
    243 
    244     const Double_t F = atof(arr[n++]->GetName());
    245 
    246     const TString val = arr[n++]->GetName();
    247 
    248     const Double_t psf = val.IsFloat() ? val.Atof() : -1;
    249 
    250     n += val.IsFloat() ? 1 : 0;
    251 
    252     TString type = arr[n-1]->GetName();
    253     type.Prepend("MMirror");
    254 
    255     for (UInt_t i=0; i<n; i++)
    256         delete arr.RemoveAt(i);
    257     arr.Compress();
    258 
    259     TString msg;
    260     TClass *cls = MParContainer::GetClass(type);
    261     if (!cls)
    262     {
    263         *fLog << err << "ERROR - Class " << type << " not in dictionary." << endl;
    264         return 0;
    265     }
    266 
    267     if (!cls->InheritsFrom(MMirror::Class()))
    268     {
    269         *fLog << err << "ERROR - Cannot create new instance of class " << type << ": " << endl;
    270         *fLog << "Class doesn't inherit from MMirror." << endl;
    271         return 0;
    272     }
    273 
    274     MMirror *m = static_cast<MMirror*>(cls->New());
    275     if (!m)
    276     {
    277         *fLog << err << "ERROR - Cannot create new instance of class " << type << ": " << endl;
    278         *fLog << " - Class has no default constructor." << endl;
    279         *fLog << " - An abstract member functions of a base class is not overwritten." << endl;
    280         return 0;
    281     }
    282 
    283     m->SetPosition(pos);
    284     m->SetNorm(norm);
    285     m->SetShape(shape);
    286     m->SetFocalLength(F);
    287     m->SetSigmaPSF(psf>=0 ? psf : defpsf);
    288 
    289     if (m->ReadM(arr)>=0)
    290         return m;
    291 
    292     *fLog << err << "ERROR - " << type << "::ReadM failed." << endl;
    293 
    294     delete m;
     155        (sqrtBac-B) / a,
     156        (sqrtBac+B) / a
     157    };
     158
     159    // Only one dz[i] can be within the lens (due to geometrical reasons
     160
     161    if (dz[0]<0 && dz[0]>fThickness)
     162        return dz[0];
     163
     164    if (dz[1]<0 && dz[1]>fThickness)
     165        return dz[1];
     166
    295167    return 0;
    296168}
    297169
    298 // --------------------------------------------------------------------------
    299 //
    300 // Read a reflector setup from a file. This needs improvemtn.
    301 //
    302 // The file structur is like:
    303 //
    304 //     x y z nx ny nz F [psf] Type ...
    305 //
    306 //  x:         x-coordinate of the mirror's center
    307 //  y:         y-coordinate of the mirror's center
    308 //  z:         z-coordinate of the mirror's center
    309 //  nx:        x-component of the normal vecor in the mirror center
    310 //  ny:        y-component of the normal vecor in the mirror center
    311 //  nz:        z-component of the normal vecor in the mirror center
    312 //  [shape]:   S for spherical <default>, P for parabolic
    313 //  F:         Focal distance of a spherical mirror
    314 //  [psf]:     This number is the psf given in the units of x,y,z and
    315 //             defined at the focal distance F. It can be used to overwrite
    316 //             the second argument given in ReadFile for individual mirrors.
    317 //  Type:      A instance of a mirrot of the class Type MMirrorType is created
    318 //             (Type can be, for example, Hex for for MMirrorHex).
    319 //  ...:       Additional arguments as defined in MMirrorType::ReadM
    320 //
    321 //
    322 // Coordinate System:
    323 //  The coordinate system is local in the reflectors frame.
    324 //  It is defined viewing from the back side of the reflector
    325 //  towards the camera. (x "right", y "up", z from reflector to camera)
    326 //  Note, that it is left-handed!
    327 //
    328 Bool_t MLens::ReadFile(TString fname, Double_t defpsf)
    329 {
    330     SetTitle(fname);
    331     fMirrors.Delete();
    332 
    333     gSystem->ExpandPathName(fname);
    334 
    335     ifstream fin(fname);
    336     if (!fin)
    337     {
    338         *fLog << err << "Cannot open file " << fname << ": ";
    339         *fLog << (errno!=0?strerror(errno):"Insufficient memory for decompression") << endl;
    340         return kFALSE;
    341     }
    342 
    343 /*
    344     Int_t idx[964];
    345     Int_t srt[964];
    346     for (int i=0; i<964; i++)
    347         srt[i] = gRandom->Integer(964);
    348 
    349     TMath::Sort(964, srt, idx);
    350     */
    351 
    352     Int_t cnt = 0;
    353 
    354     while (1)
    355     {
    356         TString line;
    357         line.ReadLine(fin);
    358         if (!fin)
    359             break;
    360 
    361         // Count lines
    362         cnt++;
    363 
    364         // Skip comments
    365         if (line.BeginsWith("#"))
    366             continue;
    367 
    368         // Remove leading and trailing whitespaces
    369         line=line.Strip(TString::kBoth);
    370 
    371         // Skip empty lines
    372         if (line.IsNull())
    373             continue;
    374 
    375         // Tokenize line
    376         TObjArray *arr = line.Tokenize(' ');
    377 
    378         // Evaluate tokens
    379         MMirror *m = EvalTokens(*arr, defpsf);
    380 
    381         // Delete now obsolete array
    382         delete arr;
    383 
    384         // Check if a new mirror could be created successfully
    385         if (!m)
    386         {
    387             *fLog << err << "Error in line " << cnt << ": " << line << endl;
    388             return kFALSE;
    389         }
    390 
    391         // Add new mirror to array
    392         fMirrors.Add(m);
    393     }
    394 
    395     InitMaxR();
    396 
    397     return kTRUE;
    398 
    399 //    fMirrors.Sort();
    400 /*
    401     for (int i=0; i<964; i++)
    402     {
    403         MMirror &ref = (MMirror&)*fMirrors[i];
    404 
    405         TArrayD dist(964);
    406         for (int j=0; j<964; j++)
    407         {
    408             const MMirror &mir = (MMirror&)*fMirrors[j];
    409             dist[j] = (ref-mir).Mod();
    410         }
    411 
    412         TArrayI idx(964);
    413         TMath::Sort(964, dist.GetArray(), idx.GetArray(), kFALSE);
    414 
    415         for (int j=0; j<964; j++)
    416         {
    417             ref.fNeighbors.Add(fMirrors[idx[j]]);
    418         }
    419     }*/
    420 }
    421 
    422 // ------------------------------------------------------------------------
    423 //
    424 Bool_t MLens::WriteFile(TString fname) const
    425 {
    426     gSystem->ExpandPathName(fname);
    427 
    428     ofstream fout(fname);
    429     if (!fout)
    430     {
    431         *fLog << err << "Cannot open file " << fname << ": ";
    432         *fLog << (errno!=0?strerror(errno):"Insufficient memory for decompression") << endl;
    433         return kFALSE;
    434     }
    435 
    436     TIter Next(&fMirrors);
    437     MMirror *m = 0;
    438     while ((m=static_cast<MMirror*>(Next())))
    439     {
    440         //  x:         x-coordinate of the mirror's center
    441         //  y:         y-coordinate of the mirror's center
    442         //  z:         z-coordinate of the mirror's center
    443         //  nx:        x-component of the normal vecor in the mirror center
    444         //  ny:        y-component of the normal vecor in the mirror center
    445         //  nz:        z-component of the normal vecor in the mirror center
    446         //  [shape]:   S for spherical <default>, P for parabolic
    447         //  F:         Focal distance of a spherical mirror
    448         //  [psf]:     This number is the psf given in the units of x,y,z and
    449         //             defined at the focal distance F. It can be used to overwrite
    450         //             the second argument given in ReadFile for individual mirrors.
    451         //  Type:      A instance of a mirrot of the class Type MMirrorType is created
    452         //             (Type can be, for example, Hex for for MMirrorHex).
    453         //  ...:       Additional arguments as defined in MMirrorType::ReadM
    454 
    455         fout << setw(12) << m->X() << " ";
    456         fout << setw(12) << m->Y() << " ";
    457         fout << setw(12) << m->Z() << "  ";
    458         fout << setw(12) << m->Nx() << " ";
    459         fout << setw(12) << m->Ny() << " ";
    460         fout << setw(12) << m->Nz() << " ";
    461         if (m->GetShape()==1)
    462             fout << "P ";
    463         fout << m->GetFocalLength() << " ";
    464         // cout << m->GetSigmaPSF() << " ";
    465         fout << m->ClassName()+7 << " ";
    466         m->WriteM(fout);
    467         fout << endl;
    468     }
    469     return kTRUE;
    470 }
    471 
    472 // --------------------------------------------------------------------------
    473 //
    474 // Print the collection of mirrors
    475 //
    476 void MLens::Print(Option_t *o) const
    477 {
    478     *fLog << all << GetDescriptor() << " (" << GetNumMirrors() << "): " << endl;
    479 
    480     fMirrors.Print(o);
    481 }
    482 
    483 // --------------------------------------------------------------------------
    484 //
    485 // Paint the collection of mirrors
    486 //
    487 void MLens::Paint(Option_t *o)
    488 {
    489     if (!TString(o).Contains("same", TString::kIgnoreCase))
    490         MH::SetPadRange(-fMaxR*1.01, -fMaxR*1.01, fMaxR*1.01, fMaxR*1.01);
    491 
    492     fMirrors.Paint(o);
    493 
    494     TEllipse e;
    495     e.SetFillStyle(0);
    496     e.SetLineColor(17);
    497     e.PaintEllipse(0, 0, fMaxR, fMaxR, 0, 360, 0);
    498 }
    499 
    500 // --------------------------------------------------------------------------
    501 //
    502 // SigmaPSF: -1
    503 // FileName: reflector.txt
    504 //
    505 // SigmaPSF can be used to set a default for the psf of the mirrors
    506 // read from the file. Note, that this can be overwritten for individual
    507 // mirrors in the file. The SigmaPSF is the sigma of a 1D-Gauss fitted
    508 // to the radial profile of the light distribution.
    509 //
    510 // For details on the file structure see MLens::ReadFile
    511 //
    512 Int_t MLens::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
    513 {
    514     Bool_t rc = kFALSE;
    515 
    516     Double_t psf = -1;
    517     if (IsEnvDefined(env, prefix, "SetSigmaPSF", print))
    518     {
    519         rc = kTRUE;
    520         psf = GetEnvValue(env, prefix, "SetSigmaPSF", -1);
    521     }
    522 
    523     if (IsEnvDefined(env, prefix, "FileName", print))
    524     {
    525         rc = kTRUE;
    526         if (!ReadFile(GetEnvValue(env, prefix, "FileName", ""), psf))
    527             return kERROR;
    528     }
    529 
    530     return rc;
    531 }
     170double CalcRefractiveIndex(const double &lambda)
     171{
     172    // https://refractiveindex.info/?shelf=organic&book=poly(methyl_methacrylate)&page=Szczurowski
     173
     174    // n^2-1=\frac{0.99654λ^2}{λ^2-0.00787}+\frac{0.18964λ^2}{λ^2-0.02191}+\frac{0.00411λ^2}{λ^2-3.85727}
     175
     176    const double l2 = lambda*lambda;
     177
     178    const double c0 = 0.99654/(1-0.00787/l2);
     179    const double c1 = 0.18964/(1-0.02191/l2);
     180    const double c2 = 0.00411/(1-3.85727/l2);
     181
     182    return sqrt(1+c0+c1+c2);
     183}
     184
     185void ApplyRefraction(MQuaternion &u, const TVector3 &n, const double &n1, const double &n2)
     186{
     187    // u: incoming direction
     188    // n: normal vector of surface
     189    // n2: refractive index of new(?) medium
     190    // n1: refractive index of old(?) medium
     191
     192    const double r = n2/n1;
     193    const double c = n*u.fVectorPart;
     194
     195    const double R = 1 - r*r*(1-c*c);
     196
     197    const auto v = r*c + sqrt(R<0 ? 0 : R);
     198
     199    u = r*u - n*v;
     200}
     201
     202Int_t MLens::ExecuteOptics(MQuaternion &p, MQuaternion &u, const Short_t &wavelength) const
     203{
     204    const double F = 50;    // [mm] focal length
     205    const double R = 25;    // [mm] radius of lens
     206    const int    N = 25;    // [cnt] Nuber of grooves within radius
     207
     208    const double w = R/N;   // [mm] Width of a single groove
     209
     210    const double n = CalcRefractiveIndex(wavelength==0 ? 450 : wavelength);
     211
     212    const int nth = TMath::FloorNint(p.R()/w); // Ray will hit the nth groove
     213
     214    const double r0  =  nth*w;      // Lower radius of nth groove
     215    const double r1  = (nth+1)*w;   // Upper radius of nth groove
     216    const double cen = (nth+0.5)*w; // Center of nth groove
     217
     218    // FIXME: Do we have to check the nth-1 and nth+1 groove as well?
     219
     220    // Angle to the center of the groove
     221    const double phi = acos(cen / (F/1.125));
     222
     223    //        theta                peak_z
     224    const Groove g( TMath::Pi()/2 - phi, r0*tan(phi));
     225
     226    // Calculate the insersection between the ray and the cone
     227    float dz = CalcIntersection(p, u, g, -3);
     228
     229    // No groove was hit
     230    if (dz>=0)
     231        return -1;
     232
     233    // Propagate to the hit point
     234    p.PropagateZ(u, dz);
     235
     236    // Check where the ray has hit
     237    const double pr = p.R();
     238
     239    // If the ray has not hit within the right radius.. reject
     240    if (pr<=r0 || pr>r1)
     241        return -1;
     242
     243    // Normal vector of the surface at the hit point
     244    // FIXME: Surface roughness?
     245    TVector3 norm;
     246    norm.SetMagThetaPhi(1, g.fTheta+TMath::Pi()/2, p.XYvector().Phi());
     247
     248    // Apply refraction at lens entrance (change directional vector)
     249    // FIXME: Surace roughness
     250    ApplyRefraction(u, norm, 1, n);
     251
     252    // Propagate ray until bottom of lens
     253    const double v = u.fRealPart; // c changes in the medium
     254    u.fRealPart = n*v;
     255    p.PropagateZ(u, -3);
     256    u.fRealPart = v;  // Set back to c
     257
     258    // New normal surface (bottom of lens)
     259    norm.SetMagThetaPhi(1, 0, 0);
     260
     261    // Apply refraction at lens exit (change directional vector)
     262    // FIXME: Surace roughness
     263    ApplyRefraction(u, norm, n, 1);
     264
     265    // To make this consistent with a mirror system, we now change our coordinate system
     266    // Rays from the lens to the camera are up-going (positive sign)
     267    u.fVectorPart.SetZ(-u.Z());
     268
     269    // This is only correct if the focal distance is calculated from the inner lens surface!
     270    p.fVectorPart.SetZ(0);
     271
     272    return nth; // Keep track of groove index
     273}
  • trunk/Mars/msimreflector/MLens.h

    r19539 r19592  
    66#endif
    77
    8 #ifndef ROOT_TObjArray
    9 #include <TObjArray.h>
    10 #endif
    11 
    128class MQuaternion;
    13 class MMirror;
    149
    1510class MLens : public MOptics
    1611{
    1712private:
    18     //Simple Array               // 5.1s
    19     TObjArray fMirrors;          // 6.1s   (Pointer)
    20     //TObjArray fMirrors;        // 6.1s   (GetObjectRef)
    21     //TObjArray fMirrors;        // 8.3s   (Next)
    22     //TObjArray fMirrors;        // 10.1s  (UncheckedAt)
    23     //TList fMirrors;            // 10.7s
    24     //TOrdCollection fMirrors;   // 23.4s
    25 
    2613    Double_t fMaxR;
    2714
    2815    void InitMaxR();
    2916
    30     // Helper for I/O
    31     MMirror *EvalTokens(TObjArray &arr, Double_t defpsf) const;
    32 
    33     // MParContainer
    34     Int_t ReadEnv(const TEnv &env, TString prefix, Bool_t print);
    35 
    3617public:
    3718    MLens(const char *name=NULL, const char *title=NULL);
    38 
    39     const MMirror **GetFirstPtr() const;
    40     const UInt_t GetNumMirrors() const;
    41 
    42     const MMirror *GetMirror(UInt_t idx) const { return idx>=GetNumMirrors()?0:*(GetFirstPtr()+idx); }
    43 
    44     Bool_t ReadFile(TString fname, Double_t defpsf=-1);
    45     Bool_t WriteFile(TString fname) const;
    4619
    4720    Double_t GetMaxR() const { return fMaxR; }
     
    5023    virtual Bool_t CanHit(const MQuaternion &p) const;
    5124
    52     Int_t ExecuteOptics(MQuaternion &p, MQuaternion &u) const;
     25    Int_t ExecuteOptics(MQuaternion &p, MQuaternion &u, const Short_t &) const;
    5326
    54     void SetSigmaPSF(Double_t psf);
    55 
    56     virtual Bool_t IsValid() const { return fMirrors.GetEntries(); }
    57 
    58 
    59     // TObject
    60     void Paint(Option_t *o);
    61     void Print(Option_t *o) const;
     27    Bool_t IsValid() const { return kTRUE; }
    6228
    6329    ClassDef(MLens, 1) // Parameter container storing the description of a lens
Note: See TracChangeset for help on using the changeset viewer.