190 Atmosphere , & ! Input
192 AtmOptics , & ! Input
193 SfcOptics , & ! Input
194 GeometryInfo, & ! Input
195 SensorIndex , & ! Input
196 ChannelIndex, & ! Input
197 RTSolution , & ! Output
200 result( error_status )
207 INTEGER,
INTENT(IN) :: sensorindex
208 INTEGER,
INTENT(IN) :: channelindex
210 INTEGER,
INTENT(OUT) :: nz
211 TYPE(
rtv_type),
INTENT(IN OUT) :: rtv
213 INTEGER :: error_status
215 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Input' 217 CHARACTER(256) :: message
218 INTEGER :: no, na, nt
221 REAL(fp) :: user_emissivity, direct_reflectivity
226 geometryinfo%Cosine_Sensor_Zenith =
one / geometryinfo%Secant_Sensor_Zenith
228 rtv%n_Added_Layers = atmosphere%n_Added_Layers
229 rtv%n_Layers = atmosphere%n_Layers
230 rtv%n_Streams = rtsolution%n_Full_Streams/2
235 no = rtsolution%n_Layers
236 na = rtv%n_Added_Layers
240 rtsolution%Layer_Optical_Depth(1:no) = atmoptics%Optical_Depth(na+1:nt)
243 rtv%Cosmic_Background_Radiance =
sc(sensorindex)%Cosmic_Background_Radiance(channelindex)
244 rtv%Sensor_Type =
sc(sensorindex)%Sensor_Type
245 rtv%Is_Solar_Channel = spccoeff_issolar(
sc(sensorindex), channelindex=channelindex )
246 rtv%COS_SUN = cos(geometryinfo%Source_Zenith_Radian)
247 rtv%Solar_Irradiance =
sc(sensorindex)%Solar_Irradiance(channelindex) * geometryinfo%AU_ratio2
251 rtv%Diffuse_Surface = .false.
257 rtv%Scattering_RT = .false.
265 rtv%Scattering_RT = .true.
279 DO i = 1, rtv%n_Streams
280 IF( abs( rtv%COS_Angle(i) - geometryinfo%Cosine_Sensor_Zenith ) <
angle_threshold )
THEN 281 sfcoptics%Index_Sat_Ang = i
282 rtv%n_Angles = rtv%n_Streams
292 rtv%n_Angles = rtv%n_Streams + 1
293 sfcoptics%Index_Sat_Ang = rtv%n_Angles
294 rtv%COS_Angle( rtv%n_Angles ) = geometryinfo%Cosine_Sensor_Zenith
295 rtv%COS_Weight( rtv%n_Angles ) =
zero 305 ELSE IF( rtv%Diffuse_Surface)
THEN 313 rtv%COS_Weight( 1 ) =
one 315 rtv%n_Angles = rtv%n_Streams + 1
316 sfcoptics%Index_Sat_Ang = rtv%n_Angles
317 rtv%COS_Angle( rtv%n_Angles ) = geometryinfo%Cosine_Sensor_Zenith
318 rtv%COS_Weight( rtv%n_Angles ) =
zero 327 rtv%COS_Angle( rtv%n_Angles ) = geometryinfo%Cosine_Sensor_Zenith
328 rtv%COS_Weight( rtv%n_Angles ) =
one 329 sfcoptics%Index_Sat_Ang = rtv%n_Angles
331 END IF determine_stream_angles
343 sfcoptics%n_Stokes = 1
348 sfcoptics%n_Angles = rtv%n_Angles
352 sfcoptics%Weight(i) = rtv%COS_Angle(i) * rtv%COS_Weight(i)
354 factor = factor + sfcoptics%Weight(i)
357 IF( factor >
zero) sfcoptics%Weight(1:nz) = sfcoptics%Weight(1:nz) / factor
362 IF ( sfcoptics%Compute )
THEN 370 IF ( error_status /=
success )
THEN 371 WRITE( message,
'("Error computing SfcOptics for ",a," channel ",i0)' ) &
372 trim(
sc(sensorindex)%Sensor_Id), &
373 sc(sensorindex)%Sensor_Channel(channelindex)
380 IF( rtv%Scattering_RT )
THEN 382 sfcoptics%Reflectivity =
zero 383 user_emissivity = sfcoptics%Emissivity(1,1)
384 sfcoptics%Emissivity(1,1) =
zero 385 direct_reflectivity = sfcoptics%Direct_Reflectivity(1,1)/
pi 386 sfcoptics%Emissivity(1:nz,1) = user_emissivity
388 sfcoptics%Direct_Reflectivity(1:nz,1) = direct_reflectivity
389 IF( rtv%Diffuse_Surface)
THEN 391 sfcoptics%Reflectivity(1:nz, 1, i, 1) = (
one-sfcoptics%Emissivity(i,1))*sfcoptics%Weight(i)
395 sfcoptics%Reflectivity(i, 1, i, 1) = (
one-sfcoptics%Emissivity(i,1))
399 user_emissivity = sfcoptics%Emissivity(1,1)
400 sfcoptics%Emissivity( sfcoptics%Index_Sat_Ang,1 ) = user_emissivity
401 sfcoptics%Reflectivity(1,1,1,1) =
one - user_emissivity
409 rtsolution%Surface_Emissivity = sfcoptics%Emissivity( sfcoptics%Index_Sat_Ang, 1 )
410 rtsolution%Surface_Reflectivity = sfcoptics%Reflectivity( sfcoptics%Index_Sat_Ang, 1, sfcoptics%Index_Sat_Ang, 1 )
415 IF( rtv%mth_Azi == 0 )
THEN 416 DO k = 1, atmosphere%n_Layers
420 atmosphere%Temperature(k), &
421 rtv%Planck_Atmosphere(k) )
427 sfcoptics%Surface_Temperature, &
430 rtv%Planck_Atmosphere =
zero 431 rtv%Planck_Surface =
zero 622 Atmosphere , & ! FWD Input
623 Surface , & ! FWD Input
624 AtmOptics , & ! FWD Input
625 SfcOptics , & ! FWD Input
626 Atmosphere_TL , & ! TL Input
627 Surface_TL , & ! TL Input
628 AtmOptics_TL , & ! TL Input
629 SfcOptics_TL , & ! TL Input/Output
630 GeometryInfo , & ! Input
631 SensorIndex , & ! Input
632 ChannelIndex , & ! Input
633 RTSolution_TL , & ! TL Output
635 User_Emissivity_TL , & ! Output
636 Direct_Reflectivity_TL, & ! Output
637 Planck_Surface_TL , & ! Output
638 Planck_Atmosphere_TL , & ! Output
642 result( error_status )
653 INTEGER,
INTENT(IN) :: sensorindex
654 INTEGER,
INTENT(IN) :: channelindex
656 INTEGER,
INTENT(OUT) :: nz
657 REAL(fp),
INTENT(OUT) :: user_emissivity_tl
658 REAL(fp),
INTENT(OUT) :: direct_reflectivity_tl
659 REAL(fp),
INTENT(OUT) :: planck_surface_tl
660 REAL(fp),
INTENT(OUT) :: planck_atmosphere_tl(0:)
661 REAL(fp),
INTENT(OUT) :: pff_tl(:,:,:)
662 REAL(fp),
INTENT(OUT) :: pbb_tl(:,:,:)
665 INTEGER :: error_status
667 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Input_TL' 669 CHARACTER(256) :: message
671 INTEGER :: no, na, nt
680 sfcoptics_tl%Index_Sat_Ang = sfcoptics%Index_Sat_Ang
684 no = rtsolution_tl%n_Layers
685 na = rtv%n_Added_Layers
689 rtsolution_tl%Layer_Optical_Depth(1:no) = atmoptics_tl%Optical_Depth(na+1:nt)
695 IF ( rtv%Scattering_RT )
THEN 711 sfcoptics_tl%n_Stokes = 1
717 sfcoptics_tl%n_Angles = sfcoptics%n_Angles
718 sfcoptics_tl%Angle = sfcoptics%Angle
719 sfcoptics_tl%Weight = sfcoptics%Weight
726 IF ( sfcoptics%Compute )
THEN 736 IF ( error_status /=
success )
THEN 737 WRITE( message,
'("Error computing SfcOptics_TL for ",a," channel ",i0)' ) &
738 trim(
sc(sensorindex)%Sensor_Id), &
739 sc(sensorindex)%Sensor_Channel(channelindex)
744 IF( rtv%Scattering_RT )
THEN 746 sfcoptics_tl%Reflectivity =
zero 747 user_emissivity_tl = sfcoptics_tl%Emissivity(1,1)
748 sfcoptics_tl%Emissivity(1,1) =
zero 749 direct_reflectivity_tl = sfcoptics_tl%Direct_Reflectivity(1,1)/
pi 750 sfcoptics_tl%Emissivity(1:nz,1) = user_emissivity_tl
752 sfcoptics_tl%Direct_Reflectivity(1:nz,1) = direct_reflectivity_tl
753 IF( rtv%Diffuse_Surface)
THEN 755 sfcoptics_tl%Reflectivity(1:nz, 1, i, 1) = -sfcoptics_tl%Emissivity(i,1)*sfcoptics%Weight(i)
759 sfcoptics_tl%Reflectivity(i, 1, i, 1) = -sfcoptics_tl%Emissivity(i,1)
763 user_emissivity_tl = sfcoptics_tl%Emissivity(1,1)
764 sfcoptics_tl%Emissivity( sfcoptics%Index_Sat_Ang,1 ) = user_emissivity_tl
765 sfcoptics_tl%Reflectivity(1,1,1,1) = - user_emissivity_tl
772 rtsolution_tl%Surface_Emissivity = sfcoptics_tl%Emissivity( sfcoptics_tl%Index_Sat_Ang, 1 )
773 rtsolution_tl%Surface_Reflectivity = sfcoptics_tl%Reflectivity( sfcoptics_tl%Index_Sat_Ang, 1, &
774 sfcoptics_tl%Index_Sat_Ang, 1 )
779 IF( rtv%mth_Azi == 0 )
THEN 781 DO k = 1, atmosphere%n_Layers
785 atmosphere%Temperature(k) , &
786 atmosphere_tl%Temperature(k), &
787 planck_atmosphere_tl(k) )
793 sfcoptics%Surface_Temperature , &
794 sfcoptics_tl%Surface_Temperature, &
797 planck_atmosphere_tl =
zero 798 planck_surface_tl =
zero 932 SfcOptics , & ! FWD Input
933 RTSolution , & ! FWD Input
934 GeometryInfo , & ! Input
935 SensorIndex , & ! Input
936 ChannelIndex , & ! Input
937 RTSolution_AD , & ! AD Output/Input
938 SfcOptics_AD , & ! AD Output
939 Planck_Surface_AD , & ! AD Output
940 Planck_Atmosphere_AD , & ! AD Output
941 Radiance_AD , & ! AD Output
944 result( error_status )
949 INTEGER,
INTENT(IN) :: sensorindex
950 INTEGER,
INTENT(IN) :: channelindex
953 REAL(fp),
INTENT(OUT) :: planck_surface_ad
954 REAL(fp),
INTENT(OUT) :: planck_atmosphere_ad(0:)
955 REAL(fp),
INTENT(OUT) :: radiance_ad
956 INTEGER,
INTENT(OUT) :: nz
959 INTEGER :: error_status
961 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Input_AD' 970 sfcoptics_ad%Index_Sat_Ang = sfcoptics%Index_Sat_Ang
973 planck_surface_ad =
zero 974 planck_atmosphere_ad =
zero 980 IF ( spccoeff_isinfraredsensor(
sc(sensorindex) ) .OR. &
981 spccoeff_ismicrowavesensor(
sc(sensorindex) ) )
THEN 982 IF( rtv%mth_Azi == 0 )
THEN 986 rtsolution%Radiance , &
987 rtsolution_ad%Brightness_Temperature, &
989 rtsolution_ad%Brightness_Temperature =
zero 994 radiance_ad = radiance_ad + rtsolution_ad%Radiance * &
995 cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1092 Atmosphere , & ! Input
1093 SfcOptics , & ! Input
1094 GeometryInfo , & ! Input
1095 SensorIndex , & ! Input
1096 ChannelIndex , & ! Input
1099 result( error_status )
1104 INTEGER ,
INTENT(IN) :: sensorindex
1105 INTEGER ,
INTENT(IN) :: channelindex
1110 INTEGER :: error_status
1112 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Output' 1114 INTEGER :: no, na, nt
1115 REAL(fp) :: radiance
1120 IF( rtv%Scattering_RT )
THEN 1123 IF ( rtv%aircraft%rt )
THEN 1124 radiance = rtv%s_Level_Rad_UP(sfcoptics%Index_Sat_Ang, rtv%aircraft%idx)
1126 radiance = rtv%s_Level_Rad_UP(sfcoptics%Index_Sat_Ang, 0)
1133 IF ( rtv%aircraft%rt )
THEN 1134 radiance = rtv%e_Level_Rad_UP(rtv%aircraft%idx)
1136 radiance = rtv%e_Level_Rad_UP(0)
1140 rtsolution%Up_Radiance = rtv%Up_Radiance
1141 rtsolution%Down_Radiance = rtv%e_Level_Rad_DOWN(atmosphere%n_Layers)
1142 rtsolution%Down_Solar_Radiance = rtv%Down_Solar_Radiance
1143 rtsolution%Surface_Planck_Radiance = rtv%Planck_Surface
1146 no = rtsolution%n_Layers
1147 na = rtv%n_Added_Layers
1151 rtsolution%Upwelling_Radiance(1:no) = rtv%e_Level_Rad_UP(na+1:nt)
1152 rtsolution%Upwelling_Overcast_Radiance(1:no) = rtv%e_Cloud_Radiance_UP(na+1:nt)
1157 rtsolution%Radiance = rtsolution%Radiance + radiance* &
1158 cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1163 IF( rtv%mth_Azi == 0 ) &
1167 rtsolution%Radiance, &
1168 rtsolution%Brightness_Temperature )
1275 SfcOptics , & ! Input
1276 RTSolution , & ! Input
1277 GeometryInfo , & ! Input
1278 Radiance_TL , & ! Input
1279 Scattering_Radiance_TL , & ! Input
1280 SensorIndex , & ! Input
1281 ChannelIndex , & ! Input
1284 result( error_status )
1289 REAL(fp) ,
INTENT(IN) :: radiance_tl
1290 REAL(fp) ,
INTENT(IN) :: scattering_radiance_tl(:)
1291 INTEGER ,
INTENT(IN) :: sensorindex
1292 INTEGER ,
INTENT(IN) :: channelindex
1297 INTEGER :: error_status
1299 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Output_TL' 1304 IF( rtv%Scattering_RT )
THEN 1306 rtsolution_tl%Radiance = rtsolution_tl%Radiance + scattering_radiance_tl( sfcoptics%Index_Sat_Ang )* &
1307 cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1311 rtsolution_tl%Radiance = rtsolution_tl%Radiance + radiance_tl* &
1312 cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1318 IF( rtv%mth_Azi == 0 ) &
1322 rtsolution%Radiance , &
1323 rtsolution_tl%Radiance , &
1324 rtsolution_tl%Brightness_Temperature )
1510 Atmosphere , & ! Input
1512 AtmOptics , & ! Input
1513 SfcOptics , & ! Input
1516 GeometryInfo , & ! Input
1517 SensorIndex , & ! Input
1518 ChannelIndex , & ! Input
1520 AtmOptics_AD , & ! Output
1521 SfcOptics_AD , & ! Output
1522 Planck_Surface_AD , & ! Output
1523 Planck_Atmosphere_AD , & ! Output
1524 User_Emissivity_AD , & ! Output
1525 Atmosphere_AD , & ! Output
1526 Surface_AD , & ! Output
1527 RTSolution_AD , & ! Output
1529 result( error_status )
1535 REAL(fp) ,
INTENT(IN OUT) :: pff_ad(:,:,:)
1536 REAL(fp) ,
INTENT(IN OUT) :: pbb_ad(:,:,:)
1538 INTEGER ,
INTENT(IN) :: sensorindex
1539 INTEGER ,
INTENT(IN) :: channelindex
1540 INTEGER ,
INTENT(IN) :: nz
1543 REAL(fp) ,
INTENT(IN OUT) :: planck_surface_ad
1544 REAL(fp) ,
INTENT(IN OUT) :: planck_atmosphere_ad(0:)
1545 REAL(fp) ,
INTENT(OUT) :: user_emissivity_ad
1552 INTEGER :: error_status
1554 CHARACTER(*),
PARAMETER :: routine_name =
'Assign_Common_Output_AD' 1556 CHARACTER(256) :: message
1557 INTEGER :: no, na, nt
1565 IF( rtv%mth_Azi == 0 )
THEN 1570 sfcoptics%Surface_Temperature , &
1571 planck_surface_ad , &
1572 sfcoptics_ad%Surface_Temperature )
1573 planck_surface_ad =
zero 1575 DO k = 1, atmosphere%n_Layers
1579 atmosphere%Temperature(k) , &
1580 planck_atmosphere_ad(k) , &
1581 atmosphere_ad%Temperature(k) )
1582 planck_atmosphere_ad(k) =
zero 1585 planck_surface_ad =
zero 1586 planck_atmosphere_ad =
zero 1596 sfcoptics_ad%n_Stokes = 1
1602 sfcoptics_ad%n_Angles = sfcoptics%n_Angles
1603 sfcoptics_ad%Angle = sfcoptics%Angle
1604 sfcoptics_ad%Weight = sfcoptics%Weight
1616 IF( rtv%Scattering_RT )
THEN 1617 user_emissivity_ad =
zero 1618 IF( rtv%Diffuse_Surface)
THEN 1620 user_emissivity_ad = user_emissivity_ad - &
1621 (sum(sfcoptics_ad%Reflectivity(1:nz,1,i,1))*sfcoptics%Weight(i))
1625 user_emissivity_ad = user_emissivity_ad - sfcoptics_ad%Reflectivity(i,1,i,1)
1631 rtsolution_ad%Surface_Emissivity = user_emissivity_ad
1633 rtsolution_ad%Surface_Emissivity = sfcoptics_ad%Emissivity(sfcoptics_ad%Index_Sat_Ang,1) - &
1634 sfcoptics_ad%Reflectivity(1,1,1,1)
1641 IF ( sfcoptics%Compute )
THEN 1651 IF ( error_status /=
success )
THEN 1652 WRITE( message,
'("Error computing SfcOptics_AD for ",a," channel ",i0)' ) &
1653 trim(
sc(sensorindex)%Sensor_Id), &
1654 sc(sensorindex)%Sensor_Channel(channelindex)
1663 IF ( rtv%Scattering_RT )
THEN 1677 no = rtsolution_ad%n_Layers
1678 na = rtv%n_Added_Layers
1682 rtsolution_ad%Layer_Optical_Depth(1:no) = atmoptics_ad%Optical_Depth(na+1:nt)
1740 AtmOptics, & ! Input
1743 TYPE(CRTM_AtmOptics_type),
INTENT(IN) :: AtmOptics
1744 TYPE(RTV_type),
INTENT(IN OUT) :: RTV
1746 INTEGER :: i, j, k, l, ifac, jn
1754 DO i = 1, rtv%n_Angles
1756 atmoptics%n_Legendre_Terms, &
1757 rtv%COS_Angle(i) , &
1761 IF(rtv%Solar_Flag_true) &
1763 atmoptics%n_Legendre_Terms , &
1765 rtv%Pleg(0:,rtv%n_Angles+1) )
1767 IF( rtv%Solar_Flag_true )
THEN 1768 jn = rtv%n_Angles + 1
1777 layer_loop:
DO k = 1, rtv%n_Layers
1789 DO i = 1, rtv%n_Angles
1794 DO l = rtv%mth_Azi, atmoptics%n_Legendre_Terms - 1
1795 ifac = (-1) ** (l - rtv%mth_Azi)
1796 rtv%Off(i,j,k) = rtv%Off(i,j,k) + ( atmoptics%Phase_Coefficient(l,1,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j) )
1797 rtv%Obb(i,j,k) = rtv%Obb(i,j,k) + ( atmoptics%Phase_Coefficient(l,1,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j)*ifac )
1800 rtv%Pff(i,j,k) = rtv%Off(i,j,k)
1801 rtv%Pbb(i,j,k) = rtv%Obb(i,j,k)
1804 IF ( rtv%mth_Azi == 0 )
THEN 1815 END IF significant_scattering
1879 AtmOptics, & ! FWD Input
1880 AtmOptics_TL, & ! TL Input
1881 Pff_TL, & ! TL Output
1882 Pbb_TL, & ! TL Output
1885 TYPE(CRTM_AtmOptics_type) ,
INTENT(IN) :: AtmOptics
1886 TYPE(CRTM_AtmOptics_type) ,
INTENT(IN) :: AtmOptics_TL
1887 REAL(fp) ,
INTENT(OUT) :: Pff_TL(:,:,:)
1888 REAL(fp) ,
INTENT(OUT) :: Pbb_TL(:,:,:)
1889 TYPE(RTV_type) ,
INTENT(IN) :: RTV
1891 INTEGER :: i, j, k, l, nZ, ifac, jn
1892 REAL(fp),
DIMENSION(RTV%n_Angles,RTV%n_Angles+1) :: Lff, Lbb
1901 IF( rtv%Solar_Flag_true )
THEN 1902 jn = rtv%n_Angles + 1
1907 layer_loop:
DO k = 1, rtv%n_Layers
1916 lff(1:nz,1:jn) = rtv%Off(1:nz,1:jn,k)
1917 lbb(1:nz,1:jn) = rtv%Obb(1:nz,1:jn,k)
1921 DO i = 1, rtv%n_Angles
1923 pff_tl(i,j,k) =
zero 1924 pbb_tl(i,j,k) =
zero 1926 DO l = rtv%mth_Azi, atmoptics%n_Legendre_Terms - 1
1927 ifac = (-1) ** (l - rtv%mth_Azi)
1928 pff_tl(i,j,k) = pff_tl(i,j,k) + ( atmoptics_tl%Phase_Coefficient(l,1,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j) )
1929 pbb_tl(i,j,k) = pbb_tl(i,j,k) + ( atmoptics_tl%Phase_Coefficient(l,1,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j)*ifac )
1934 IF ( rtv%mth_Azi == 0 )
THEN 1935 IF ( rtv%Off(i,j,k) <
zero )
THEN 1936 pff_tl(i,j,k) =
zero 1940 IF ( rtv%Obb(i,j,k) <
zero )
THEN 1941 pbb_tl(i,j,k) =
zero 1949 IF ( rtv%mth_Azi == 0 )
THEN 1958 END IF significant_scattering
2033 AtmOptics, & ! FWD Input
2034 Pff_AD, & ! AD Input
2035 Pbb_AD, & ! AD Input
2036 AtmOptics_AD, & ! AD Output
2039 TYPE(CRTM_AtmOptics_type),
INTENT(IN) :: AtmOptics
2040 REAL(fp) ,
INTENT(IN OUT) :: Pff_AD(:,:,:)
2041 REAL(fp) ,
INTENT(IN OUT) :: Pbb_AD(:,:,:)
2042 TYPE(CRTM_AtmOptics_type),
INTENT(IN OUT) :: AtmOptics_AD
2043 TYPE(RTV_type) ,
INTENT(IN) :: RTV
2045 INTEGER :: i, j, k, l, nZ, ifac, jn
2046 REAL(fp),
DIMENSION(RTV%n_Angles,RTV%n_Angles+1) :: Lff, Lbb
2055 IF( rtv%Solar_Flag_true )
THEN 2056 jn = rtv%n_Angles + 1
2061 layer_loop:
DO k = 1, rtv%n_Layers
2075 lff(1:nz,1:jn) = rtv%Off(1:nz,1:jn,k)
2076 lbb(1:nz,1:jn) = rtv%Obb(1:nz,1:jn,k)
2078 IF ( rtv%mth_Azi == 0 )
THEN 2080 DO i = 1, rtv%n_Angles
2099 DO i = 1, rtv%n_Angles
2103 IF ( rtv%mth_Azi == 0 )
THEN 2104 IF ( rtv%Off(i,j,k) <
zero) pff_ad(i,j,k) =
zero 2105 IF ( rtv%Obb(i,j,k) <
zero) pbb_ad(i,j,k) =
zero 2108 DO l = rtv%mth_Azi, atmoptics%n_Legendre_Terms - 1
2109 ifac = (-1) ** (l - rtv%mth_Azi)
2110 atmoptics_ad%Phase_Coefficient(l,1,k) = atmoptics_ad%Phase_Coefficient(l,1,k) + &
2111 ( pff_ad(i,j,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j) )
2112 atmoptics_ad%Phase_Coefficient(l,1,k) = atmoptics_ad%Phase_Coefficient(l,1,k) + &
2113 ( pbb_ad(i,j,k)*rtv%Pleg(l,i)*rtv%Pleg(l,j)*ifac )
2116 pff_ad(i,j,k) =
zero 2117 pbb_ad(i,j,k) =
zero 2122 END IF significant_scattering
2134 INTEGER,
INTENT(IN) :: k
2135 TYPE(RTV_type),
INTENT(IN OUT) :: RTV
2142 rtv%Sum_Fac(0,k)=
zero 2143 DO i = 1, rtv%n_Streams
2144 rtv%n_Factor(i,k)=
zero 2146 rtv%n_Factor(i,k)=rtv%n_Factor(i,k)+(rtv%Pff(i,j,k)+rtv%Pbb(i,j,k))*rtv%COS_Weight(j)
2149 rtv%Pff(i,j,k)=rtv%Pff(i,j,k)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k))
2150 rtv%Pbb(i,j,k)=rtv%Pbb(i,j,k)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k))
2152 rtv%Sum_Fac(i,k)=
zero 2155 rtv%Pff(j,i,k)=rtv%Pff(i,j,k)
2156 rtv%Pbb(j,i,k)=rtv%Pbb(i,j,k)
2159 rtv%Sum_Fac(i,k)=rtv%Sum_Fac(i,k) + (rtv%Pff(i+1,j,k)+rtv%Pbb(i+1,j,k))*rtv%COS_Weight(j)
2164 IF( rtv%n_Streams < nz )
THEN 2166 rtv%n_Factor(nz,k) = rtv%Sum_Fac(nz-1,k)
2168 rtv%Pff(j,nz,k) = rtv%Pff(j,nz,k)/rtv%n_Factor(nz,k)
2169 rtv%Pbb(j,nz,k) = rtv%Pbb(j,nz,k)/rtv%n_Factor(nz,k)
2172 rtv%Pff(nz,j,k) = rtv%Pff(j,nz,k)
2173 rtv%Pbb(nz,j,k) = rtv%Pbb(j,nz,k)
2182 INTEGER ,
INTENT(IN) :: k
2183 TYPE(RTV_type),
INTENT(IN) :: RTV
2184 REAL(fp) ,
INTENT(IN) :: Pff(:,:)
2185 REAL(fp) ,
INTENT(IN) :: Pbb(:,:)
2186 REAL(fp) ,
INTENT(IN OUT) :: Pff_TL(:,:)
2187 REAL(fp) ,
INTENT(IN OUT) :: Pbb_TL(:,:)
2189 REAL(fp) :: n_Factor_TL
2190 REAL(fp) :: Sum_Fac_TL(0:RTV%n_Angles)
2196 sum_fac_tl(0) =
zero 2197 DO i = 1, rtv%n_Streams
2200 n_factor_tl=n_factor_tl+(pff_tl(i,j)+pbb_tl(i,j))*rtv%COS_Weight(j)
2203 pff_tl(i,j) = pff_tl(i,j)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k)) - &
2204 pff(i,j)/rtv%n_Factor(i,k)/rtv%n_Factor(i,k)*n_factor_tl*(
one-rtv%Sum_Fac(i-1,k)) - &
2205 pff(i,j)/rtv%n_Factor(i,k)*sum_fac_tl(i-1)
2207 pbb_tl(i,j) = pbb_tl(i,j)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k)) - &
2208 pbb(i,j)/rtv%n_Factor(i,k)/rtv%n_Factor(i,k)*n_factor_tl*(
one-rtv%Sum_Fac(i-1,k)) - &
2209 pbb(i,j)/rtv%n_Factor(i,k)*sum_fac_tl(i-1)
2215 pff_tl(j,i) = pff_tl(i,j)
2216 pbb_tl(j,i) = pbb_tl(i,j)
2219 sum_fac_tl(i)=sum_fac_tl(i) + (pff_tl(j,i+1)+pbb_tl(j,i+1))*rtv%COS_Weight(j)
2224 IF( rtv%n_Streams < nz )
THEN 2226 n_factor_tl = sum_fac_tl(nz-1)
2228 pff_tl(j,nz) = pff_tl(j,nz)/rtv%n_Factor(nz,k) - &
2229 pff(j,nz)/rtv%n_Factor(nz,k)/rtv%n_Factor(nz,k)*n_factor_tl
2231 pbb_tl(j,nz) = pbb_tl(j,nz)/rtv%n_Factor(nz,k) - &
2232 pbb(j,nz)/rtv%n_Factor(nz,k)/rtv%n_Factor(nz,k)*n_factor_tl
2235 pff_tl(nz,j) = pff_tl(j,nz)
2236 pbb_tl(nz,j) = pbb_tl(j,nz)
2245 INTEGER ,
INTENT(IN) :: k
2246 TYPE(RTV_type),
INTENT(IN) :: RTV
2247 REAL(fp) ,
INTENT(IN) :: Pff(:,:)
2248 REAL(fp) ,
INTENT(IN) :: Pbb(:,:)
2249 REAL(fp) ,
INTENT(IN OUT) :: Pff_AD(:,:)
2250 REAL(fp) ,
INTENT(IN OUT) :: Pbb_AD(:,:)
2253 REAL(fp) :: n_Factor_AD
2254 REAL(fp) :: Sum_Fac_AD(0:RTV%n_Angles)
2261 IF( rtv%n_Streams < nz )
THEN 2266 pff_ad(j,nz) = pff_ad(j,nz) + pff_ad(nz,j)
2268 pbb_ad(j,nz) = pbb_ad(j,nz) + pbb_ad(nz,j)
2272 n_factor_ad = n_factor_ad - pff(j,nz)/rtv%n_Factor(nz,k)/rtv%n_Factor(nz,k)*pff_ad(j,nz)
2273 pff_ad(j,nz) = pff_ad(j,nz)/rtv%n_Factor(nz,k)
2275 n_factor_ad = n_factor_ad - pbb(j,nz)/rtv%n_Factor(nz,k)/rtv%n_Factor(nz,k)*pbb_ad(j,nz)
2276 pbb_ad(j,nz) = pbb_ad(j,nz)/rtv%n_Factor(nz,k)
2278 sum_fac_ad(nz-1) = n_factor_ad
2282 DO i = rtv%n_Streams, 1, -1
2286 pbb_ad(j,i+1) = pbb_ad(j,i+1) + sum_fac_ad(i)*rtv%COS_Weight(j)
2287 pff_ad(j,i+1) = pff_ad(j,i+1) + sum_fac_ad(i)*rtv%COS_Weight(j)
2290 pff_ad(i,j) = pff_ad(i,j) + pff_ad(j,i)
2292 pbb_ad(i,j) = pbb_ad(i,j) + pbb_ad(j,i)
2296 sum_fac_ad(i) =
zero 2298 sum_fac_ad(i-1) = sum_fac_ad(i-1) - pbb(i,j)/rtv%n_Factor(i,k)*pbb_ad(i,j)
2299 n_factor_ad = n_factor_ad -pbb(i,j)/rtv%n_Factor(i,k)/rtv%n_Factor(i,k) * &
2300 pbb_ad(i,j)*(
one-rtv%Sum_Fac(i-1,k))
2301 pbb_ad(i,j) = pbb_ad(i,j)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k))
2303 sum_fac_ad(i-1) = sum_fac_ad(i-1) - pff(i,j)/rtv%n_Factor(i,k)*pff_ad(i,j)
2304 n_factor_ad = n_factor_ad -pff(i,j)/rtv%n_Factor(i,k)/rtv%n_Factor(i,k) * &
2305 pff_ad(i,j)*(
one-rtv%Sum_Fac(i-1,k))
2306 pff_ad(i,j) = pff_ad(i,j)/rtv%n_Factor(i,k)*(
one-rtv%Sum_Fac(i-1,k))
2309 pbb_ad(i,j) = pbb_ad(i,j) + n_factor_ad*rtv%COS_Weight(j)
2310 pff_ad(i,j) = pff_ad(i,j) + n_factor_ad*rtv%COS_Weight(j)
2314 sum_fac_ad(0) =
zero integer function, public crtm_compute_sfcoptics_tl(Surface, SfcOptics, Surface_TL, GeometryInfo, SensorIndex, ChannelIndex, SfcOptics_TL, iVar)
integer function, public assign_common_output(Atmosphere, SfcOptics, GeometryInfo, SensorIndex, ChannelIndex, RTV, RTSolution)
integer function, public assign_common_output_ad(Atmosphere, Surface, AtmOptics, SfcOptics, Pff_AD, Pbb_AD, GeometryInfo, SensorIndex, ChannelIndex, nZ, AtmOptics_AD, SfcOptics_AD, Planck_Surface_AD, Planck_Atmosphere_AD, User_Emissivity_AD, Atmosphere_AD, Surface_AD, RTSolution_AD, RTV)
real(fp), parameter, public zero
subroutine normalize_phase(k, RTV)
subroutine, public legendre_m(MF, N, U, PL)
subroutine normalize_phase_ad(k, RTV, Pff, Pbb, Pff_AD, Pbb_AD)
integer, parameter, public fp
integer function, public assign_common_output_tl(SfcOptics, RTSolution, GeometryInfo, Radiance_TL, Scattering_Radiance_TL, SensorIndex, ChannelIndex, RTV, RTSolution_TL)
character(*), parameter module_version_id
real(fp), parameter, public angle_threshold
real(fp), parameter, public scattering_albedo_threshold
subroutine, public crtm_planck_radiance_ad(n, l, Temperature, Radiance_AD, Temperature_AD)
pure logical function, public crtm_include_scattering(atmoptics)
subroutine crtm_phase_matrix_ad(AtmOptics, Pff_AD, Pbb_AD, AtmOptics_AD, RTV)
integer function, public crtm_compute_sfcoptics(Surface, GeometryInfo, SensorIndex, ChannelIndex, SfcOptics, iVar)
real(fp), parameter, public secant_diffusivity
subroutine normalize_phase_tl(k, RTV, Pff, Pbb, Pff_TL, Pbb_TL)
subroutine crtm_phase_matrix_tl(AtmOptics, AtmOptics_TL, Pff_TL, Pbb_TL, RTV)
subroutine, public crtm_planck_temperature_ad(n, l, Radiance, Temperature_AD, Radiance_AD)
elemental logical function, public crtm_rtsolution_associated(RTSolution)
real(fp), parameter, public one
recursive subroutine, public display_message(Routine_Name, Message, Error_State, Message_Log)
subroutine, public double_gauss_quadrature(NUM, ABSCISSAS, WEIGHTS)
subroutine, public crtm_planck_temperature_tl(n, l, Radiance, Radiance_TL, Temperature_TL)
integer function, public assign_common_input_ad(SfcOptics, RTSolution, GeometryInfo, SensorIndex, ChannelIndex, RTSolution_AD, SfcOptics_AD, Planck_Surface_AD, Planck_Atmosphere_AD, Radiance_AD, nz, RTV)
real(fp), parameter, public degrees_to_radians
type(spccoeff_type), dimension(:), allocatable, save, public sc
integer function, public assign_common_input(Atmosphere, Surface, AtmOptics, SfcOptics, GeometryInfo, SensorIndex, ChannelIndex, RTSolution, nz, RTV)
integer function, public crtm_compute_sfcoptics_ad(Surface, SfcOptics, SfcOptics_AD, GeometryInfo, SensorIndex, ChannelIndex, Surface_AD, iVar)
subroutine crtm_phase_matrix(AtmOptics, RTV)
integer function, public assign_common_input_tl(Atmosphere, Surface, AtmOptics, SfcOptics, Atmosphere_TL, Surface_TL, AtmOptics_TL, SfcOptics_TL, GeometryInfo, SensorIndex, ChannelIndex, RTSolution_TL, nz, User_Emissivity_TL, Direct_Reflectivity_TL, Planck_Surface_TL, Planck_Atmosphere_TL, Pff_TL, Pbb_TL, RTV)
real(fp), parameter, public phase_threshold
subroutine, public crtm_planck_temperature(n, l, Radiance, Temperature)
subroutine, public crtm_planck_radiance_tl(n, l, Temperature, Temperature_TL, Radiance_TL)
integer, parameter, public success
subroutine, public crtm_planck_radiance(n, l, Temperature, Radiance)
real(fp), parameter, public pi