rl_raylib  1.1.10
rl_Polarization.cc
1 // File: rl_Polarization.cc
2 // Author: Terry Gaetz
3 
4 // --8<--8<--8<--8<--
5 //
6 // Copyright (C) 2006, 2007 Smithsonian Astrophysical Observatory
7 //
8 // This file is part of rl_raylib
9 //
10 // rl_raylib is free software; you can redistribute it and/or
11 // modify it under the terms of the GNU General Public License
12 // as published by the Free Software Foundation; either version 2
13 // of the License, or (at your option) any later version.
14 //
15 // rl_raylib is distributed in the hope that it will be useful,
16 // but WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with this program; if not, write to the
22 // Free Software Foundation, Inc.
23 // 51 Franklin Street, Fifth Floor
24 // Boston, MA 02110-1301, USA
25 //
26 // -->8-->8-->8-->8--
27 
28 #include <rl_Polarization.h> // rl_Polarization rl_PolCSPOD
29 
30 #include <cmath> // sqrt
31 #include <cfloat> // DBL_EPSILON
32 #include <mathconst/mathconst.h> // M_DEG2RAD
33 #include <iostream> // ostream
34 
35 #include <rl_raylib/rl_ReflectionCoefPOD.h> // rl_ReflectionCoefPOD
36 
37 using namespace std;
38 
39 /* mean graze angles for shells */
40 #ifdef CONST_SINGRAZE
41  static double alpha[] = {51.2, 46.2, 41.2, 36.4, 31.7, 27.0};
42 #endif
43 
44 #ifdef POLREFL_DEBUG
45  static FILE *dbgw; /* --- tjg debug --- */
46 #endif /* POLREFL_DEBUG */
47 
48 //#########################################################################
49 // forward declaration of statics...
50 //#########################################################################
51 
52 //-------------------------------------------------------------------------
53 static void
54 make_Triad( dvm3_Vector& x,
55  dvm3_Vector& y,
56  dvm3_Vector const& s );
57 
58 //-------------------------------------------------------------------------
59 static void
60 pol_perp_par(
61  dvm3_Vector& perp, // out: perpendicular direction vector
62  dvm3_Vector& parv1, // out: parallel direction vector ("v1")
63  dvm3_Vector& parv2, // out: parallel direction vector ("v2")
64  dvm3_Vector const& norm, // in: normal vector
65  dvm3_Vector const& v1, // in: initial ray direction vector
66  dvm3_Vector const& v2 // in: final ray direction vector
67 );
68 
69 
70 //#########################################################################
71 // rl_PolCSPOD methods
72 //#########################################################################
73 
74 //-------------------------------------------------------------------------
75 void rl_PolCSPOD::
76 init( double b_over_a, // polarization ellipse minor/major
77  double psi // polarization ellipse angle
78  )
79 {
80  //--- initialize ray polarization ----
81 
82  double cpsi = cos(psi * M_DEG2RAD);
83  double spsi = sin(psi * M_DEG2RAD);
84 
85  double sqre = sqrt(1.0 + b_over_a * b_over_a );
86 
87  c2_[0] = rl_Traits::complex( cpsi / sqre, (double)0.0 );
88  c2_[1] = rl_Traits::complex( spsi / sqre, (double)0.0 );
89 
90  rl_Traits::complex ctmp( -spsi / sqre, (double)0.0 );
91  s2_[0] = b_over_a * ctmp;
92 
93  ctmp = rl_Traits::complex( cpsi / sqre, (double)0.0 );
94  s2_[1] = b_over_a * ctmp;
95 }
96 
97 //-------------------------------------------------------------------------
98 double rl_PolCSPOD::
99 intensity() const
100 {
101  rl_Traits::complex crossterms( c2_[0] * conj( s2_[0] )
102  + c2_[1] * conj( s2_[1] ) );
103  double rfl = norm(c2_[0]) + norm(s2_[0]) + norm(c2_[1]) + norm(s2_[1])
104  + 2.0 * crossterms.real();
105 
106  #ifdef OSAC_WEIGHT_KLUDGE
107  /*
108  * for exact agreement with OSAC - degrade weight to 8 places.
109  *
110  * 'Make w agree with its integer-equivalent output value'
111  * (see vflect.f)
112  */
113  return 1.e-8 * ((int)(1.0e8 * rfl + 0.5));
114  #else
115  return rfl;
116  #endif
117 }
118 
119 //-------------------------------------------------------------------------
120 void rl_PolCSPOD::
121 attenuate( double factor )
122 {
123  const double fac = sqrt(factor);
124  c2_[0] *= fac;
125  c2_[1] *= fac;
126  s2_[0] *= fac;
127  s2_[1] *= fac;
128 }
129 
130 //-------------------------------------------------------------------------
131 // write polarization amplitude to stream os. The state is
132 // written as "[(r,i)(r,i)][(r,i)(r,i)]".
133 
134 std::ostream& rl_PolCSPOD::
135 print_on( std::ostream& os ) const
136 {
137  os << "[(" << c2_[0].real() << "," << c2_[0].imag() << ")("
138  << c2_[1].real() << "," << c2_[1].imag() << ")][("
139  << s2_[0].real() << "," << s2_[0].imag() << ")("
140  << s2_[1].real() << "," << s2_[1].imag() << ")]";
141  return os;
142 }
143 
144 //-------------------------------------------------------------------------
145 // write polarization amplitude to stream os. The state is
146 // written as "[(r,i)(r,i)]\n[(r,i)(r,i)]\nintensity".
147 
148 void rl_PolCSPOD::
149 cprint_on( FILE* of ) const
150 {
151  fprintf(of, "\nc2: [(%.15e, %.15e)(%.15e, %.15e)]",
152  c2_[0].real(), c2_[0].imag(),
153  c2_[1].real(), c2_[1].imag() );
154  fprintf(of, "\ns2: [(%.15e, %.15e)(%.15e, %.15e)]",
155  s2_[0].real(), s2_[0].imag(),
156  s2_[1].real(), s2_[1].imag() );
157  fprintf(of, "\n%.15e\n", intensity());
158 }
159 
160 
161 
162 //#########################################################################
163 // rl_Polarization methods
164 //#########################################################################
165 
166 //=========================================================================
167 // ctors/dtors...
168 
169 //-------------------------------------------------------------------------
171 init( rl_PolCSPOD const& cs, dvm3_Vector const& dir )
172 {
173  dvm3_Vector x1;
174  dvm3_Vector y1;
175 
176  // define x1 and y1 unit vectors...
177  //
178  // x1 is the normalized y1 <cross-product> cs,
179  // y1 is v <cross-product> x1
180  //
181  // the set [x1, y1, dir] forms an orthonormal basis.
182 
183  make_Triad( x1, y1, dir );
184 
185  // Construct real & imag parts of C and S vectors from
186  // c2_.q1, c2_.q2, s1_.q1, and s2_.q2 amplitudes.
187  //
188  // Note that c2_ and s2_ are complex coefficients
189  // multiplying the (real) polarization unit vectors.
190  //
191  // C_r_ = Re[ c2_.q1 * x1 + c2_.q2 * y1 ]
192  // C_i_ = Im[ c2_.q1 * x1 + c2_.q2 * y1 ]
193  //
194  // S_r_ = Re[ s2_.q1 * x1 + s2_.q2 * y1 ]
195  // S_i_ = Im[ s2_.q1 * x1 + s2_.q2 * y1 ]
196 
197  C_r_.lincomb( cs.c2_[0].real(), x1, cs.c2_[1].real(), y1 );
198  C_i_.lincomb( cs.c2_[0].imag(), x1, cs.c2_[1].imag(), y1 );
199  S_r_.lincomb( cs.s2_[0].real(), x1, cs.s2_[1].real(), y1 );
200  S_i_.lincomb( cs.s2_[0].imag(), x1, cs.s2_[1].imag(), y1 );
201 }
202 
203 //=========================================================================
204 // accessors...
205 
206 //-------------------------------------------------------------------------
207 double rl_Polarization::
208 intensity() const
209 {
210  return dot(C_r_, C_r_) + dot(C_i_, C_i_) + dot(S_r_, S_r_) + dot(S_i_, S_i_)
211  + 2.0 * ( dot(C_i_, S_r_) - dot(C_r_, S_i_) );
212 }
213 
214 //-------------------------------------------------------------------------
215 
217 get_PolCSPOD( rl_PolCSPOD& cs, dvm3_Vector const& dir ) const
218 {
219  // define x2 and y2 unit vectors...
220  // x2 is the normalized y2 <cross-product> cs,
221  // y2 is v <cross-product> x2
222 
223  dvm3_Vector x2;
224  dvm3_Vector y2;
225  make_Triad( x2, y2, dir );
226 
227  // Convert the C, S complex vectors back into c2_ and s2_
228  // complex amplitudes
229 
230  double d1 = dot( x2, C_r_);
231  double d2 = dot( x2, C_i_);
232  cs.c2_[0] = rl_Traits::complex( d1, d2 );
233  cs.c2_[1] = rl_Traits::complex( dot( y2, C_r_), dot( y2, C_i_) );
234  cs.s2_[0] = rl_Traits::complex( dot( x2, S_r_), dot( x2, S_i_) );
235  cs.s2_[1] = rl_Traits::complex( dot( y2, S_r_), dot( y2, S_i_) );
236 }
237 
238 //=========================================================================
239 // mutators...
240 
241 //-------------------------------------------------------------------------
243 attenuate( double factor )
244 {
245  const double fac = sqrt(factor);
246  C_r_ *= fac;
247  C_i_ *= fac;
248  S_r_ *= fac;
249  S_i_ *= fac;
250 }
251 
252 //-------------------------------------------------------------------------
253 // Description: Reflect the polarization state of a ray:
254 //
255 // Given surface normal, and incident and reflected ray direction vectors,
256 // and the parallel and perpendicular reflection coefficients,
257 // evaluate the reflected polarization vectors.
258 //
259 // Notes : based in part on OSAC routine "pzrefl"
260 
263  dvm3_Vector const& normal,
264  dvm3_Vector const& dir_in,
265  dvm3_Vector const& dir_out,
266  rl_ReflectionCoefPOD const& rflcoef
267  )
268 {
269  // calculate sine of the grazing angle
270  // double sg = dot( dir_in, normal );
271 
272  #ifdef CONST_SINGRAZE
273  sg = sin(alpha[CONST_SINGRAZE-1] * 0.0002908882);
274  #endif
275 
276  //x // reflect the direction vector:
277  //x // v_o = (1 - 2 norm <diadprod> norm) <dotprod> v_i
278  //x
279  //x dir_out.lincomb( 1.0, dir_in, -2.0 * sg, &srf->norm );
280 
281  // form the incoming/outgoing polarization basis vectors
282 
283  dvm3_Vector perp; // perpendicular to vo-vn plane
284  dvm3_Vector parvi; // normal to vi, in vi-vn plane
285  dvm3_Vector parvo; // normal to vo, in vo-vn plane
286 
287  pol_perp_par( perp, parvi, parvo, normal,
288  dir_in, dir_out );
289 
290  // build components of T matrix:
291  //
292  // T_par = parvo <dyadproduct> parvi
293  // T_prp = prpo <dyadproduct> prpi
294  //
295  // Note that (for now) prpo = prpi = perp
296 
297  dvm3_Matrix T_prp; // projector for perp vector
298  dvm3_Matrix T_par; // projector for par vector
299 
300  T_par.dyad_product( parvo, parvi );
301  T_prp.dyad_product( perp, perp );
302 
303  // combine with complex reflectivities; form real and imag parts of T
304 
305  dvm3_Matrix T_re; // real part of T matrix
306  dvm3_Matrix T_im; // imag part of T matrix
307 
308  T_re.lincomb( rflcoef.para().real(), T_par,
309  rflcoef.perp().real(), T_prp );
310  T_im.lincomb( rflcoef.para().imag(), T_par,
311  rflcoef.perp().imag(), T_prp );
312 
313  /*
314  * Apply Re, Im part of T matrix to Re, Im polarization vectors.
315  */
316  dvm3_Vector Tr_Cr; T_re.mvmult( Tr_Cr, C_r_ ); // T_re . C_r_
317  dvm3_Vector Tr_Ci; T_re.mvmult( Tr_Ci, C_i_ ); // T_re . C_i_
318  dvm3_Vector Ti_Cr; T_im.mvmult( Ti_Cr, C_r_ ); // T_im . C_r_
319  dvm3_Vector Ti_Ci; T_im.mvmult( Ti_Ci, C_i_ ); // T_im . C_i_
320 
321  dvm3_Vector Tr_Sr; T_re.mvmult( Tr_Sr, S_r_ ); // T_re . S_r_
322  dvm3_Vector Tr_Si; T_re.mvmult( Tr_Si, S_i_ ); // T_re . S_i_
323  dvm3_Vector Ti_Sr; T_im.mvmult( Ti_Sr, S_r_ ); // T_im . S_r_
324  dvm3_Vector Ti_Si; T_im.mvmult( Ti_Si, S_i_ ); // T_im . S_i_
325 
326  // collect and reassemble into outgoing PolVec...
327 
328  C_r_ = Tr_Cr - Ti_Ci;
329  C_i_ = Tr_Ci + Ti_Cr;
330  S_r_ = Tr_Sr - Ti_Si;
331  S_i_ = Tr_Si + Ti_Sr;
332 }
333 
334 //-------------------------------------------------------------------------
335 // write polarization vectors to stream os. The state is
336 // written as "[x y z][x y z][x y z][x y z]".
337 
338 std::ostream& rl_Polarization::
339 print_on( std::ostream& os ) const
340 {
341  os << C_r_ << C_i_ << S_r_ << S_i_;
342  return os;
343 }
344 
345 //-------------------------------------------------------------------------
346 // write polarization vectors to stream os. The state is
347 // written as "[x y z][x y z][x y z][x y z]".
348 
350 cprint_on( FILE* of ) const
351 {
352  double* vec = new double[3];
353  C_r_.copy_to_cvec( vec );
354  fprintf(of, "\nC_r_: [%.15e, %.15e %.15e]\n",
355  vec[0], vec[1], vec[2]);
356  C_i_.copy_to_cvec( vec );
357  fprintf(of, "\nC_i_: [%.15e, %.15e %.15e]\n",
358  vec[0], vec[1], vec[2]);
359  S_r_.copy_to_cvec( vec );
360  fprintf(of, "\nS_r_: [%.15e, %.15e %.15e]\n",
361  vec[0], vec[1], vec[2]);
362  S_i_.copy_to_cvec( vec );
363  fprintf(of, "\nS_i_: [%.15e, %.15e %.15e]\n",
364  vec[0], vec[1], vec[2]);
365  delete [] vec;
366 }
367 
368 #define TINY_VALUE (100.0 * DBL_EPSILON)
369 
370 //#########################################################################
371 // file static definitions...
372 //#########################################################################
373 
374 //-------------------------------------------------------------------------
375 void
376 pol_perp_par(
377  dvm3_Vector& perp, // out: perpendicular direction vector
378  dvm3_Vector& parv1, // out: parallel direction vector ("v1")
379  dvm3_Vector& parv2, // out: parallel direction vector ("v2")
380  dvm3_Vector const& norm, // in: normal vector
381  dvm3_Vector const& v1, // in: initial ray direction vector
382  dvm3_Vector const& v2 // in: final ray direction vector
383 )
384 {
385  // find perpendicular unit vector
386 
387  if ( fabs(dot(v1, v2)) < 0.99999 )
388  {
389  perp.cross( v2, v1 ); // use v2 cross v1
390  }
391  else
392  {
393  if ( fabs(dot(norm,v1)) < 0.99999 )
394  {
395  perp.cross( norm, v1 ); // use norm cross v1
396  }
397  else
398  {
399  if ( fabs(1.0 - fabs( norm[0] )) < TINY_VALUE )
400  {
401  perp.init( norm[2], 0.0, -norm[0] ); // use norm cross e_y
402  }
403  else
404  {
405  perp.init( 0.0, -norm[2], -norm[1] ); // use norm cross e_x
406  }
407  }
408  }
409  perp.unitize();
410 
411  // now find the parallel unit vectors
412 
413  parv1.cross( v1, perp );
414  parv2.cross( v2, perp );
415 }
416 
417 //-------------------------------------------------------------------------
418 void
419 make_Triad( dvm3_Vector& x,
420  dvm3_Vector& y,
421  dvm3_Vector const& s )
422 {
423  // form the "x" and "y" unit vectors
424 
425  // make "x" unit vector
426 
427  if ( fabs(1.0 - fabs( s[1] )) > TINY_VALUE )
428  {
429  x.init( s[2], 0.0, -s[0] ); // x = e_x <crossprod> s
430  }
431  else
432  {
433  x.init( 1.0, 0.0, 0.0 ); // x = e_x is close enough
434  }
435  x.unitize(); // make x a unit vector
436 
437  // y = s <crossprod> x; s and x are orthogonal
438  // unit vectors so normalization isn't needed.
439 
440  y.cross( s, x ); // y = s <crossprod> x
441 
442  // we now have an orthonormal triad.
443 }
444 #undef TINY_VALUE
445 
446 #undef CONST_SINGRAZE
std::complex< double > complex
Typedef for the complex type.
Definition: rl_Traits.h:61
void cprint_on(std::FILE *of) const
Write the polarization amplitude to FILE* of.
rl_Traits::complex c2_[2]
polarization amplitude (cosine component)
void init(double b_over_a=1.0, double psi=0.0)
Initialization method.
void cprint_on(FILE *of) const
Write the polarization amplitude to FILE* of.
double intensity() const
Evaluate the intensity.
std::ostream & print_on(std::ostream &os) const
Write the polarization vectors to stream os.
std::ostream & print_on(std::ostream &os) const
Write the polarization amplitude to stream os.
rl_Traits::complex s2_[2]
polarization amplitude (sine component)
A Plain Ol' Data class representing complex reflection coefficients.
A Plain Ol' Data struct (POD) encapsulating the OSAC-style complex polarization amplitudes.
void attenuate(double factor)
Attenuate intensity by a factor.
void attenuate(double factor)
Attenuate reflectivity by a factor.
void init(rl_PolCSPOD const &cs, dvm3_Vector const &dir)
Initialize the polarization state.
void get_PolCSPOD(rl_PolCSPOD &cs, dvm3_Vector const &dir) const
Evaluate rl_PolCSPOD corresponding to this polarization state.
double intensity() const
Evaluate intensity.
void reflect(dvm3_Vector const &normal, dvm3_Vector const &dir_in, dvm3_Vector const &dir_out, rl_ReflectionCoefPOD const &rflcoef)
Evaluate the reflected polarization vectors.
rl_Traits::complex para() const
rl_Traits::complex perp() const