rl_Polarization.cc

00001 // File:   rl_Polarization.cc
00002 // Author: Terry Gaetz
00003 
00004 // --8<--8<--8<--8<--
00005 //
00006 // Copyright (C) 2006, 2007 Smithsonian Astrophysical Observatory
00007 //
00008 // This file is part of rl_raylib
00009 //
00010 // rl_raylib is free software; you can redistribute it and/or
00011 // modify it under the terms of the GNU General Public License
00012 // as published by the Free Software Foundation; either version 2
00013 // of the License, or (at your option) any later version.
00014 //
00015 // rl_raylib is distributed in the hope that it will be useful,
00016 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00018 // GNU General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU General Public License
00021 // along with this program; if not, write to the 
00022 //       Free Software Foundation, Inc. 
00023 //       51 Franklin Street, Fifth Floor
00024 //       Boston, MA  02110-1301, USA
00025 //
00026 // -->8-->8-->8-->8--
00027 
00028 #include <rl_Polarization.h>                 // rl_Polarization rl_PolCSPOD
00029 
00030 #include <cmath>                             // sqrt
00031 #include <cfloat>                            // DBL_EPSILON
00032 #include <mathconst/mathconst.h>             // M_DEG2RAD
00033 #include <iostream>                          // ostream
00034 
00035 #include <rl_raylib/rl_ReflectionCoefPOD.h>  // rl_ReflectionCoefPOD
00036 
00037 using namespace std;
00038 
00039 /* mean graze angles for shells */
00040 #ifdef CONST_SINGRAZE
00041   static double alpha[] = {51.2, 46.2, 41.2, 36.4, 31.7, 27.0};
00042 #endif
00043 
00044 #ifdef    POLREFL_DEBUG
00045   static FILE *dbgw; /* --- tjg debug --- */
00046 #endif /* POLREFL_DEBUG */
00047 
00048 //#########################################################################
00049 // forward declaration of statics...
00050 //#########################################################################
00051 
00052 //-------------------------------------------------------------------------
00053 static void
00054 make_Triad( dvm3_Vector&       x,
00055             dvm3_Vector&       y,
00056             dvm3_Vector const& s );
00057 
00058 //-------------------------------------------------------------------------
00059 static void
00060 pol_perp_par(
00061         dvm3_Vector&       perp,  // out: perpendicular direction vector
00062         dvm3_Vector&       parv1, // out: parallel direction vector ("v1")
00063         dvm3_Vector&       parv2, // out: parallel direction vector ("v2")
00064         dvm3_Vector const& norm,  // in:  normal vector
00065         dvm3_Vector const& v1,    // in:  initial  ray direction vector
00066         dvm3_Vector const& v2     // in:  final    ray direction vector
00067 );
00068 
00069 
00070 //#########################################################################
00071 // rl_PolCSPOD methods
00072 //#########################################################################
00073 
00074 //-------------------------------------------------------------------------
00075 void rl_PolCSPOD::
00076 init( double  b_over_a,  // polarization ellipse minor/major
00077       double  psi        // polarization ellipse angle
00078     )
00079 {
00080   //--- initialize ray polarization ----
00081 
00082   double cpsi = cos(psi * M_DEG2RAD);
00083   double spsi = sin(psi * M_DEG2RAD);
00084 
00085   double sqre = sqrt(1.0 + b_over_a * b_over_a );
00086 
00087   c2_[0] =  rl_Traits::complex( cpsi / sqre, (double)0.0 );
00088   c2_[1] =  rl_Traits::complex( spsi / sqre, (double)0.0 );
00089 
00090   rl_Traits::complex ctmp( -spsi / sqre, (double)0.0 );
00091   s2_[0] =  b_over_a * ctmp;
00092 
00093   ctmp   =  rl_Traits::complex( cpsi / sqre, (double)0.0 );
00094   s2_[1] =  b_over_a * ctmp;
00095 }
00096 
00097 //-------------------------------------------------------------------------
00098 double rl_PolCSPOD::
00099 intensity() const
00100 {
00101   rl_Traits::complex crossterms( c2_[0] * conj( s2_[0] )
00102                                + c2_[1] * conj( s2_[1] ) );
00103   double rfl = norm(c2_[0]) + norm(s2_[0]) + norm(c2_[1]) + norm(s2_[1])
00104              + 2.0 * crossterms.real();
00105 
00106   #ifdef OSAC_WEIGHT_KLUDGE
00107     /*
00108      * for exact agreement with OSAC - degrade weight to 8 places.
00109      *
00110      * 'Make w agree with its integer-equivalent output value'
00111      *    (see vflect.f)
00112      */
00113     return 1.e-8 * ((int)(1.0e8 * rfl + 0.5));
00114   #else
00115      return rfl;
00116   #endif
00117 }
00118 
00119 //-------------------------------------------------------------------------
00120 void rl_PolCSPOD::
00121 attenuate( double factor )
00122 {
00123   const double fac = sqrt(factor);
00124   c2_[0] *= fac;
00125   c2_[1] *= fac;
00126   s2_[0] *= fac;
00127   s2_[1] *= fac;
00128 }
00129 
00130 //-------------------------------------------------------------------------
00131 // write polarization amplitude to stream os.  The state is
00132 // written as "[(r,i)(r,i)][(r,i)(r,i)]".
00133 
00134 std::ostream& rl_PolCSPOD::
00135 print_on( std::ostream& os ) const
00136 {
00137   os << "[(" << c2_[0].real() << "," << c2_[0].imag() << ")("
00138              << c2_[1].real() << "," << c2_[1].imag() << ")][("
00139              << s2_[0].real() << "," << s2_[0].imag() << ")("
00140              << s2_[1].real() << "," << s2_[1].imag() << ")]";
00141   return os;
00142 }
00143 
00144 //-------------------------------------------------------------------------
00145 // write polarization amplitude to stream os.  The state is
00146 // written as "[(r,i)(r,i)]\n[(r,i)(r,i)]\nintensity".
00147 
00148 void rl_PolCSPOD::
00149 cprint_on( FILE* of ) const
00150 {
00151   fprintf(of, "\nc2: [(%.15e, %.15e)(%.15e, %.15e)]",
00152               c2_[0].real(), c2_[0].imag(),
00153               c2_[1].real(), c2_[1].imag() );
00154   fprintf(of, "\ns2: [(%.15e, %.15e)(%.15e, %.15e)]",
00155               s2_[0].real(), s2_[0].imag(),
00156               s2_[1].real(), s2_[1].imag() );
00157   fprintf(of, "\n%.15e\n", intensity());
00158 }
00159 
00160 
00161 
00162 //#########################################################################
00163 // rl_Polarization methods
00164 //#########################################################################
00165 
00166 //=========================================================================
00167 // ctors/dtors...
00168 
00169 //-------------------------------------------------------------------------
00170 void rl_Polarization::
00171 init( rl_PolCSPOD const& cs, dvm3_Vector const& dir )
00172 {
00173   dvm3_Vector x1;
00174   dvm3_Vector y1;
00175 
00176   // define x1 and y1 unit vectors...
00177   //
00178   // x1 is the normalized y1 <cross-product> cs,
00179   // y1 is v <cross-product> x1
00180   //
00181   // the set [x1, y1, dir] forms an orthonormal basis.
00182 
00183   make_Triad( x1, y1, dir );
00184 
00185   // Construct real & imag parts of C and S vectors from
00186   // c2_.q1, c2_.q2, s1_.q1, and s2_.q2 amplitudes.
00187   //
00188   // Note that c2_ and s2_ are complex coefficients
00189   // multiplying the (real) polarization unit vectors.
00190   //
00191   // C_r_ = Re[ c2_.q1 * x1 + c2_.q2 * y1 ]
00192   // C_i_ = Im[ c2_.q1 * x1 + c2_.q2 * y1 ]
00193   //
00194   // S_r_ = Re[ s2_.q1 * x1 + s2_.q2 * y1 ]
00195   // S_i_ = Im[ s2_.q1 * x1 + s2_.q2 * y1 ]
00196 
00197   C_r_.lincomb( cs.c2_[0].real(),  x1,  cs.c2_[1].real(),  y1 );
00198   C_i_.lincomb( cs.c2_[0].imag(),  x1,  cs.c2_[1].imag(),  y1 );
00199   S_r_.lincomb( cs.s2_[0].real(),  x1,  cs.s2_[1].real(),  y1 );
00200   S_i_.lincomb( cs.s2_[0].imag(),  x1,  cs.s2_[1].imag(),  y1 );
00201 }
00202 
00203 //=========================================================================
00204 // accessors...
00205 
00206 //-------------------------------------------------------------------------
00207 double rl_Polarization::
00208 intensity() const
00209 {
00210   return dot(C_r_, C_r_) + dot(C_i_, C_i_) + dot(S_r_, S_r_) + dot(S_i_, S_i_)
00211                  + 2.0 * ( dot(C_i_, S_r_) - dot(C_r_, S_i_) );
00212 }
00213 
00214 //-------------------------------------------------------------------------
00215 
00216 void rl_Polarization::
00217 get_PolCSPOD( rl_PolCSPOD& cs, dvm3_Vector const& dir ) const
00218 {
00219   // define x2 and y2 unit vectors...
00220   // x2 is the normalized y2 <cross-product> cs,
00221   // y2 is v <cross-product> x2
00222 
00223   dvm3_Vector x2;
00224   dvm3_Vector y2;
00225   make_Triad( x2, y2, dir );
00226 
00227   // Convert the C, S complex vectors back into c2_ and s2_
00228   // complex amplitudes
00229 
00230   double d1 = dot( x2, C_r_);
00231   double d2 = dot( x2, C_i_);
00232   cs.c2_[0] = rl_Traits::complex( d1, d2 );
00233   cs.c2_[1] = rl_Traits::complex( dot( y2, C_r_), dot( y2, C_i_) );
00234   cs.s2_[0] = rl_Traits::complex( dot( x2, S_r_), dot( x2, S_i_) );
00235   cs.s2_[1] = rl_Traits::complex( dot( y2, S_r_), dot( y2, S_i_) );
00236 }
00237 
00238 //=========================================================================
00239 // mutators...
00240 
00241 //-------------------------------------------------------------------------
00242 void rl_Polarization::
00243 attenuate( double factor )
00244 {
00245   const double fac = sqrt(factor);
00246   C_r_ *= fac;
00247   C_i_ *= fac;
00248   S_r_ *= fac;
00249   S_i_ *= fac;
00250 }
00251 
00252 //-------------------------------------------------------------------------
00253 // Description: Reflect the polarization state of a ray:
00254 //
00255 // Given surface normal, and incident and reflected ray direction vectors,
00256 //  and the parallel and perpendicular reflection coefficients,
00257 //  evaluate the reflected polarization vectors.
00258 //
00259 // Notes      : based in part on OSAC routine "pzrefl"
00260 
00261 void rl_Polarization::
00262 reflect(
00263          dvm3_Vector          const& normal,
00264          dvm3_Vector          const& dir_in,
00265          dvm3_Vector          const& dir_out,
00266          rl_ReflectionCoefPOD const& rflcoef
00267        )
00268 {
00269   // calculate sine of the grazing angle
00270   // double sg = dot( dir_in, normal );
00271 
00272   #ifdef CONST_SINGRAZE
00273     sg = sin(alpha[CONST_SINGRAZE-1] * 0.0002908882);
00274   #endif
00275 
00276   //x // reflect the direction vector:
00277   //x //    v_o = (1 - 2 norm <diadprod> norm) <dotprod> v_i
00278   //x 
00279   //x dir_out.lincomb( 1.0, dir_in, -2.0 * sg, &srf->norm );
00280 
00281   // form the incoming/outgoing polarization basis vectors 
00282 
00283   dvm3_Vector  perp;        // perpendicular to vo-vn plane
00284   dvm3_Vector  parvi;       // normal to vi, in vi-vn plane
00285   dvm3_Vector  parvo;       // normal to vo, in vo-vn plane
00286 
00287   pol_perp_par( perp, parvi, parvo, normal, 
00288                 dir_in, dir_out );
00289   
00290   // build components of T matrix:
00291   //
00292   //   T_par = parvo <dyadproduct> parvi
00293   //   T_prp = prpo  <dyadproduct> prpi
00294   //
00295   //   Note that (for now) prpo = prpi = perp
00296 
00297   dvm3_Matrix  T_prp;       // projector for perp vector
00298   dvm3_Matrix  T_par;       // projector for par  vector
00299 
00300   T_par.dyad_product( parvo, parvi );
00301   T_prp.dyad_product( perp,  perp  );
00302 
00303   // combine with complex reflectivities; form real and imag parts of T 
00304 
00305   dvm3_Matrix  T_re;        // real part of T matrix
00306   dvm3_Matrix  T_im;        // imag part of T matrix
00307 
00308   T_re.lincomb( rflcoef.para().real(), T_par, 
00309                 rflcoef.perp().real(), T_prp );
00310   T_im.lincomb( rflcoef.para().imag(), T_par, 
00311                 rflcoef.perp().imag(), T_prp );
00312 
00313   /* 
00314    * Apply Re, Im part of T matrix to Re, Im polarization vectors.
00315    */ 
00316   dvm3_Vector  Tr_Cr;    T_re.mvmult( Tr_Cr, C_r_ );              // T_re . C_r_
00317   dvm3_Vector  Tr_Ci;    T_re.mvmult( Tr_Ci, C_i_ );              // T_re . C_i_
00318   dvm3_Vector  Ti_Cr;    T_im.mvmult( Ti_Cr, C_r_ );              // T_im . C_r_
00319   dvm3_Vector  Ti_Ci;    T_im.mvmult( Ti_Ci, C_i_ );              // T_im . C_i_
00320 
00321   dvm3_Vector  Tr_Sr;    T_re.mvmult( Tr_Sr, S_r_ );              // T_re . S_r_
00322   dvm3_Vector  Tr_Si;    T_re.mvmult( Tr_Si, S_i_ );              // T_re . S_i_
00323   dvm3_Vector  Ti_Sr;    T_im.mvmult( Ti_Sr, S_r_ );              // T_im . S_r_
00324   dvm3_Vector  Ti_Si;    T_im.mvmult( Ti_Si, S_i_ );              // T_im . S_i_
00325 
00326   // collect and reassemble into outgoing PolVec... 
00327 
00328   C_r_ = Tr_Cr - Ti_Ci;
00329   C_i_ = Tr_Ci + Ti_Cr;
00330   S_r_ = Tr_Sr - Ti_Si;
00331   S_i_ = Tr_Si + Ti_Sr;
00332 }
00333 
00334 //-------------------------------------------------------------------------
00335 // write polarization vectors to stream os.  The state is
00336 // written as "[x y z][x y z][x y z][x y z]".
00337 
00338 std::ostream& rl_Polarization::
00339 print_on( std::ostream& os ) const
00340 {
00341   os << C_r_ << C_i_ << S_r_ << S_i_;
00342   return os;
00343 }
00344 
00345 //-------------------------------------------------------------------------
00346 // write polarization vectors to stream os.  The state is
00347 // written as "[x y z][x y z][x y z][x y z]".
00348 
00349 void rl_Polarization::
00350 cprint_on( FILE* of ) const
00351 {
00352   double* vec = new double[3];
00353   C_r_.copy_to_cvec( vec );
00354   fprintf(of, "\nC_r_: [%.15e, %.15e %.15e]\n",
00355               vec[0], vec[1], vec[2]);
00356   C_i_.copy_to_cvec( vec );
00357   fprintf(of, "\nC_i_: [%.15e, %.15e %.15e]\n",
00358               vec[0], vec[1], vec[2]);
00359   S_r_.copy_to_cvec( vec );
00360   fprintf(of, "\nS_r_: [%.15e, %.15e %.15e]\n",
00361               vec[0], vec[1], vec[2]);
00362   S_i_.copy_to_cvec( vec );
00363   fprintf(of, "\nS_i_: [%.15e, %.15e %.15e]\n",
00364               vec[0], vec[1], vec[2]);
00365   delete [] vec;
00366 }
00367 
00368 #define TINY_VALUE (100.0 * DBL_EPSILON)
00369 
00370 //#########################################################################
00371 // file static definitions...
00372 //#########################################################################
00373 
00374 //-------------------------------------------------------------------------
00375 void 
00376 pol_perp_par(
00377         dvm3_Vector&       perp,  // out: perpendicular direction vector
00378         dvm3_Vector&       parv1, // out: parallel direction vector ("v1")
00379         dvm3_Vector&       parv2, // out: parallel direction vector ("v2")
00380         dvm3_Vector const& norm,  // in:  normal vector
00381         dvm3_Vector const& v1,    // in:  initial  ray direction vector
00382         dvm3_Vector const& v2     // in:  final    ray direction vector
00383 )
00384 {
00385   // find perpendicular unit vector
00386 
00387   if ( fabs(dot(v1, v2)) < 0.99999 )
00388   {
00389     perp.cross( v2, v1 );       // use v2 cross v1
00390   }
00391   else
00392   {
00393     if ( fabs(dot(norm,v1)) < 0.99999 )
00394     {
00395       perp.cross( norm, v1 );   // use norm cross v1
00396     }
00397     else
00398     {
00399       if ( fabs(1.0 - fabs( norm[0] )) < TINY_VALUE )
00400       {
00401         perp.init( norm[2], 0.0, -norm[0] );    // use norm cross e_y
00402       }
00403       else
00404       {
00405         perp.init( 0.0, -norm[2], -norm[1] );   // use norm cross e_x
00406       }
00407     }
00408   }
00409   perp.unitize();
00410 
00411   // now find the parallel unit vectors
00412 
00413   parv1.cross( v1, perp );
00414   parv2.cross( v2, perp );
00415 }
00416 
00417 //-------------------------------------------------------------------------
00418 void
00419 make_Triad( dvm3_Vector&       x,
00420             dvm3_Vector&       y,
00421             dvm3_Vector const& s )
00422 {
00423   // form the "x" and "y" unit vectors
00424 
00425   // make "x" unit vector
00426 
00427   if ( fabs(1.0 - fabs( s[1] )) > TINY_VALUE )
00428   {
00429     x.init( s[2], 0.0, -s[0] );       // x = e_x <crossprod> s
00430   }
00431   else
00432   {
00433     x.init( 1.0, 0.0, 0.0 );          // x = e_x is close enough
00434   }
00435   x.unitize();                        // make x a unit vector
00436 
00437   // y = s <crossprod> x; s and x are orthogonal
00438   // unit vectors so normalization isn't needed.
00439 
00440   y.cross( s, x );                    // y = s <crossprod> x
00441 
00442   // we now have an orthonormal triad.
00443 }
00444 #undef TINY_VALUE
00445 
00446 #undef CONST_SINGRAZE

Generated on Mon Nov 3 18:15:05 2008 for rl_raylib by  doxygen 1.5.6