NASA Logo
Ocean Color Science Software

ocssw V2022
interpolate_libraries.f90
Go to the documentation of this file.
2 
3  implicit none
4 
5  private
6 
7  public :: libraryinterpolate, scatangle, interpolate_wind_speed, lib_clean, lib_init ! GetReflForGivenWindSpeedIce, GetReflForGivenWindSpeedWater
8 
9 ! REAL*8 :: CONST_PI, RAD_PER_DEG , DEG_PER_RAD
10  real, parameter, dimension(3) :: wind_speed = (/3.0, 7.0, 15.0/)
11  integer :: NUM_OF_TAUS, NUM_OF_WAVELENGTHS_RFTF
12 
13  real :: RLphase
14  real, dimension(:), allocatable :: aeroPhase
15  real, dimension(:,:), allocatable :: fOmega_water, fOmega_ice
16  real, dimension(:,:,:), allocatable :: tauC_source_water, tauC_source_ice
17 
18  REAL,ALLOCATABLE::multiScatRefl_water_lamb(:,:,:)
19  REAL,ALLOCATABLE::multiScatRefl_ice_lamb(:,:,:)
20 
21  REAL,ALLOCATABLE::multiScatRefl_water_CM(:,:,:,:)
22  REAL,ALLOCATABLE::multiScatRefl_ice_CM(:,:,:,:)
23 
24 ! single-scattering arrays
25  real, allocatable :: fPhase_lamb(:,:), tauC(:,:)
26  REAL, ALLOCATABLE :: gsquare(:), aero_denom(:)
27 
28  integer :: lib_call_counter
29 
30  REAL,ALLOCATABLE::phaseFunVals_water(:,:), phaseFunVals_ice(:,:), ssRefl_water(:,:,:), ssRefl_ice(:,:,:)
31  real, dimension(:), allocatable :: LUT_TAU_VALUES
32 
33  real:: COS_SCAT, COS_VZA, COS_SZA
34 
35 
36 ! ** extra SS stuff
37 
38  real:: prev_ratio
39 
40  REAL, ALLOCATABLE :: fPhase_liq(:,:),tauLR_liq(:,:),EXP0_liq(:,:),EXP1_liq(:,:)
41  REAL, ALLOCATABLE :: fPhase_ice(:,:),tauLR_ice(:,:),EXP0_ice(:,:),EXP1_ice(:,:)
42  REAL, dimension(:), allocatable ::ocean_aero_Omega, ocean_aero_phase
43 
44 ! ** end extra SS stuff
45 
46 
47  contains
48 
49 
50  subroutine lib_init
51 
52  use libraryarrays
53 
54  integer :: nwl, iwl, itau, ntau, i
55  real, dimension(:,:), allocatable :: fqe_water, fqe_ice
56 
57  lib_call_counter = 0
58 
59  nwl = number_wavelengths
60  ntau = number_taus
61  allocate(aerophase(nwl))
62 
63  allocate(fqe_water(number_wavelengths, number_waterradii), fomega_water(number_wavelengths, number_waterradii))
65 
66 ! ** extra SS stuff
67 
68  ALLOCATE(fphase_liq(1:nwl,1:number_waterradii),taulr_liq(1:nwl,1:number_waterradii), &
69  exp0_liq(1:nwl,1:number_waterradii),exp1_liq(1:nwl,1:number_waterradii))
70  ALLOCATE(fphase_ice(1:nwl,1:number_iceradii),taulr_ice(1:nwl,1:number_iceradii), &
71  exp0_ice(1:nwl,1:number_iceradii),exp1_ice(1:nwl,1:number_iceradii))
72 
73  allocate(ocean_aero_omega(nwl), ocean_aero_phase(nwl))
74  ! WDR make re-entrant
75  if( .not. allocated(aero_denom) ) allocate(aero_denom(nwl))
76 
77  do i=1, nwl
78  if (aerosol_tau(i) > 0.) then
79  ocean_aero_omega(i) = ((aerosol_tau(i) -0.1) * 1.0 + 0.1 * aerosol_ssa(i))/aerosol_tau(i)
80  else
81  ocean_aero_omega(i) = 0.
82  endif
83  end do
84 
85  aero_denom = (aerosol_tau -0.1) + 0.1 * aerosol_ssa
86 
87  prev_ratio = -999.
88 
89 
90 ! ** end extra SS stuff
91 
92 
93 
94  fqe_water(1,:) = 1.0
95  do iwl = 2,nwl
96  fqe_water(iwl,:) = extinction_water(iwl,:)/extinction_water(1,:)
97  enddo
98 
99  fomega_water = truncation_factor_water * singlescattering_water !ssalbedo * ftrunc= (omega*f)
100 
101  fqe_ice(1,:) = 1.0
102  do iwl = 2,nwl
103  fqe_ice(iwl,:) = extinction_ice(iwl,:)/extinction_ice(1,:)
104  enddo
105 
106  fomega_ice = truncation_factor_ice * singlescattering_ice !ssalbedo * ftrunc= (omega*f)
107 
108  allocate(tauc_source_water(number_taus+1, number_wavelengths, number_waterradii))
109  allocate(tauc_source_ice(number_taus+1, number_wavelengths, number_iceradii))
110 
111  tauc_source_water(1,:,:) = 0.0
112  tauc_source_ice(1,:,:) = 0.0
113 
114  do itau = 2, number_taus+1
115  tauc_source_water(itau,:,:) = (library_taus(itau-1) * fqe_water) * (1.0 - fomega_water) !Qe sacling, tau scaling, phase func truncation
116  tauc_source_ice(itau,:,:) = (library_taus(itau-1) * fqe_ice) * (1.0 - fomega_ice) !Qe sacling, tau scaling, phase func truncation
117  end do
118 
119  deallocate(fqe_water, fqe_ice)
120 
121 ! ice_int(ntau+1,number_wavelengths,number_iceradii) )
122 
123  allocate( multiscatrefl_water_lamb(ntau+1,number_wavelengths,number_waterradii), &
124  multiscatrefl_ice_lamb(ntau+1,number_wavelengths,number_iceradii), &
125  multiscatrefl_water_cm(ntau+1,number_wavelengths,number_waterradii, 3), &
126  multiscatrefl_ice_cm(ntau+1,number_wavelengths,number_iceradii, 3) )
127 
128  allocate( fphase_lamb(nwl, number_waterradii))
129 
130  allocate( gsquare(nwl))
131 
132 
133 
134  gsquare = aerosol_asym*aerosol_asym
135  aerophase = 0.0
136 
137  num_of_taus = number_taus + 1
138  num_of_wavelengths_rftf = 2
139 
140 
141  allocate(lut_tau_values(num_of_taus))
142  lut_tau_values(1) = 0.
143  lut_tau_values(2:num_of_taus) = library_taus(1:number_taus)
144 
145  ALLOCATE(phasefunvals_ice(number_wavelengths,number_iceradii) )
146  !phaseFunVals_ice = 0. ! WDR UIV again?
147 ! ALLOCATE(ssRefl_water(NUM_OF_TAUS,number_wavelengths,number_waterradii), &
148 ! ssRefl_ice(NUM_OF_TAUS,number_wavelengths,number_iceradii), &
149 ! phaseFunVals_water(number_wavelengths,number_waterradii), &
150 ! phaseFunVals_ice(number_wavelengths,number_iceradii))
151  ALLOCATE(ssrefl_water(num_of_taus,number_wavelengths,number_waterradii))
152  ALLOCATE(phasefunvals_water(number_wavelengths,number_waterradii) )
153  !phaseFunVals_water = 999. ! WDR UIV again?
154  ALLOCATE(ssrefl_ice(num_of_taus,number_wavelengths,number_iceradii) )
155 
156 
157  end subroutine lib_init
158 
159 
160  subroutine lib_clean
161 
162  deallocate(aerophase)
163  deallocate(fomega_water, fomega_ice)
164 
165  deallocate(tauc_source_water, tauc_source_ice)
166 
167  deallocate(multiscatrefl_water_lamb, multiscatrefl_ice_lamb)
168  deallocate(multiscatrefl_water_cm, multiscatrefl_ice_cm)
169 
170  deallocate(fphase_lamb)
171  deallocate(gsquare)
172  deallocate(lut_tau_values)
173  DEALLOCATE(phasefunvals_water, phasefunvals_ice, ssrefl_water, ssrefl_ice)
174 
175 ! ** extra SS stuff
176  deallocate(fphase_liq, taulr_liq, exp0_liq, exp1_liq)
177  deallocate(fphase_ice, taulr_ice, exp0_ice, exp1_ice)
178  deallocate(ocean_aero_omega, ocean_aero_phase)
179 ! ** end extra SS stuff
180 
181 
182  end subroutine lib_clean
183 
184 
185 
186  subroutine libraryinterpolate(local_solarzenith, &
187  local_sensorzenith, &
188  local_relativeazimuth, &
189  local_scatangle, &
190  local_wind_speed, &
191  wind_speed_only, interp_MS, interp_SS, &
192  debug, &
193  status, i, j)
194  !
195  ! so, this interpolates library values to the current conditions of
196  ! 1st 5 inputs
197  ! What it interpolates?
198  use core_arrays
199  use generalauxtype
200  use libraryarrays
203 
204  implicit none
205 
206  real, intent(in) :: local_solarzenith, local_sensorzenith, local_relativeazimuth, local_wind_speed
207  real, intent(in) :: local_scatangle
208  logical,intent(in) :: debug
209  integer, intent(in) :: i, j
210  integer,intent(out) :: status
211  logical, intent(in) :: wind_speed_only, interp_ms, interp_ss
212 
213  real :: dtheta
214  integer :: ianghi, ianglow
215 
216  status = 0
217 
218 
219 ! CONST_PI = 4.d0 * ATAN(1.d0) !value of PI
220 ! RAD_PER_DEG = CONST_PI/180.d0 !radians per degree
221 ! DEG_PER_RAD = 180.d0/CONST_PI !degrees per radian
222 
223  cos_vza = cos(local_sensorzenith * d2r)
224  cos_sza = cos(local_solarzenith * d2r)
225  cos_scat = cos(local_scatangle * d2r)
226 ! NUM_OF_TAUS = number_taus + 1
227 ! NUM_OF_WAVELENGTHS_RFTF = 2
228 
229 ! For Lambertian libraries: land, also ocean when 0.86um saturates.
230  if (.not. cox_munk) then
231 
234 
235  call getrefl(cos_sza, cos_vza, local_relativeazimuth, local_scatangle, &
236  int_reflectance_water(:,:,:), int_reflectance_ice(:,:,:), interp_ms, interp_ss, status)
237  if( status == 5 ) return ! outside table
238 
239  if (interp_ms) then
240 
243 
244  call getsdevrefllamb(cos_sza, cos_vza, local_relativeazimuth,&
246  int_reflectance_ice_sdev(:,:,:), &
247  status)
248  if( status == 5 ) return ! out of table space
249  call setup_emissivity_flux(cos_sza, library_fluxsolarzenith, ianghi, ianglow, dtheta)
250 
252  ianghi, ianglow, dtheta, status)
254  ianghi, ianglow, dtheta, status)
255 
256  call setup_emissivity_flux(cos_vza, library_fluxsensorzenith, ianghi, ianglow, dtheta)
257 
259  ianghi, ianglow, dtheta, status)
261  ianghi, ianglow, dtheta, status)
262 
264  ianghi, ianglow, dtheta, status)
266  ianghi, ianglow, dtheta, status)
267 
268  endif
269 
270  else
271 
274 
275  call getreflforgivenwindspeed(cos_sza, cos_vza, local_relativeazimuth, &
276  local_scatangle, cos_scat, local_wind_speed, &
278  wind_speed_only, interp_ms, interp_ss, status)
279  if( status == 5 ) return
280 
281  if (interp_ms) then
282 
285 
286  call getsdevreflforgivenwindspeed(cos_sza, cos_vza, &
287  local_relativeazimuth, local_wind_speed, &
288  int_refl_water_sdev_wspeed(:,:,:,:), &
290  int_refl_ice_sdev_wspeed(:,:,:,:), &
291  int_reflectance_ice_sdev(:,:,:), &
292  wind_speed_only, status)
293  if( status == 5 ) return ! point beyond table range
294 
295  call setup_emissivity_flux(cos_vza, library_sensor_zenith, ianghi, ianglow, dtheta)
296 
297  call interpemissforagivenwspeed(cos_vza, local_wind_speed, surface_emissivity_water,&
299  wind_speed_only, ianghi, ianglow, dtheta, status)
300  call interpemissforagivenwspeed(cos_vza, local_wind_speed, surface_emissivity_ice, &
302  wind_speed_only, ianghi, ianglow, dtheta, status)
303  call interpemissforagivenwspeed(cos_vza, local_wind_speed, cloud_emissivity_water, &
305  wind_speed_only, ianghi, ianglow, dtheta, status)
306  call interpemissforagivenwspeed(cos_vza, local_wind_speed, cloud_emissivity_ice, &
308  wind_speed_only, ianghi, ianglow, dtheta, status)
309 
310 ! These quantities are not used for Collection 6 calculations
311  call interpemissforagivenwspeed(cos_vza, local_wind_speed, surface_emissivity_water_sdev,&
313  wind_speed_only, ianghi, ianglow, dtheta, status)
314  call interpemissforagivenwspeed(cos_vza, local_wind_speed, surface_emissivity_ice_sdev, &
316  wind_speed_only, ianghi, ianglow, dtheta, status)
317 
318  call interpemissforagivenwspeed(cos_vza, local_wind_speed, cloud_emissivity_water_sdev, &
320  wind_speed_only, ianghi, ianglow, dtheta, status)
321  call interpemissforagivenwspeed(cos_vza, local_wind_speed, cloud_emissivity_ice_sdev, &
323  wind_speed_only, ianghi, ianglow, dtheta, status)
324  endif
325 
326  endif
327 
328 
329 
330 
331 end subroutine libraryinterpolate
332 
333 REAL FUNCTION scatangle(solarAng,viewAng,relAzm)
334 !................................................................................
335 !
336 !................................................................................
337 ! !F90
338 !
339 ! !DESCRIPTION:
340 ! This function computes the scattering angle for a given sun-satellite geometry!
341 !
342 ! !INPUT PARAMETERS:
343 ! solarAng : solar angle in degrees (scalar)
344 ! viewAng : view angle in degrees (scalar)
345 ! relAzm : rel azimuth (from the principle plane of the SUN) (scalar)
346 !
347 !
348 ! !OUTPUT PARAMETERS:
349 ! scatAngle : scattering angle in degrees (scalar)
350 !
351 !................................................................................
352 ! !Revision History:
353 ! Revision 1.0 2009/06/03 GW
354 ! Initial integration of the standalone version. No changes from the standalone.
355 !
356 !................................................................................
357 ! !Team-unique Header:
358 ! Cloud Retrieval Group, NASA/GSFC
359 !
360 ! !PROGRAMMER:
361 !
362 ! Nandana Amarasinghe (SSAI)
363 ! Climate and Radiation Branch
364 ! NASA Goddard Space Flight Center
365 ! Greenbelt, Maryland 20771, U.S.A.
366 !
367 ! Gala Wind (SSAI)
368 ! Climate and Radiation Branch
369 ! NASA Goddard Space Flight Center
370 ! Greenbelt, Maryland 20771, U.S.A.
371 !
372 !*******************************************************************************
373 ! !END
374  use science_parameters, only:d2r
375 
376  IMPLICIT NONE
377  !CALLING ARGUMENTS
378  REAL, INTENT(IN):: solarang,viewang,relazm
379 
380  !LOCAL VARIABLES
381  REAL::solarmu,viewmu !,pi,torad,todeg
382 
383 ! pi = acos(-1.d0 )
384 ! torad = pi / 180.d0
385 ! todeg = 180.d0 / pi
386 
387  solarmu = cos(solarang * d2r)
388  viewmu = cos(viewang * d2r)
389 
390  scatangle= -solarmu*viewmu + sqrt( ( 1.- solarmu**2 )* &
391  ( 1.- viewmu**2 ) )*cos( relazm*d2r)
392 
393 
394  IF(scatangle < -1.0)scatangle = -1.0 !guard against round off error
395  IF(scatangle > 1.0)scatangle = 1.0
396 
397 
398  scatangle = acos(scatangle) / d2r !* todeg
399 
400 END FUNCTION scatangle
401 
402 
403 SUBROUTINE getphasefunctionvalues(scatAngle, scatAngleArray, num_angles, phaseFunArray, phaseNormConst, &
404  phaseFunVals, ierror)
405 !................................................................................
406 !
407 !................................................................................
408 ! !F90
409 !
410 ! !DESCRIPTION:
411 ! This program interpolates the phase fun value array for a given scattering angle an returns
412 ! the corresponding values for particular lambda and re (see the inputs and outputs)!
413 ! !INPUT PARAMETERS:
414  ! scatAngle : scattering angle in degrees (scalar)
415  ! scatAngleArray : 1-D array of sacattering angles for which phase fun vals are precomputed
416  ! phaseFunArray : 3-D array of phase fun vals: 1st DIM = num of angles,
417  ! 2nd DIM = num of wavelengths, 3rd DIM = num of REs
418  ! phaseNormConst : Normalization factor for the phase function. If it is normalized to 1
419  ! then it is OK. But, this number (= pmom(0)) from the dfit routines has
420  ! been saved and included in the library file in adition to d_fit coeffs.
421 !
422 !
423 ! !OUTPUT PARAMETERS:
424  ! phaseFunVals : 2-D array of interpolated phase fun values at scatAngle
425  ! 1st dim = num of wavelengths, 2nd dim = num of rE's
426  ! ierror : integer value
427  ! ierror = 0 success,
428  ! = 3 invalid data
429 !................................................................................
430 ! !Revision History:
431 ! Revision 1.0 2009/06/03 GW
432 ! Initial integration of the standalone version. No changes from the standalone.
433 !
434 !................................................................................
435 ! !Team-unique Header:
436 ! Cloud Retrieval Group, NASA/GSFC
437 !
438 ! !PROGRAMMER:
439 !
440 ! Nandana Amarasinghe (SSAI)
441 ! Climate and Radiation Branch
442 ! NASA Goddard Space Flight Center
443 ! Greenbelt, Maryland 20771, U.S.A.
444 !
445 ! Gala Wind (SSAI)
446 ! Climate and Radiation Branch
447 ! NASA Goddard Space Flight Center
448 ! Greenbelt, Maryland 20771, U.S.A.
449 !
450 !*******************************************************************************
451 ! !END
452 
454  use core_arrays, only: platform_name
455  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
456  IMPLICIT NONE
457  !Arguments
458  REAL,INTENT(IN)::scatAngle, scatAngleArray(:),phaseFunArray(:,:,:),phaseNormConst(:,:)
459  REAL,INTENT(OUT)::phaseFunVals(:,:)
460  integer, intent(in):: num_angles
461  INTEGER,INTENT(OUT)::ierror
462 
463  !local variables
464  INTEGER:: iAng, iAngLow, iAngHi, nwl
465  REAL::dtheta
466  CHARACTER(150)::tempstr
467 
468  !start execution
469  ! WDR 13 Jun 2022 add the OCI platform_name in if()
470  if (platform_name == "AVIRIS" .or. platform_name == "RSP" .or. &
471  platform_name == "AHI" .or. platform_name == "EPIC" .or. &
472  platform_name == "SSFR" .or. c2_sensor_id == oci_id .or. &
473  c2_sensor_id == ocis_id ) then
474  nwl = number_wavelengths
475  else
476  nwl = number_wavelengths-1
477  endif
478 
479 
480  ierror = 0 !start with no error
481  ianglow = 1
482  ianghi = num_angles
483  IF(scatangle > scatanglearray(ianghi))THEN
484  ianglow = ianghi-1 !no need to do the search; do extrapolation
485 
486  ELSEIF(scatangle < scatanglearray(ianglow))THEN
487  ianghi = ianglow+1
488  ELSE
489  DO !do a bisection search
490  IF((ianghi - ianglow) == 1)EXIT
491  iang = (ianghi + ianglow)/2
492  IF(scatanglearray(iang) < scatangle)THEN
493  ianglow = iang
494  ELSE
495  ianghi = iang
496  ENDIF
497  ENDDO
498  ENDIF
499 
500 
501  dtheta = scatanglearray(ianghi) - scatanglearray(ianglow)
502  !print*,iAngHi,iAngLow,dtheta,scatAngleArray(iAngHi)-scatAngle,scatAngle - scatAngleArray(iAngLow)
503 ! IF(dtheta ==0)THEN
504 ! ierror = 3 !invalid data read
505 ! RETURN
506 ! ENDIF
507 
508  phasefunvals(1:nwl,:) = ( phasefunarray(ianglow,1:nwl,:) * (scatanglearray(ianghi)-scatangle) + &
509  phasefunarray(ianghi,1:nwl,:) * (scatangle - scatanglearray(ianglow) ))/dtheta
510 
511 
512  phasefunvals(1:nwl,:) = phasefunvals(1:nwl,:)/phasenormconst(1:nwl,:)
513 
514  RETURN
515 
516 END SUBROUTINE getphasefunctionvalues
517 
518 
519 subroutine get_aero_params(cos_scatAngle, aeroG)
521  use core_arrays, only: platform_name
522  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
523 
524  real, intent(in) :: aeroG(:), cos_scatAngle
525 
526  INTEGER:: nwl
527 
528  if (platform_name == "AVIRIS" .or. platform_name == "RSP" .or. &
529  platform_name == "AHI" .or. platform_name == "EPIC" .or. &
530  platform_name == "SSFR" .or. c2_sensor_id == oci_id .or. &
531  c2_sensor_id == ocis_id ) then
532  nwl = number_wavelengths
533  else
534  nwl = number_wavelengths-2
535  endif
536 
537  rlphase = 0.75 * (1 + cos_scatangle*cos_scatangle)
538 
539 ! it is generally advisable to avoid a real power whenever possible in order to make the code run faster
540 ! GW: 10.25.13
541  aerophase(1:nwl) = (1.0 - gsquare(1:nwl))/ &
542  (sqrt(1 + gsquare(1:nwl) &
543  - 2.0 * aerog(1:nwl) * cos_scatangle)**3)
544 
545 end subroutine get_aero_params
546 
547 
548 SUBROUTINE interpolatefluxes(solarOrViewAng, solarAngMuArray, inFluxArray, outFluxArray, iAngHi, iAngLow, dtheta, ierror)
549 !................................................................................
550 ! !F90
551 !
552 ! !DESCRIPTION:
553  !This program interpolates(linear) the RF or TF array for a given solar/oview angle and returns
554  !the corresponding values for particular tau, lambda and re (see the inputs and outputs)
555 ! !INPUT PARAMETERS:
556  ! solarOrViewAng : input angle in degrees (scalar)
557  ! solarAngleMuArray : 1-D array of solar angles for which RF/TF vals are precomputed
558  ! inFluxArray : 4-D array of flux vals: 1st DIM = num of angles, 2nd DIM = ntaus
559  ! 3rdd DIM = num of wavelengths, 4th DIM = num of REs
560 !
561 !
562 ! !OUTPUT PARAMETERS:
563  ! outFluxArray : 3-D array of interpolated RF or TF values at the input angle
564  ! 1st dim = ntaus, 2nd dim = num of wavelengths, 3rdd dim = num of rE's
565  ! ierror : integer value
566  ! ierror = 0 success,
567  ! = 3 invalid data
568 !................................................................................
569 ! !Revision History:
570 ! Revision 1.0 2009/06/03 GW
571 ! Initial integration of the standalone version. No changes from the standalone.
572 !
573 !................................................................................
574 ! !Team-unique Header:
575 ! Cloud Retrieval Group, NASA/GSFC
576 !
577 ! !PROGRAMMER:
578 !
579 ! Nandana Amarasinghe (SSAI)
580 ! Climate and Radiation Branch
581 ! NASA Goddard Space Flight Center
582 ! Greenbelt, Maryland 20771, U.S.A.
583 !
584 ! Gala Wind (SSAI)
585 ! Climate and Radiation Branch
586 ! NASA Goddard Space Flight Center
587 ! Greenbelt, Maryland 20771, U.S.A.
588 !
589 !*******************************************************************************
590 ! !END
591  IMPLICIT NONE
592  !Arguments
593  REAL,INTENT(IN):: solarOrViewAng, solarAngMuArray(:),inFluxArray(:,:,:,:)
594  REAL,INTENT(OUT)::outFluxArray(:,:,:)
595  INTEGER,INTENT(OUT)::ierror
596  integer, intent(in) :: iAngLow, iAngHi
597  real, intent(in) :: dtheta
598 
599  !local variables
600  real*8::umu
601 
602  real :: diff_one_way, diff_other_way
603 
604  !start esecution
605  ierror = 0
606 
607  umu = solarorviewang
608 
609  diff_one_way = solarangmuarray(ianghi)- umu
610  diff_other_way = umu - solarangmuarray(ianglow)
611 
612  outfluxarray(:,:,:) = 0.
613  outfluxarray(:,:,:) = ( influxarray(ianglow,:,:,:) * diff_one_way + &
614  influxarray(ianghi,:,:,:) * diff_other_way ) / dtheta
615 
616  RETURN
617 
618 END SUBROUTINE interpolatefluxes
619 
620  SUBROUTINE calcsingscatpartofreflectance(phaseFunVals_water, phaseFunVals_ice, &
621  theta, theta0, ssRefl_water, ssRefl_ice, ierror)
622 
623 !................................................................................
624 ! !F90
625 !
626 ! !DESCRIPTION:
627  !Calculates the single scattering part of the reflectance according to Nakajima/Tanaka method
628  !implemented in DISORT Version2. We had to implement this way, because of the way single
629  !scattering part is takent out in DISORT. Please refer to the DISORT2 documentation
630  !Formula Used here: refl = omega * PC * {1.0 - EXP[-tauc *(1/mu + 1/mu0)]}
631  ! refl = refl/(4.0 * (mu + mu0))
632  ! where tauc = tau * [Qe(lambda,re)/Qe(CH1,re)] * (1.0 - ftrunc * omega)
633  ! PC = Ptrue / (1 - ftrunc * omega)
634  ! also remember refl_lUT = (I_disort) * pi/( mu0* solarflux)
635  ! so that multiplication is already taken care of here, no need to do it
636  ! outside of this routine
637 ! !INPUT PARAMETERS:
638  ! PhaseFunVals : phase function calculated for all the wavelengths and effective radii
639  ! DIMS[1:nlambda,1:nradii] (look at the program getPhaseFunValues.f90)
640  ! ssAlbedo : single scattering albedo DIMS[1:nlambda,1:nradii]
641  ! extCoeff : extinction coefficient DIMS[1:nlambda,1:nradii]
642  ! dfitF : phase function truncation factor in dfit method, DIMS[1:nlambda,1:nradii]
643  ! tauVals : opticat thickness values in the LUT, DIMS[1:ntau]
644  ! theta : view angle in degrees
645  ! theta0 : solar angle in degrees
646 !
647 !
648 ! !OUTPUT PARAMETERS:
649  ! ssRefl : single scattering part of the reflectance , here no need to input the beam
650  ! strength, since we calculate the reflectance.
651  ! ierror : integer value
652  ! ierror = 0 success,
653  ! = 2 memory allocation error
654  ! = 3 invalid data
655 !................................................................................
656 ! !Revision History:
657 ! Revision 1.0 2009/06/03 GW
658 ! Initial integration of the standalone version. No changes from the standalone.
659 !
660 !................................................................................
661 ! !Team-unique Header:
662 ! Cloud Retrieval Group, NASA/GSFC
663 !
664 ! !PROGRAMMER:
665 !
666 ! Nandana Amarasinghe (SSAI)
667 ! Climate and Radiation Branch
668 ! NASA Goddard Space Flight Center
669 ! Greenbelt, Maryland 20771, U.S.A.
670 !
671 ! Gala Wind (SSAI)
672 ! Climate and Radiation Branch
673 ! NASA Goddard Space Flight Center
674 ! Greenbelt, Maryland 20771, U.S.A.
675 !
676 !*******************************************************************************
677 ! !END
678  use libraryarrays
679  use core_arrays, only: platform_name
680  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
681 
682  REAL, INTENT(IN) :: phaseFunVals_water(:,:), phaseFunVals_ice(:,:)
683  REAL, INTENT(IN) :: theta, theta0
684  REAL, INTENT(OUT) :: ssRefl_water(:,:,:), ssRefl_ice(:,:,:)
685  INtEGER,INTENT(OUT)::ierror
686 
687 
688  INTEGER:: ntau, nwl, nre, itau,iwl, allocateStatus, ii,jj
689  REAL :: MUPLUSMU0,ratio, mytau
690 
691  ntau = number_taus !size(tauVals)
692 
693  if( platform_name == "AVIRIS" .or. platform_name == "RSP" .or. &
694  platform_name == "AHI" .or. platform_name == "EPIC" .or. &
695  platform_name == "SSFR" .or. c2_sensor_id == oci_id .or. &
696  c2_sensor_id == ocis_id ) then
697  nwl = number_wavelengths + 1 ! we have to make sure not to skip the last channel
698  else
699  nwl = number_wavelengths
700  endif
701 
702 
703  muplusmu0 = theta + theta0
704  ratio = muplusmu0/(theta * theta0)
705 
706 
707  fphase_lamb = phasefunvals_water/ (1.0 - fomega_water) !add something here to safeguard div by zero
708 
709  ssrefl_water = 0.0
710  ssrefl_ice = 0.0
711 
712  do itau = 1, ntau
713 
714  do ii=1, nwl-1
715  do jj=1, number_waterradii
716 
717  mytau = tauc_source_water(itau+1,ii,jj)*ratio
718 
719  if (mytau == 0) then
720  exp1_liq(ii,jj) = 1.
721  else if (mytau > 10.) then
722  exp1_liq(ii,jj) = 0.
723  else
724  exp1_liq(ii,jj) = exp(-mytau)
725 
726 ! EXP1_liq(ii,jj) = 1. + (mytau ) + (mytau )**2/2. + (mytau )**3/6. + &
727 ! (mytau )**4/24. + (mytau )**5/120. + mytau**6/720.
728 ! EXP1_liq(ii,jj) = 1. / EXP1_liq(ii,jj)
729  endif
730 
731  ssrefl_water(itau,ii,jj) = (singlescattering_water(ii,jj)* fphase_lamb(ii,jj) * (1 - exp1_liq(ii,jj))) / (4.0*muplusmu0)
732 
733  end do
734  end do
735 
736 
737  enddo
738 
739  fphase_lamb(:, 1:number_iceradii) = phasefunvals_ice/ (1.0 - fomega_ice) !add something here to safeguard div by zero
740 
741  do itau = 1, ntau
742 
743 
744  do ii=1, nwl-1
745  do jj=1, number_iceradii
746 
747 
748  mytau = tauc_source_ice(itau+1,ii,jj)*ratio
749 
750  if (mytau == 0) then
751  exp1_ice(ii,jj) = 1.
752  else if (mytau > 10.) then
753  exp1_ice(ii,jj) = 0.
754  else
755  exp1_ice(ii,jj) = exp(-mytau)
756  endif
757 
758  ssrefl_ice(itau,ii,jj) = (singlescattering_ice(ii,jj)* fphase_lamb(ii,jj) * (1 - exp1_ice(ii,jj))) / (4.0*muplusmu0)
759 
760  end do
761  end do
762 
763 
764  enddo
765 
766  END SUBROUTINE calcsingscatpartofreflectance
767 
768 
769  SUBROUTINE single_scattering_calcs_ocean(phaseFunVals_liq, phaseFunVals_ice, ssAlbedo_liq, ssAlbedo_ice, &
770  RLphase, aeroPhase, RLTau,aeroTau, aeroOmega,&
771  theta, theta0,ssRefl_liq, ssRefl_ice)
772 
774  use core_arrays, only: platform_name
775  use ch_xfr, only : c2_sensor_id, oci_id, ocis_id
776 
777  REAL, INTENT(IN) :: phaseFunVals_liq(:,:), ssAlbedo_liq(:,:), &
778  phaseFunVals_ice(:,:), ssAlbedo_ice(:,:), &
779  RLTau(:), aeroTau(:), &
780  aeroOmega(:), RLphase, aeroPhase(:)
781  REAL, INTENT(IN) :: theta, theta0
782  REAL, INTENT(OUT) :: ssRefl_liq(:,:,:), ssRefl_ice(:,:,:)
783 
784 
785  INTEGER:: ntau, nwl, itau, iwl, wave_limit
786  REAL :: MUPLUSMU0,ratio
787 
788  real :: ocean_mul(20)
789  logical:: changed_ratio
790  integer :: ii, jj
791 
792  ntau = num_of_taus
793 
794  if (platform_name == "AVIRIS" .or. platform_name == "RSP" .or. &
795  platform_name == "AHI" .or. platform_name == "EPIC" .or. &
796  platform_name == "SSFR" .or. c2_sensor_id == oci_id .or. &
797  c2_sensor_id == ocis_id ) then
798  nwl = number_wavelengths+1 ! capture the last channel, it is not to be skipped
799  else
800  nwl = number_wavelengths
801  endif
802 
803  muplusmu0 = theta + theta0
804  ratio = muplusmu0/(theta * theta0)
805 
806  fphase_liq = phasefunvals_liq/ (1.0 - fomega_water) !add something here to safeguard div by zero
807  fphase_ice = phasefunvals_ice/ (1.0 - fomega_ice)
808 
809 
810 #if AVIRIS_INST | RSP_INST | AHI_INST | EPIC_INST | SSFR_INST
811  do iwl = 1, number_wavelengths
812 #else
813  if( ( c2_sensor_id == oci_id ) .or. ( c2_sensor_id == ocis_id ) ) then
814  wave_limit = number_wavelengths - 3
815  else
816  wave_limit = number_wavelengths - 2
817  endif
818  do iwl = 1, wave_limit
819 #endif
820  ocean_aero_phase(iwl) = ( (aerotau(iwl) -0.1) * rlphase + 0.1 * aeroomega(iwl) * aerophase(iwl) ) / aero_denom(iwl)
821  ocean_mul(iwl) = ocean_aero_phase(iwl) * ocean_aero_omega(iwl)
822  end do
823 
824  exp1_liq = 0.
825  exp1_ice = 0.
826 
827  DO itau = 1, ntau
828 
829  taulr_liq = tauc_source_water(itau,:,:) * ratio
830  taulr_ice = tauc_source_ice(itau,:,:) * ratio
831 
832 
833  do ii=1, nwl-1
834  do jj=1, number_waterradii
835 
836  if (taulr_liq(ii,jj) == 0) then
837  exp1_liq(ii,jj) = 1.
838  else if (taulr_liq(ii,jj) > 10.) then
839  exp1_liq(ii,jj) = 0.
840  else
841  exp1_liq(ii,jj) = exp(-taulr_liq(ii,jj))
842  endif
843 
844  end do
845  end do
846 
847  do ii=1, nwl-1
848  do jj=1, number_iceradii
849 
850  if (taulr_ice(ii,jj) == 0) then
851  exp1_ice(ii,jj) = 1.
852  else if (taulr_ice(ii,jj) > 10.) then
853  exp1_ice(ii,jj) = 0.
854  else
855  exp1_ice(ii,jj) = exp(-taulr_ice(ii,jj))
856  endif
857 
858  end do
859  end do
860 
861 
862  ssrefl_liq(itau,:,:) = ssalbedo_liq * fphase_liq * (1. - exp1_liq)
863  ssrefl_ice(itau,:,:) = ssalbedo_ice * fphase_ice * (1. - exp1_ice)
864 
865  exp0_liq = exp1_liq
866  exp0_ice = exp1_ice
867 
868 #if AVIRIS_INST | RSP_INST | AHI_INST | EPIC_INST | SSFR_INST
869  DO iwl =1, number_wavelengths
870 #else
871  if( ( c2_sensor_id == oci_id ) .or. ( c2_sensor_id == ocis_id ) ) then
872  wave_limit = number_wavelengths - 3
873  else
874  wave_limit = number_wavelengths - 2
875  endif
876  do iwl = 1, wave_limit
877 #endif
878  taulr_liq(iwl,:) = (tauc_source_water(itau,iwl,:) + rltau(iwl))*ratio
879  taulr_ice(iwl,:) = (tauc_source_ice(itau,iwl,:) + rltau(iwl))*ratio
880 
881  ii = iwl
882  do jj=1, number_waterradii
883 
884  if (taulr_liq(ii,jj) == 0) then
885  exp1_liq(ii,jj) = 1.
886  else if (taulr_liq(ii,jj) > 10.) then
887  exp1_liq(ii,jj) = 0.
888  else
889  exp1_liq(ii,jj) = exp(-taulr_liq(ii,jj))
890  endif
891 
892  end do
893 
894  ii = iwl
895  do jj=1, number_iceradii
896 
897  if (taulr_ice(ii,jj) == 0) then
898  exp1_ice(ii,jj) = 1.
899  else if (taulr_ice(ii,jj) > 10.) then
900  exp1_ice(ii,jj) = 0.
901  else
902  exp1_ice(ii,jj) = exp(-taulr_ice(ii,jj))
903  endif
904 
905  end do
906 
907 
908  ssrefl_liq(itau,iwl,:) = ssrefl_liq(itau,iwl,:) + &
909  rlphase * 1.0 * (exp0_liq(iwl,:) - exp1_liq(iwl,:))
910  ssrefl_ice(itau,iwl,:) = ssrefl_ice(itau,iwl,:) + &
911  rlphase * 1.0 * (exp0_ice(iwl,:) - exp1_ice(iwl,:))
912 
913  exp0_liq(iwl,:) = exp1_liq(iwl,:)
914  exp0_ice(iwl,:) = exp1_ice(iwl,:)
915 
916  taulr_liq(iwl,:) = (tauc_source_water(itau,iwl,:) + rltau(iwl) + aerotau(iwl)) * ratio
917  taulr_ice(iwl,:) = (tauc_source_ice(itau,iwl,:) + rltau(iwl) + aerotau(iwl)) * ratio
918 
919  !coming in aeroTau is the (rayleightau + 0.1) from DISORT run
920 
921  ii = iwl
922  do jj=1, number_waterradii
923 
924  if (taulr_liq(ii,jj) == 0) then
925  exp1_liq(ii,jj) = 1.
926  else if (taulr_liq(ii,jj) > 10.) then
927  exp1_liq(ii,jj) = 0.
928  else
929  exp1_liq(ii,jj) = exp(-taulr_liq(ii,jj))
930  endif
931 
932  end do
933 
934  ii = iwl
935  do jj=1, number_iceradii
936 
937  if (taulr_ice(ii,jj) == 0) then
938  exp1_ice(ii,jj) = 1.
939  else if (taulr_ice(ii,jj) > 10.) then
940  exp1_ice(ii,jj) = 0.
941  else
942  exp1_ice(ii,jj) = exp(-taulr_ice(ii,jj))
943  endif
944 
945  end do
946 
947 
948 
949  ssrefl_liq(itau,iwl,:) = ssrefl_liq(itau,iwl,:) + ocean_mul(iwl) * &
950  (exp0_liq(iwl,:) - exp1_liq(iwl,:))
951  ssrefl_ice(itau,iwl,:) = ssrefl_ice(itau,iwl,:) + ocean_mul(iwl) * &
952  (exp0_ice(iwl,:) - exp1_ice(iwl,:))
953  ENDDO
954 
955  ENDDO
956 
957 
958  ssrefl_liq = ssrefl_liq/(4.0 * muplusmu0)
959  ssrefl_ice = ssrefl_ice/(4.0 * muplusmu0)
960 
961 
962 
963  END SUBROUTINE single_scattering_calcs_ocean
964 
965 
966 
967 
968 
969  !----------------------SUBROUTINE CalcSingScatPartOfOceanReflectance---------------------------------
970  !Calculates the single scattering part of the reflectance according to Nakajima/Tanaka method
971  !implemented in DISORT Version2. We had to implement this way, because of the way single
972  !scattering part is takent out in DISORT. Please refer to the DISORT2 documentation
973  !Formula Used here: refl = omega * PC * {1.0 - EXP[-tauc *(1/mu + 1/mu0)]}
974  ! refl = refl/(4.0 * (mu + mu0))
975  ! where tauc = tau * [Qe(lambda,re)/Qe(CH1,re)] * (1.0 - ftrunc * omega)
976  ! PC = Ptrue / (1 - ftrunc * omega)
977  ! also remember refl_lUT = (I_disort) * pi/( mu0* solarflux)
978  ! so that multiplication is already taken care of here, no need to do it
979  ! outside of this routine
980  !INPUTS:
981  ! PhaseFunVals : phase function calculated for all the wavelengths and effective radii
982  ! DIMS[1:nlambda,1:nradii] (look at the program getPhaseFunValues.f90)
983  ! RLphase : Rayleigh phase func value (fn of the scattering angle)
984  ! aeroPhse : mixed phase function of rayleigh layer and the aerosol layer
985  ! ssAlbedo : single scattering albedo DIMS[1:nlambda,1:nradii]
986  ! extCoeff : extinction coefficient DIMS[1:nlambda,1:nradii]
987  ! dfitF : phase function truncation factor in dfit method, DIMS[1:nlambda,1:nradii]
988  ! tauVals : opticat thickness values in the LUT, DIMS[1:ntau]
989  ! RLTau : Rayleigh optical thickness (combined several layers upto 8km)
990  ! aeroTau : aerosol optical thickness(0.1) + rayleigh layer(lowest)
991  ! aeroOmega : ss albedo of aerosol (not the mixed, mixing will be done in SS calc routine)
992  ! theta : view angle in degrees
993  ! theta0 : solar angle in degrees
994  !OUTPUTS:
995  ! ssRefl : single scattering part of the reflectance , here no need to input the beam
996  ! strength, since we calculate the reflectance.
997  ! ierror : integer value
998  ! ierror = 0 success,
999  ! = 2 memory allocation error
1000  ! = 3 invalid data
1001  !-----------------------------------------------------------------------------------------------
1002 
1003  subroutine setup_emissivity_flux(angle, angle_array, idx_hi, idx_lo, dtheta)
1005  real, intent(in) :: angle, angle_array(:)
1006  real, intent(out) :: dtheta
1007  integer, intent(out) :: idx_hi, idx_lo
1008 
1009 
1010  real :: UMU
1011  integer :: iAng
1012 
1013  umu = angle !COS(solarOrViewAng * RAD_PER_DEG)
1014  idx_lo = 1
1015  idx_hi = size(angle_array)
1016  IF(umu < angle_array(idx_lo))THEN
1017  idx_hi = idx_lo+1 !no need to do the search; do extrapolation
1018  if (idx_hi > size(angle_array)) idx_hi = idx_lo ! we don't have enough points
1019  ELSE
1020  DO !do a bisection search
1021  IF((idx_hi - idx_lo) == 1)EXIT
1022  iang = (idx_hi + idx_lo)/2
1023  IF(angle_array(iang) < umu)THEN
1024  idx_lo = iang
1025  ELSE
1026  idx_hi = iang
1027  ENDIF
1028  ENDDO
1029  ENDIF
1030 
1031 
1032  dtheta = angle_array(idx_hi) - angle_array(idx_lo)
1033 
1034  end subroutine setup_emissivity_flux
1035 
1036 
1037  !--------------------SUBROUTINE InterpCloudEmissForAGivenWSpeed ---------------------------------
1038  !This program interpolates(linear) the cloud emissivity arrays for a given wind speed and
1039  !returns the corresponding values for particular tau, lambda and re(see the inputs and outputs)
1040  !INPUTS:
1041  ! iphase : integer 1 =Water Cloud, 2 =Ice Cloud
1042  ! solarOrViewAng : input angle in degrees (scalar)
1043  ! wspeed : real scalar, wind speed in m/s
1044  !OUTPUTS:
1045  ! outFluxArray : 3-D array ofinterpolated cloud emissivity values at the input ang
1046  ! 1st dim = ntaus, 2nd dim = num of wavelengths, 3rdd dim = num of rE's
1047  ! ierror : integer value
1048  ! ierror = 0 success,
1049  ! = 3 invalid data
1050  !-----------------------------------------------------------------------------------------------
1051  !
1052  SUBROUTINE interpemissforagivenwspeed(solarOrViewAng,wspeed, inEmissArray, outEmissArray_ws, outEmissArray, &
1053  wspeed_only, iAngHi, iAngLow, dtheta,ierror)
1054 
1055  use libraryarrays
1056 
1057  IMPLICIT NONE
1058  !Arguments
1059  REAL,INTENT(IN):: solarOrViewAng,wspeed, dtheta
1060  real, intent(in) :: inEmissArray(:,:,:,:,:)
1061  REAL,INTENT(OUT)::outEmissArray(:,:,:)
1062  real, intent(inout) :: outEmissArray_ws(:,:,:,:)
1063  INTEGER,INTENT(OUT)::ierror
1064  logical, intent(in) :: wspeed_only
1065  integer, intent(in) :: iAngHi, iAngLow
1066 
1067 
1068  real :: diff_one_way, diff_other_way
1069 
1070 
1071  diff_one_way = library_sensor_zenith(ianghi)- solarorviewang
1072  diff_other_way = solarorviewang - library_sensor_zenith(ianglow)
1073 
1074  outemissarray_ws(:,:,:,1) = (inemissarray(ianglow,:,:,:,1) * diff_one_way + &
1075  inemissarray(ianghi,:,:,:,1) * diff_other_way ) / dtheta
1076  outemissarray_ws(:,:,:,2) = (inemissarray(ianglow,:,:,:,2) * diff_one_way + &
1077  inemissarray(ianghi,:,:,:,2) * diff_other_way) / dtheta
1078  outemissarray_ws(:,:,:,3) = (inemissarray(ianglow,:,:,:,3) * diff_one_way + &
1079  inemissarray(ianghi,:,:,:,3) * diff_other_way) / dtheta
1080 
1081  outemissarray = 0.
1082  call interpolate_wind_speed(wspeed, outemissarray_ws, outemissarray)
1083 
1084  RETURN
1085 
1086  END SUBROUTINE interpemissforagivenwspeed
1087 
1088 
1089 ! SUBROUTINE TriLinearInterpolationRefl(solarAng,viewAng,azmAng, &
1090 ! solarMuArray, viewMuArray, relAzmArray, &
1091 ! multiScatReflLUT_water, multiScatRefl_water, &
1092 ! multiScatReflLUT_ice, multiScatRefl_ice,ierror)
1093 !!................................................................................
1094 !! !F90
1095 !!
1096 !! !DESCRIPTION:
1097 ! !perform 3-D linear interpolation on multiple scattering part of the reflectance. Instead of
1098 ! !indexing, a search is done to find the the indices of the cube that enclose the interpolation
1099 ! !point.
1100 !! !INPUT PARAMETERS:
1101 ! ! solarAng : solar zenith angle in degrees (scalar)
1102 ! ! viewAng : sensor zenith angle in degrees (scalar)
1103 ! ! azmAng : rel azm ang in degrees (scalar)]
1104 ! ! solarMuArray : array of solar zenith angles (from the LUT) 1-D array
1105 ! ! viewMuArray : array of sensor zenith angles (from the LUT) 1D array
1106 ! ! relAzmArray : array of rel azimuth angles (from the LUT) 1-D array
1107 ! ! multiScatReflLUT: multi scat part of the reflectance (from LUT) 6D array
1108 ! ! DIMS[1:nview, 1:nsolar,1:nazm, 1:ntau, 1:nlambda, 1:nre]
1109 !!
1110 !!
1111 !! !OUTPUT PARAMETERS:
1112 ! ! multiScatRefl : interpolated MS reflectance for the input geometry 3D array
1113 ! ! DIMS[1:ntau, 1:nlambda, 1:nre]
1114 ! ! ierror : integer value
1115 ! ! ierror = 0 success,
1116 ! ! = 1 I/O error
1117 ! ! = 2 memory allocation error
1118 ! ! = 3 invalid data
1119 !!................................................................................
1120 !! !Revision History:
1121 !! Revision 1.0 2009/06/03 GW
1122 !! Initial integration of the standalone version. No changes from the standalone.
1123 !!
1124 !!................................................................................
1125 !! !Team-unique Header:
1126 !! Cloud Retrieval Group, NASA/GSFC
1127 !!
1128 !! !PROGRAMMER:
1129 !!
1130 !! Nandana Amarasinghe (SSAI)
1131 !! Climate and Radiation Branch
1132 !! NASA Goddard Space Flight Center
1133 !! Greenbelt, Maryland 20771, U.S.A.
1134 !!
1135 !! Gala Wind (SSAI)
1136 !! Climate and Radiation Branch
1137 !! NASA Goddard Space Flight Center
1138 !! Greenbelt, Maryland 20771, U.S.A.
1139 !!
1140 !!*******************************************************************************
1141 !! !END
1142 ! use libraryarrays
1143 !
1144 ! IMPLICIT NONE
1145 !
1146 ! REAL,INTENT(IN):: solarAng, viewAng, azmAng
1147 ! REAL,INTENT(IN):: solarMuArray(:), viewMuArray(:), relAzmArray(:)
1148 ! REAL,INTENT(IN):: multiScatReflLUT_water(:,:,:,:,:,:)
1149 ! REAL,INTENT(IN):: multiScatReflLUT_ice(:,:,:,:,:,:)
1150 !
1151 ! REAL,INTENT(OUT)::multiScatRefl_water(:,:,:)
1152 ! REAL,INTENT(OUT)::multiScatRefl_ice(:,:,:)
1153 ! INTEGER,INTENT(OUT)::ierror
1154 !
1155 !
1156 ! REAL*8::prod11,prod12,prod21,prod22,vol1,vol2,vol3,vol4,vol5,vol6,vol7,vol8,volume
1157 ! REAL*4::solarmu,viewmu
1158 ! REAL*8::dmu0,dmu,dmu0_1,dmu0_2,dmu_1,dmu_2,dazm,dazm_1,dazm_2
1159 ! INTEGER:: ishi,islow, ivhi,ivlow, iazhi,iazlow,n1,n2,n3
1160 !
1161 ! n1 = number_solarzenith !size(solarMuArray)
1162 ! n2 = number_sensorzenith !size(viewMuArray)
1163 ! n3 = number_relazimuth !size(relAzmArray)
1164 !
1165 !
1166 ! ierror = 0 !no error
1167 ! solarmu = solarAng !COS(solarAng * RAD_PER_DEG)
1168 ! viewmu = viewAng !COS(viewAng * RAD_PER_DEG)
1169 !
1170 ! CALL BisectionSearch(solarmu,solarMuArray(1:n1),n1,islow,ishi)
1171 !
1172 ! CALL BisectionSearch(viewmu,viewMuArray(1:n2),n2,ivlow,ivhi)
1173 !
1174 ! CALL BisectionSearch(azmAng,relAzmArray(1:n3),n3,iazlow,iazhi)
1175 !
1176 ! dmu = solarMuArray(ishi) - solarMuArray(islow)
1177 ! dmu0 = viewMuArray(ivhi) - viewMuArray(ivlow)
1178 ! dazm = relAzmArray(iazhi) - relAzmArray(iazlow)
1179 !
1180 ! dmu_1 = solarMu - solarMuArray(islow) ; dmu_2 = solarMuArray(ishi) - solarMu
1181 ! dmu0_1 = viewMu - viewMuArray(ivlow) ; dmu0_2 = viewMuArray(ivhi) - viewMu
1182 ! dazm_1 = azmAng - relAzmArray(iazlow); dazm_2 = relAzmArray(iazhi) - azmAng
1183 !
1184 ! volume = dmu0 * dmu * dazm
1185 ! IF(volume == 0.0)THEN
1186 ! ierror = 3 !invalid data
1187 ! RETURN
1188 ! ENDIF
1189 ! !print*,ivhi,ishi,iazlow
1190 ! !print*,dmu0_1,dmu0_2,dmu_1,dmu_2,dazm_1,dazm_2
1191 !
1192 ! prod22 = dmu_2 * dmu0_2
1193 ! prod12 = dmu_1 * dmu0_2
1194 ! prod11 = dmu_1 * dmu0_1
1195 ! prod21 = dmu_2 * dmu0_1
1196 !
1197 ! vol8 = abs(prod22 * dazm_2)
1198 ! vol7 = abs(prod22 * dazm_1)
1199 ! vol6 = abs(prod12 * dazm_1)
1200 ! vol5 = abs(prod12 * dazm_2)
1201 !
1202 ! vol1 = abs(prod11 * dazm_1)
1203 ! vol2 = abs(prod21 * dazm_1)
1204 ! vol3 = abs(prod21 * dazm_2)
1205 ! vol4 = abs(prod11 * dazm_2)
1206 ! !volume = vol1 + vol2 + vol3 + vol4 + vol5 + vol6 + vol7 + vol8
1207 !
1208 !
1209 ! multiScatRefl_water(:,:,:) = ( multiScatReflLUT_water(ivlow,islow,iazlow,:,:,:) * vol8 + &
1210 ! multiScatReflLUT_water(ivhi,islow,iazlow,:,:,:) * vol3 + &
1211 ! multiScatReflLUT_water(ivhi,islow,iazhi,:,:,:) * vol2 + &
1212 ! multiScatReflLUT_water(ivlow,islow,iazhi,:,:,:) * vol7 + &
1213 ! multiScatReflLUT_water(ivlow,ishi,iazlow,:,:,:) * vol5 + &
1214 ! multiScatReflLUT_water(ivhi,ishi,iazlow,:,:,:) * vol4 + &
1215 ! multiScatReflLUT_water(ivhi,ishi,iazhi,:,:,:) * vol1 + &
1216 ! multiScatReflLUT_water(ivlow,ishi,iazhi,:,:,:) * vol6 )/volume
1217 !
1218 ! multiScatRefl_ice(:,:,:) = ( multiScatReflLUT_ice(ivlow,islow,iazlow,:,:,:) * vol8 + &
1219 ! multiScatReflLUT_ice(ivhi,islow,iazlow,:,:,:) * vol3 + &
1220 ! multiScatReflLUT_ice(ivhi,islow,iazhi,:,:,:) * vol2 + &
1221 ! multiScatReflLUT_ice(ivlow,islow,iazhi,:,:,:) * vol7 + &
1222 ! multiScatReflLUT_ice(ivlow,ishi,iazlow,:,:,:) * vol5 + &
1223 ! multiScatReflLUT_ice(ivhi,ishi,iazlow,:,:,:) * vol4 + &
1224 ! multiScatReflLUT_ice(ivhi,ishi,iazhi,:,:,:) * vol1 + &
1225 ! multiScatReflLUT_ice(ivlow,ishi,iazhi,:,:,:) * vol6 )/volume
1226 !
1227 !
1228 !
1229 ! RETURN
1230 !
1231 !
1232 ! CONTAINS
1233 !
1234 ! SUBROUTINE BisectionSearch(theta,thetaArray,ntheta,iAngLow, iAngHi)
1235 !
1236 ! INTEGER,INTENT(IN)::ntheta
1237 ! REAL,INTENT(IN)::theta,thetaArray(1:ntheta)
1238 ! INTEGER,INTENT(OUT)::iAngLow, iAngHi
1239 !
1240 ! INTEGER:: iAng
1241 !
1242 !
1243 ! iAngLow = 1
1244 ! iAngHi = ntheta
1245 !
1246 ! if (ntheta == 1) then
1247 ! iAngLow = 1
1248 ! iAngHi = 1
1249 ! return
1250 ! endif
1251 !
1252 ! DO
1253 ! IF((iAngHi - iAngLow) == 1)EXIT
1254 ! iAng = (iAngHi + iAngLow)/2
1255 ! IF(thetaArray(iAng) < theta)THEN
1256 ! iAngLow = iAng
1257 ! ELSE
1258 ! iAngHi = iAng
1259 ! ENDIF
1260 ! ENDDO
1261 !
1262 ! if (iAngLow < 1) iAngLow = 1
1263 ! if (ianghi < 1) ianghi = 1
1264 ! if (iAngHi > ntheta) iAngHi = ntheta
1265 ! if (ianglow > ntheta) ianglow = ntheta
1266 !
1267 !
1268 ! RETURN
1269 !
1270 ! END SUBROUTINE BisectionSearch
1271 !
1272 ! END SUBROUTINE TriLinearInterpolationRefl
1273 
1274 
1275  SUBROUTINE getrefl(solarAng,viewAng,azmAng,in_scat, refl_water, refl_ice, interp_MS, interp_SS, ierror)
1276 !................................................................................
1277 ! !F90
1278 !
1279 ! !DESCRIPTION:
1280  !This subroutine outouts the Reflectance table (R(tau,lambda, re)) for a given sun-satellite
1281  !geometry
1282 ! !INPUT PARAMETERS:
1283  ! solarAng : solar zenith in degrees (scalar)
1284  ! viewAng : viewing angle in degrees (scalar)
1285  ! azmAng : relative azimuth between sat and solar planes
1286 !
1287 !
1288 ! !OUTPUT PARAMETERS:
1289  ! refl : reflectance as a function of tau,lambda,re : 3-D Array
1290  ! DIM-1 = ntau, DIM-2 = nlambda, DIM-3 = nre
1291  ! ierror : integer value
1292  ! ierror = 0 success,
1293  ! = 1 I/O error
1294  ! = 2 memory allocation error
1295  ! = 3 invalid data
1296 !................................................................................
1297 ! !Revision History:
1298 ! Revision 1.0 2009/06/03 GW
1299 ! Initial integration of the standalone version. Changed from the standalone
1300 ! delivered version to not apply the albedo inside this routine. That is done in a
1301 ! different place later in corescience.
1302 !
1303 !................................................................................
1304 ! !Team-unique Header:
1305 ! Cloud Retrieval Group, NASA/GSFC
1306 !
1307 ! !PROGRAMMER:
1308 !
1309 ! Nandana Amarasinghe (SSAI)
1310 ! Climate and Radiation Branch
1311 ! NASA Goddard Space Flight Center
1312 ! Greenbelt, Maryland 20771, U.S.A.
1313 !
1314 ! Gala Wind (SSAI)
1315 ! Climate and Radiation Branch
1316 ! NASA Goddard Space Flight Center
1317 ! Greenbelt, Maryland 20771, U.S.A.
1318 !
1319 !*******************************************************************************
1320 ! !END
1321 
1322  use libraryarrays
1323  use generalauxtype
1324  ! WDR to get the scan line
1325  use ch_xfr, only : c2_scan
1326  !
1327  REAL, INTENT(IN) :: solarAng,viewAng,azmAng, in_scat
1328  REAL,INTENT(OUT) :: refl_water(:,:,:), refl_ice(:,:,:)
1329  INTEGER,INTENT(OUT)::ierror
1330  logical, intent(in) :: interp_MS, interp_SS
1331 
1332  INTEGER::allocateStatus
1333  REAL :: scat_angle
1334  INTEGER::do_mgr
1335 
1336  !WDR TEMP set some loop controls for refl table output
1337  !INTEGER::ksen,klam,ktau,kre
1338 
1339  scat_angle = in_scat ! scatAngle(solarAng,viewAng,azmAng)
1340 
1341  if (interp_ms) then
1342 
1343  !interpolate multi scattering part of the reflection
1344  multiscatrefl_water_lamb = 0
1345  multiscatrefl_ice_lamb = 0
1346  call mng_ms_get_lib( solarang, viewang, azmang, 0, 0, c2_scan, &
1347  multiscatrefl_water_lamb, multiscatrefl_ice_lamb, ierror )
1348  if( ierror /= 0 ) return
1349  endif
1350 
1351  if (interp_ss) then
1352 
1353  !get the phase function values
1355  phase_fun_norm_constant_water, phasefunvals_water,ierror)
1357  phase_fun_norm_constant_ice, phasefunvals_ice,ierror)
1358 
1359 
1360 ! IF(ierror /= 0)RETURN !if there is any error no point of continuing
1361 
1362  !calculate the single scattering part of the reflection
1363  CALL calcsingscatpartofreflectance(phasefunvals_water, phasefunvals_ice, &
1364  viewang, solarang, ssrefl_water, ssrefl_ice, ierror)
1365 
1366  endif
1367 
1368 
1369 ! IF(ierror /= 0)RETURN !if there is any error no point of continuing
1370 
1371 
1372  refl_water = ssrefl_water + multiscatrefl_water_lamb !total reflection with asurf = 0.0
1373  refl_ice = ssrefl_ice + multiscatrefl_ice_lamb !total reflection with asurf = 0.0
1374 
1375  RETURN
1376  END SUBROUTINE getrefl
1377 
1378 
1379 
1380  !--------------------------SUBROUTINE GetReflForGivenWindSpeedWater-----------------------------
1381  !This subroutine outputs the Reflectance table (R(tau,lambda, re)) for a given sun-satellite
1382  !geometry and a wind speed
1383  !INPUTS :
1384  ! solarAng : solar zenith in degrees (scalar)
1385  ! viewAng : viewing angle in degrees (scalar)
1386  ! azmAng : relative azimuth between sat and solar planes
1387  ! wspeed : wind speed m/s
1388  !OUTPUTS :
1389  ! reflAsurf : reflectance as a function of tau,lambda,re : 3-D Array
1390  ! DIM-1 = ntau, DIM-2 = nlambda, DIM-3 = nre
1391  ! ierror : integer value
1392  ! ierror = 0 success,
1393  ! = 1 I/O error
1394  ! = 2 memory allocation error
1395  ! = 3 invalid data
1396  !ROUTINES CAlLED :
1397  ! ReadPhaseFunctionDataHDF
1398  ! ReadWaterLibraryHDF
1399  ! getAllPhaseFunctionValues
1400  ! TriLinearInterpolationReflSearch
1401  ! CalcSingScatPartOfReflectance`
1402  !-----------------------------------------------------------------------------------------------
1403  SUBROUTINE getreflforgivenwindspeed(solarAng,viewAng,azmAng, in_scat, cos_scat, wspeed, reflAsurf_water, &
1404  reflAsurf_ice, wind_speed_only, interp_MS, interp_SS, ierror)
1406 
1407  use libraryarrays
1409  use ch_xfr, only : c2_scan
1410  !
1411  REAL, INTENT(IN) :: solarAng,viewAng,azmAng,wspeed, in_scat, cos_scat
1412  REAL,INTENT(OUT) :: reflAsurf_water(:,:,:), reflAsurf_ice(:,:,:)
1413  INTEGER,INTENT(OUT)::ierror
1414  logical, intent(in) :: wind_speed_only, interp_MS, interp_SS
1415 
1416  INTEGER::allocateStatus,iuhi,iulow
1417  REAL :: scat_angle,deltau
1418  INTEGER::do_mgr
1419 
1420  scat_angle = in_scat !scatAngle(solarAng,viewAng,azmAng)
1421 
1422  if (interp_ss) then
1423 
1424  !get the phase function values
1425  call get_aero_params(cos_scat, aerosol_asym)
1426 
1428  phase_fun_norm_constant_water, phasefunvals_water,ierror)
1429 
1431  phase_fun_norm_constant_ice, phasefunvals_ice,ierror)
1432 
1433  !calculate the single scattering part of the reflection
1434  call single_scattering_calcs_ocean(phasefunvals_water, phasefunvals_ice, singlescattering_water, singlescattering_ice, &
1435  rlphase, aerophase, rayleigh_tau, aerosol_tau, aerosol_ssa, &
1436  solarang,viewang, ssrefl_water, ssrefl_ice)
1437 
1438  endif
1439 
1440 
1441 
1442  if (interp_ms) then
1443 
1444 ! ! WDR wire up the sfc typ 1 for 3 m/s
1445  multiscatrefl_water_cm(:,:,:,1) = 0
1446  multiscatrefl_ice_cm(:,:,:,1) = 0
1447  call mng_ms_get_lib( solarang, viewang, azmang, 1, 0, c2_scan, &
1448  multiscatrefl_water_cm(:,:,:,1), multiscatrefl_ice_cm(:,:,:,1), &
1449  ierror )
1450  if( ierror /= 0 ) return
1451 
1452 !
1453 ! ! WDR wire up the sfc typ 1 for 7 m/s
1454  multiscatrefl_water_cm(:,:,:,2) = 0
1455  multiscatrefl_ice_cm(:,:,:,2) = 0
1456  call mng_ms_get_lib( solarang, viewang, azmang, 2, 0, c2_scan, &
1457  multiscatrefl_water_cm(:,:,:,2), multiscatrefl_ice_cm(:,:,:,2), &
1458  ierror )
1459  if( ierror /= 0 ) return
1460 
1461 !
1462 ! ! WDR wire up the sfc typ 1 for 15 m/s
1463  multiscatrefl_water_cm(:,:,:,3) = 0
1464  multiscatrefl_ice_cm(:,:,:,3) = 0
1465  call mng_ms_get_lib( solarang, viewang, azmang, 3, 0, c2_scan, &
1466  multiscatrefl_water_cm(:,:,:,3), multiscatrefl_ice_cm(:,:,:,3), &
1467  ierror )
1468  if( ierror /= 0 ) return
1469  endif
1470 
1471 
1472  int_reflectance_water_wspeed(:,:,:,1) = multiscatrefl_water_cm(:,:,:,1) + ssrefl_water(:,:,:)
1473  int_reflectance_water_wspeed(:,:,:,2) = multiscatrefl_water_cm(:,:,:,2) + ssrefl_water(:,:,:)
1474  int_reflectance_water_wspeed(:,:,:,3) = multiscatrefl_water_cm(:,:,:,3) + ssrefl_water(:,:,:)
1475 
1476  int_reflectance_ice_wspeed(:,:,:,1) = multiscatrefl_ice_cm(:,:,:,1) + ssrefl_ice(:,:,:)
1477  int_reflectance_ice_wspeed(:,:,:,2) = multiscatrefl_ice_cm(:,:,:,2) + ssrefl_ice(:,:,:)
1478  int_reflectance_ice_wspeed(:,:,:,3) = multiscatrefl_ice_cm(:,:,:,3) + ssrefl_ice(:,:,:)
1479 
1480 
1481 
1482  call interpolate_wind_speed(wspeed, int_reflectance_water_wspeed, reflasurf_water)
1483  call interpolate_wind_speed(wspeed, int_reflectance_ice_wspeed, reflasurf_ice)
1484 
1485  END SUBROUTINE getreflforgivenwindspeed
1486 
1487 
1488 
1489  subroutine interpolate_wind_speed(wspeed, data_in, data_out)
1491  use ch_xfr, only : c2_scan
1492  real, intent(in) :: wspeed
1493  real, dimension(:,:,:,:), intent(in) :: data_in
1494  real, dimension(:,:,:), intent(inout) :: data_out
1495 
1496  real :: deltau
1497  integer :: iuhi, iulow
1498 
1499  IF( (wspeed - wind_speed(1)) <= 0.2)THEN
1500 
1501  data_out(:,:,:) = data_in(:,:,:,1)
1502 
1503 
1504  ELSEIF( (wspeed - wind_speed(3)) >= -0.2)THEN
1505 
1506  data_out(:,:,:) = data_in(:,:,:,3)
1507 
1508 
1509  ELSEIF(abs(wspeed - wind_speed(2)) <= 0.2)THEN
1510 
1511 
1512  data_out(:,:,:) = data_in(:,:,:,2)
1513 
1514  ELSE
1515 
1516  IF(wspeed - wind_speed(2) < 0.0)THEN
1517  deltau = wind_speed(2) - wind_speed(1)
1518  iuhi = 2
1519  iulow = 1
1520 
1521  data_out(:,:,:) = ( data_in(:,:,:,1) * (wind_speed(iuhi)- wspeed) + &
1522  data_in(:,:,:,2) * (wspeed - wind_speed(iulow) ))/deltau
1523 
1524 
1525  ELSE
1526 
1527  deltau = wind_speed(3) - wind_speed(2)
1528  iuhi = 3
1529  iulow = 2
1530 
1531  data_out(:,:,:) = ( data_in(:,:,:,2) * (wind_speed(iuhi)- wspeed) + &
1532  data_in(:,:,:,3) * (wspeed - wind_speed(iulow) ))/deltau
1533 
1534 
1535 
1536  ENDIF
1537 
1538 
1539  ENDIF
1540 
1541 
1542  end subroutine interpolate_wind_speed
1543 
1544 
1545 
1546  !--------------------------SUBROUTINE GetReflForGivenWindSpeedIce-------------------------------
1547  !This subroutine outputs the Reflectance table (R(tau,lambda, re)) for a given sun-satellite
1548  !geometry and a wind speed(m/s)
1549  !INPUTS :
1550  ! solarAng : solar zenith in degrees (scalar)
1551  ! viewAng : viewing angle in degrees (scalar)
1552  ! azmAng : relative azimuth between sat and solar planes
1553  ! wspeed : wind speed m/s
1554  !OUTPUTS :
1555  ! reflAsurf : reflectance as a function of tau,lambda,re : 3-D Array
1556  ! DIM-1 = ntau, DIM-2 = nlambda, DIM-3 = nre
1557  ! ierror : integer value
1558  ! ierror = 0 success,
1559  ! = 1 I/O error
1560  ! = 2 memory allocation error
1561  ! = 3 invalid data
1562  !ROUTINES CAlLED :
1563  ! ReadPhaseFunctionDataHDF
1564  ! ReadIceLibraryHDF
1565  ! getAllPhaseFunctionValues
1566  ! TriLinearInterpolationReflSearch
1567  ! CalcSingScatPartOfReflectance
1568  !-----------------------------------------------------------------------------------------------
1569  SUBROUTINE getsdevreflforgivenwindspeed(solarAng,viewAng,azmAng,wspeed, &
1570  inrefl_ws_water, reflAsurf_water, inrefl_ws_ice, reflAsurf_ice, &
1571  wind_speed_only, ierror)
1573  use libraryarrays
1574  use ch_xfr, only : c2_scan
1575  !
1576  REAL, INTENT(IN) :: solarAng,viewAng,azmAng,wspeed
1577  real, intent(inout) :: inrefl_ws_water(:,:,:,:), inrefl_ws_ice(:,:,:,:)
1578  REAL,INTENT(OUT) :: reflAsurf_water(:,:,:), reflAsurf_ice(:,:,:)
1579  INTEGER,INTENT(OUT) :: ierror
1580  logical, intent(in) :: wind_speed_only
1581 
1582  INTEGER :: allocateStatus,iuhi,iulow
1583  REAL :: deltau
1584  integer :: num_radii
1585  INTEGER::do_mgr
1586 
1587  if (.not. wind_speed_only) then
1588 
1589 !
1590 ! ! WDR wire up the sfc typ 1 for 3 m/s
1591  inrefl_ws_water(:,:,:,1) = 0
1592  inrefl_ws_ice(:,:,:,1) = 0
1593  call mng_ms_get_lib( solarang, viewang, azmang, 1, 1, c2_scan, &
1594  inrefl_ws_water(:,:,:,1), inrefl_ws_ice(:,:,:,1), ierror )
1595  if( ierror /= 0 ) return
1596 
1597 !
1598 ! ! WDR wire up the sfc typ 2 for 7 m/s
1599  inrefl_ws_water(:,:,:,2) = 0
1600  inrefl_ws_ice(:,:,:,2) = 0
1601  call mng_ms_get_lib( solarang, viewang, azmang, 2, 1, c2_scan, &
1602  inrefl_ws_water(:,:,:,2), inrefl_ws_ice(:,:,:,2), ierror )
1603  if( ierror /= 0 ) return
1604 
1605 !
1606 ! ! WDR wire up the sfc typ 3 for 15 m/s
1607  inrefl_ws_water(:,:,:,3) = 0
1608  inrefl_ws_ice(:,:,:,3) = 0
1609  call mng_ms_get_lib( solarang, viewang, azmang, 3, 1, c2_scan, &
1610  inrefl_ws_water(:,:,:,3), inrefl_ws_ice(:,:,:,3), ierror )
1611  if( ierror /= 0 ) return
1612 
1613  !deallocate( wtrsd_int, icesd_int )
1614 
1615  endif
1616 
1617 
1618  call interpolate_wind_speed(wspeed, inrefl_ws_water, reflasurf_water)
1619  call interpolate_wind_speed(wspeed, inrefl_ws_ice, reflasurf_ice)
1620 
1621 
1622  RETURN
1623 
1624  END SUBROUTINE getsdevreflforgivenwindspeed
1625  !-----------------------------------------------------------------------------------------------
1626  SUBROUTINE getsdevrefllamb(solarAng,viewAng,azmAng, reflAsurf_water, &
1627  reflAsurf_ice, ierror)
1629 
1630  use libraryarrays
1631  ! WDR access the scan line
1632  use ch_xfr, only : c2_scan
1633  !
1634  REAL, INTENT(IN) :: solarAng,viewAng,azmAng
1635  REAL,INTENT(OUT) :: reflAsurf_water(:,:,:), reflAsurf_ice(:,:,:)
1636  INTEGER,INTENT(OUT) :: ierror
1637 
1638  INTEGER :: allocateStatus,iuhi,iulow
1639  REAL :: deltau
1640  integer :: num_radii
1641  INTEGER::do_mgr
1642 
1643  reflasurf_water = 0
1644  reflasurf_ice = 0
1645  call mng_ms_get_lib( solarang, viewang, azmang, 0, 1, c2_scan, &
1646  reflasurf_water, reflasurf_ice, ierror )
1647  if( ierror /= 0 ) return
1648 
1649 ! deallocate( wtrsd_int, icesd_int )
1650 
1651  END SUBROUTINE getsdevrefllamb
1652 
1653  end module interpolate_libraries
Definition: ch_xfr.f90:1
subroutine interpolatefluxes(solarOrViewAng, solarAngMuArray, inFluxArray, outFluxArray, iAngHi, iAngLow, dtheta, ierror)
integer number_wavelengths
real(single), dimension(:,:,:,:,:), allocatable surface_emissivity_water
real(single), dimension(:,:,:), allocatable int_surface_emissivity_water_sdev
real(single), dimension(:,:,:,:,:), allocatable cloud_emissivity_ice
real(single), dimension(:,:), allocatable phase_fun_norm_constant_ice
integer ocis_id
Definition: ch_xfr.f90:51
real(single), dimension(:,:,:,:,:), allocatable surface_emissivity_water_sdev
real(single), dimension(:,:,:,:,:), allocatable surface_emissivity_ice
subroutine getrefl(solarAng, viewAng, azmAng, in_scat, refl_water, refl_ice, interp_MS, interp_SS, ierror)
real(single), dimension(:,:,:,:), allocatable int_surface_emis_ice_wspeed
integer number_phase_angles_water
real(single), dimension(:,:,:,:), allocatable flux_down_ice_sensor
real(single), dimension(:,:), allocatable phase_fun_norm_constant_water
real function, public scatangle(solarAng, viewAng, relAzm)
real(single), dimension(:,:,:), allocatable int_reflectance_ice
subroutine, public libraryinterpolate(local_solarzenith, local_sensorzenith, local_relativeazimuth, local_scatangle, local_wind_speed, wind_speed_only, interp_MS, interp_SS, debug, status, i, j)
real(single), dimension(:,:,:,:), allocatable int_cloud_emis_water_wspeed
real(single), dimension(:,:,:), allocatable int_fluxdownwater_solar
real(single), dimension(:,:,:,:), allocatable flux_down_water_solar
character *15 platform_name
#define real
Definition: DbAlgOcean.cpp:26
real(single), dimension(:,:), allocatable singlescattering_water
real(single), dimension(:,:,:,:,:), allocatable cloud_emissivity_water_sdev
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_ice
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_water
real(single), dimension(:), allocatable aerosol_ssa
real(single), dimension(:,:,:,:), allocatable flux_down_ice_solar
real(single), dimension(:), allocatable library_fluxsolarzenith
real(single), dimension(:), allocatable library_taus
real(single), dimension(:,:,:), allocatable int_reflectance_ice_sdev
real(single), dimension(:,:,:,:), allocatable int_cloud_emis_water_sdev_wspeed
real(single), dimension(:,:,:), allocatable int_reflectance_water_sdev
real(single), dimension(:,:,:,:), allocatable int_cloud_emis_ice_sdev_wspeed
real(single), dimension(:,:,:,:,:), allocatable surface_emissivity_ice_sdev
subroutine, public interpolate_wind_speed(wspeed, data_in, data_out)
real(single), dimension(:,:,:), allocatable int_surface_emissivity_ice
real(single), dimension(:,:), allocatable singlescattering_ice
real(single), dimension(:,:), allocatable extinction_ice
real(single), dimension(:,:,:), allocatable int_fluxdownice_solar
real(single), dimension(:,:,:,:), allocatable int_refl_ice_sdev_wspeed
real(single), dimension(:,:), allocatable extinction_water
subroutine getphasefunctionvalues(scatAngle, scatAngleArray, num_angles, phaseFunArray, phaseNormConst, phaseFunVals, ierror)
real(single), dimension(:), allocatable library_sensor_zenith
real(single), dimension(:,:,:), allocatable int_surface_emissivity_water
real(single), dimension(:,:,:), allocatable phase_funcs_ice
real(single), dimension(:), allocatable aerosol_asym
real(single), dimension(:,:,:,:), allocatable flux_up_water_sensor
real(single), dimension(:,:,:,:), allocatable int_surface_emis_water_wspeed
real(single), dimension(:,:,:,:), allocatable flux_down_water_sensor
real(single), dimension(:,:,:), allocatable int_surface_emissivity_ice_sdev
integer c2_sensor_id
Definition: ch_xfr.f90:50
real(single), dimension(:), allocatable rayleigh_tau
subroutine get_aero_params(cos_scatAngle, aeroG)
integer number_phase_angles_ice
subroutine getsdevreflforgivenwindspeed(solarAng, viewAng, azmAng, wspeed, inrefl_ws_water, reflAsurf_water, inrefl_ws_ice, reflAsurf_ice, wind_speed_only, ierror)
subroutine getsdevrefllamb(solarAng, viewAng, azmAng, reflAsurf_water, reflAsurf_ice, ierror)
real(single), dimension(:,:,:,:), allocatable int_surface_emis_water_sdev_wspeed
real(single), dimension(:,:,:), allocatable int_fluxupice_sensor
real(single), dimension(:,:,:), allocatable int_reflectance_water
integer number_waterradii
subroutine setup_emissivity_flux(angle, angle_array, idx_hi, idx_lo, dtheta)
subroutine getreflforgivenwindspeed(solarAng, viewAng, azmAng, in_scat, cos_scat, wspeed, reflAsurf_water, reflAsurf_ice, wind_speed_only, interp_MS, interp_SS, ierror)
real(single), dimension(:,:,:,:,:), allocatable cloud_emissivity_ice_sdev
subroutine single_scattering_calcs_ocean(phaseFunVals_liq, phaseFunVals_ice, ssAlbedo_liq, ssAlbedo_ice, RLphase, aeroPhase, RLTau, aeroTau, aeroOmega, theta, theta0, ssRefl_liq, ssRefl_ice)
integer c2_scan
Definition: ch_xfr.f90:46
real(single), dimension(:,:,:), allocatable int_fluxdownwater_sensor
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_ice_sdev
real(single), dimension(:,:), allocatable truncation_factor_ice
integer number_iceradii
integer oci_id
Definition: ch_xfr.f90:52
real(single), dimension(:,:,:,:), allocatable int_cloud_emis_ice_wspeed
real(single), dimension(:,:,:,:,:), allocatable cloud_emissivity_water
real(single), dimension(:,:,:), allocatable phase_funcs_water
real(single), dimension(:,:,:), allocatable int_fluxdownice_sensor
real(single), dimension(:,:,:,:), allocatable int_surface_emis_ice_sdev_wspeed
real(single), dimension(:,:,:,:), allocatable int_reflectance_water_wspeed
real(single), dimension(:), allocatable phase_angles_water
real(single), dimension(:,:,:,:), allocatable flux_up_ice_sensor
real(single), dimension(:,:,:,:), allocatable int_refl_water_sdev_wspeed
real(single), dimension(:), allocatable phase_angles_ice
#define abs(a)
Definition: misc.h:90
real(single), dimension(:,:,:,:), allocatable int_reflectance_ice_wspeed
real(single), dimension(:), allocatable library_fluxsensorzenith
real(single), dimension(:), allocatable aerosol_tau
real(single), dimension(:,:), allocatable truncation_factor_water
real(single), dimension(:,:,:), allocatable int_cloud_emissivity_water_sdev
real(single), dimension(:,:,:), allocatable int_fluxupwater_sensor
integer number_taus