FV3 Bundle
Common_RTSolution.f90
Go to the documentation of this file.
1 !
2 ! Common_RTSolution
3 !
4 ! Module containing common ancillary functions for
5 ! available radiative transfer algorithms.
6 !
7 ! CREATION HISTORY:
8 ! Written by: David Neil Groff, IMSG at EMC; david.groff@noaa.gov
9 ! 09-Sep-2010
10 
11 
13 
14  ! ------------------
15  ! Environment set up
16  ! ------------------
17  ! Module use statements
18  USE type_kinds, ONLY: fp
19  USE crtm_parameters, ONLY: one, zero, pi, &
33  USE crtm_spccoeff
36  USE rtv_define
37  USE crtm_sfcoptics
39  USE crtm_utility
41 
42  ! Disable all implicit typing
43  IMPLICIT NONE
44 
45  ! ------------
46  ! Visibilities
47  ! ------------
48  ! Everything private by default
49  PRIVATE
50 
51  PUBLIC :: assign_common_input
52  PUBLIC :: assign_common_output
53  PUBLIC :: assign_common_input_tl
54  PUBLIC :: assign_common_output_tl
55  PUBLIC :: assign_common_input_ad
56  PUBLIC :: assign_common_output_ad
57 
58  ! -----------------
59  ! Module parameters
60  ! -----------------
61  ! Version Id for the module
62  CHARACTER(*), PARAMETER :: module_version_id = &
63  '$Id: $'
64 
65 CONTAINS
66 
67 !################################################################################
68 !################################################################################
69 !## ##
70 !## ## PUBLIC MODULE ROUTINES ## ##
71 !## ##
72 !################################################################################
73 !################################################################################
74 
75 !-------------------------------------------------------------------------
76 !
77 ! NAME:
78 ! Assign_Common_Input
79 !
80 ! PURPOSE:
81 ! Function to assign input that is used when calling both the ADA and
82 ! SOI RT algorithms.
83 !
84 ! CALLING SEQUENCE:
85 ! Error_Status = Assign_Common_Input( Atmosphere , & ! Input
86 ! Surface , & ! Input
87 ! AtmOptics , & ! Input
88 ! SfcOptics , & ! Input
89 ! GeometryInfo, & ! Input
90 ! SensorIndex , & ! Input
91 ! ChannelIndex, & ! Input
92 ! RTSolution , & ! Output
93 ! RTV ) ! Internal variable output
94 !
95 ! INPUT ARGUMENTS:
96 ! Atmosphere: Structure containing the atmospheric state data.
97 ! UNITS: N/A
98 ! TYPE: CRTM_Atmosphere_type
99 ! DIMENSION: Scalar
100 ! ATTRIBUTES: INTENT(IN)
101 !
102 ! Surface: Structure containing the surface state data.
103 ! UNITS: N/A
104 ! TYPE: CRTM_Surface_type
105 ! DIMENSION: Scalar
106 ! ATTRIBUTES: INTENT(IN)
107 !
108 ! AtmOptics: Structure containing the combined atmospheric
109 ! optical properties for gaseous absorption, clouds,
110 ! and aerosols.
111 ! UNITS: N/A
112 ! TYPE: CRTM_AtmOptics_type
113 ! DIMENSION: Scalar
114 ! ATTRIBUTES: INTENT(IN)
115 !
116 ! SfcOptics: Structure containing the surface optical properties
117 ! data. Argument is defined as INTENT (IN OUT ) as
118 ! different RT algorithms may compute the surface
119 ! optics properties before this routine is called.
120 ! UNITS: N/A
121 ! TYPE: CRTM_SfcOptics_type
122 ! DIMENSION: Scalar
123 ! ATTRIBUTES: INTENT(IN OUT)
124 !
125 ! GeometryInfo: Structure containing the view geometry data.
126 ! UNITS: N/A
127 ! TYPE: CRTM_GeometryInfo_type
128 ! DIMENSION: Scalar
129 ! ATTRIBUTES: INTENT(IN OUT)
130 !
131 ! SensorIndex: Sensor index id. This is a unique index associated
132 ! with a (supported) sensor used to access the
133 ! shared coefficient data for a particular sensor.
134 ! See the ChannelIndex argument.
135 ! UNITS: N/A
136 ! TYPE: INTEGER
137 ! DIMENSION: Scalar
138 ! ATTRIBUTES: INTENT(IN)
139 !
140 ! ChannelIndex: Channel index id. This is a unique index associated
141 ! with a (supported) sensor channel used to access the
142 ! shared coefficient data for a particular sensor's
143 ! channel.
144 ! See the SensorIndex argument.
145 ! UNITS: N/A
146 ! TYPE: INTEGER
147 ! DIMENSION: Scalar
148 ! ATTRIBUTES: INTENT(IN)
149 !
150 ! OUTPUT ARGUMENTS:
151 ! RTSolution: Structure containing the soluition to the RT equation
152 ! for the given inputs.
153 ! UNITS: N/A
154 ! TYPE: CRTM_RTSolution_type
155 ! DIMENSION: Scalar
156 ! ATTRIBUTES: INTENT(IN OUT)
157 !
158 ! nz Assignment of RTV%n_Angles for
159 ! existing RTSolution interfaces
160 ! UNITS: N/A
161 ! TYPE: INTEGER
162 ! DIMENSION: Scalar
163 ! ATTRIBUTES: INTENT(OUT)
164 !
165 !
166 ! RTV: Structure containing internal variables required for
167 ! tangent-linear or adjoint model calls in CRTM_RTSolution.
168 ! UNITS: N/A
169 ! TYPE: RTV_type
170 ! DIMENSION: Scalar
171 ! ATTRIBUTES: INTENT(IN OUT)
172 !
173 ! FUNCTION RESULT:
174 ! Error_Status: The return value is an integer defining the error status.
175 ! The error codes are defined in the Message_Handler module.
176 ! If == SUCCESS the computation was sucessful
177 ! == FAILURE an unrecoverable error occurred
178 ! UNITS: N/A
179 ! TYPE: INTEGER
180 ! DIMENSION: Scalar
181 !
182 ! COMMENTS:
183 ! Note the INTENT on the output RTSolution argument is IN OUT rather than
184 ! just OUT. This is necessary because the argument is defined upon
185 ! input. To prevent memory leaks, the IN OUT INTENT is a must.
186 !
187 !--------------------------------------------------------------------------------
188 
189  FUNCTION assign_common_input( &
190  Atmosphere , & ! Input
191  Surface , & ! Input
192  AtmOptics , & ! Input
193  SfcOptics , & ! Input
194  GeometryInfo, & ! Input
195  SensorIndex , & ! Input
196  ChannelIndex, & ! Input
197  RTSolution , & ! Output
198  nz , & ! Input
199  RTV ) & ! Internal variable output
200  result( error_status )
201  ! Arguments
202  TYPE(crtm_atmosphere_type), INTENT(IN) :: atmosphere
203  TYPE(crtm_surface_type), INTENT(IN) :: surface
204  TYPE(crtm_atmoptics_type), INTENT(IN) :: atmoptics
205  TYPE(crtm_sfcoptics_type), INTENT(IN OUT) :: sfcoptics
206  TYPE(crtm_geometryinfo_type), INTENT(IN OUT) :: geometryinfo
207  INTEGER, INTENT(IN) :: sensorindex
208  INTEGER, INTENT(IN) :: channelindex
209  TYPE(crtm_rtsolution_type), INTENT(IN OUT) :: rtsolution
210  INTEGER, INTENT(OUT) :: nz
211  TYPE(rtv_type), INTENT(IN OUT) :: rtv
212  ! Function result
213  INTEGER :: error_status
214  ! Local parameters
215  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Input'
216  ! Local variables
217  CHARACTER(256) :: message
218  INTEGER :: no, na, nt
219  INTEGER :: i, k
220  REAL(fp) :: factor ! SfcOptics quadrature weights normalisation factor
221  REAL(fp) :: user_emissivity, direct_reflectivity
222 
223  ! Set Up
224  error_status = success
225 
226  geometryinfo%Cosine_Sensor_Zenith = one / geometryinfo%Secant_Sensor_Zenith
227 
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
231 
232  ! Save the optical depth if required
233  IF ( crtm_rtsolution_associated( rtsolution ) ) THEN
234  ! Shorter names for indexing
235  no = rtsolution%n_Layers ! Original no. of layers
236  na = rtv%n_Added_Layers ! No. of added layers
237  nt = rtv%n_Layers ! Current total no. of layers
238  ! Assign only the optical depth profile
239  ! defined by the user input layering
240  rtsolution%Layer_Optical_Depth(1:no) = atmoptics%Optical_Depth(na+1:nt)
241  END IF
242  ! Required SpcCoeff components
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
248 
249  ! Determine the surface emission behavior
250  ! By default, surface is SPECULAR.
251  rtv%Diffuse_Surface = .false.
252 
253  ! -------------------------------------------
254  ! Determine the quadrature angles and weights
255  ! -------------------------------------------
256  ! Default is to assume scattering RT
257  rtv%Scattering_RT = .false.
258 
259  ! Check for scattering conditions
260  determine_stream_angles: IF( rtsolution%Scattering_FLAG .AND. crtm_include_scattering(atmoptics) ) THEN
261 
262  ! --------------------------
263  ! Set the scattering RT flag
264  ! --------------------------
265  rtv%Scattering_RT = .true.
266 
267  ! ---------------------------------------
268  ! Compute the Gaussian angles and weights
269  ! ---------------------------------------
270  CALL double_gauss_quadrature( rtv%n_Streams, rtv%COS_Angle, rtv%COS_Weight )
271 
272 
273  ! ----------------------------------------------------------------
274  ! Is an additional stream required for the satellite zenith angle?
275  ! That is, is the satellite zenith angle sufficiently different
276  ! from the stream angles?
277  ! ----------------------------------------------------------------
278  ! Search for a matching angle
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 ! Identify the sensor zenith angle
282  rtv%n_Angles = rtv%n_Streams
283  ! *****FLAW*****
284  ! CAN WE REPLACE THIS CONSTRUCT WITH ANOTHER
285  GO TO 100
286  ! *****FLAW*****
287  END IF
288  END DO
289 
290  ! No match found, so flag an additional
291  ! spot for the satellite zenith angle
292  rtv%n_Angles = rtv%n_Streams + 1
293  sfcoptics%Index_Sat_Ang = rtv%n_Angles ! Identify the sensor zenith angle
294  rtv%COS_Angle( rtv%n_Angles ) = geometryinfo%Cosine_Sensor_Zenith
295  rtv%COS_Weight( rtv%n_Angles ) = zero
296 
297  100 CONTINUE
298 
299  ! ------------------------
300  ! Compute the phase matrix
301  ! ------------------------
302  CALL crtm_phase_matrix( atmoptics, & ! Input
303  rtv ) ! Internal variable
304 
305  ELSE IF( rtv%Diffuse_Surface) THEN
306 
307  ! ------------------
308  ! Lambertian surface
309  ! ------------------
310  ! The default diffusivity angle
311  rtv%n_Streams = 1
312  rtv%COS_Angle( 1 ) = one/secant_diffusivity
313  rtv%COS_Weight( 1 ) = one
314  ! The satellite zenith angle
315  rtv%n_Angles = rtv%n_Streams + 1
316  sfcoptics%Index_Sat_Ang = rtv%n_Angles ! Identify the sensor zenith angle
317  rtv%COS_Angle( rtv%n_Angles ) = geometryinfo%Cosine_Sensor_Zenith
318  rtv%COS_Weight( rtv%n_Angles ) = zero
319 
320  ELSE
321 
322  ! ----------------
323  ! Specular surface
324  ! ----------------
325  rtv%n_Streams = 0
326  rtv%n_Angles = 1
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
330 
331  END IF determine_stream_angles
332  ! -------------------------------------------------
333  ! Copy the RTV angles to a short name for use below
334  ! -------------------------------------------------
335  nz = rtv%n_Angles
336 
337  ! ---------------------------------------------------------------------
338  ! Assign the number of Stokes parameters. Currently, this is set to 1
339  ! for decoupled polarization between surface and atmosphere.
340  ! Remember for polarised microwave instruments, each frequency's Stokes
341  ! parameter is treated as a separate channel.
342  ! ---------------------------------------------------------------------
343  sfcoptics%n_Stokes = 1
344 
345  ! ----------------------------------------------------
346  ! Assign the angles and weights to SfcOptics structure
347  ! ----------------------------------------------------
348  sfcoptics%n_Angles = rtv%n_Angles
349  factor = zero
350  DO i = 1, nz
351  sfcoptics%Angle(i) = acos( rtv%COS_Angle(i) ) / degrees_to_radians
352  sfcoptics%Weight(i) = rtv%COS_Angle(i) * rtv%COS_Weight(i)
353  ! Compute the weight normalisation factor
354  factor = factor + sfcoptics%Weight(i)
355  END DO
356  ! Normalise the quadrature weights
357  IF( factor > zero) sfcoptics%Weight(1:nz) = sfcoptics%Weight(1:nz) / factor
358 
359  ! --------------------------------
360  ! Populate the SfcOptics structure
361  ! --------------------------------
362  IF ( sfcoptics%Compute ) THEN
363  error_status = crtm_compute_sfcoptics( &
364  surface , & ! Input
365  geometryinfo, & ! Input
366  sensorindex , & ! Input
367  channelindex, & ! Input
368  sfcoptics , & ! In/Output
369  rtv%SOV ) ! Internal variable output
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)
374  CALL display_message( routine_name, trim(message), error_status )
375  RETURN
376  END IF
377 
378  ELSE
379 
380  IF( rtv%Scattering_RT ) THEN
381  ! Replicate the user emissivity for all angles
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
387  ! Replicate the user reflectivities for all angles
388  sfcoptics%Direct_Reflectivity(1:nz,1) = direct_reflectivity
389  IF( rtv%Diffuse_Surface) THEN
390  DO i = 1, nz
391  sfcoptics%Reflectivity(1:nz, 1, i, 1) = (one-sfcoptics%Emissivity(i,1))*sfcoptics%Weight(i)
392  END DO
393  ELSE ! Specular surface
394  DO i = 1, nz
395  sfcoptics%Reflectivity(i, 1, i, 1) = (one-sfcoptics%Emissivity(i,1))
396  END DO
397  END IF
398  ELSE
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
402  END IF
403 
404  END IF
405 
406  ! ---------------------------
407  ! Save the surface properties
408  ! ---------------------------
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 )
411 
412  ! ------------------------
413  ! Compute Planck radiances
414  ! ------------------------
415  IF( rtv%mth_Azi == 0 ) THEN
416  DO k = 1, atmosphere%n_Layers
417  CALL crtm_planck_radiance( &
418  sensorindex , & ! Input
419  channelindex , & ! Input
420  atmosphere%Temperature(k), & ! Input
421  rtv%Planck_Atmosphere(k) ) ! Output
422  END DO
423  ! Surface radiance
424  CALL crtm_planck_radiance( &
425  sensorindex , & ! Input
426  channelindex , & ! Input
427  sfcoptics%Surface_Temperature, & ! Input
428  rtv%Planck_Surface ) ! Output
429  ELSE
430  rtv%Planck_Atmosphere = zero
431  rtv%Planck_Surface = zero
432  END IF
433 
434  END FUNCTION assign_common_input
435 
436 !--------------------------------------------------------------------------------
437 !
438 ! NAME:
439 ! Assign_Common_Input_TL
440 !
441 ! PURPOSE:
442 ! Function to assign tangent-linear input common to both the ADA and
443 ! SOI RT models.
444 !
445 ! CALLING SEQUENCE:
446 ! Error_Status = Assign_Common_Input_TL( Atmosphere , & ! FWD Input
447 ! Surface , & ! FWD Input
448 ! AtmOptics , & ! FWD Input
449 ! SfcOptics , & ! FWD Input
450 ! RTSolution , & ! FWD Input
451 ! Atmosphere_TL, & ! TL Input
452 ! Surface_TL , & ! TL Input
453 ! AtmOptics_TL , & ! TL Input
454 ! SfcOptics_TL , & ! TL Input
455 ! GeometryInfo , & ! Input
456 ! SensorIndex , & ! Input
457 ! ChannelIndex , & ! Input
458 ! RTSolution_TL, & ! TL Output
459 ! RTV ) ! Internal variable input
460 !
461 ! INPUT ARGUMENTS:
462 ! Atmosphere: Structure containing the atmospheric state data.
463 ! UNITS: N/A
464 ! TYPE: CRTM_Atmosphere_type
465 ! DIMENSION: Scalar
466 ! ATTRIBUTES: INTENT(IN)
467 !
468 ! Surface: Structure containing the surface state data.
469 ! UNITS: N/A
470 ! TYPE: CRTM_Surface_type
471 ! DIMENSION: Scalar
472 ! ATTRIBUTES: INTENT(IN)
473 !
474 ! AtmOptics: Structure containing the combined atmospheric
475 ! optical properties for gaseous absorption, clouds,
476 ! and aerosols.
477 ! UNITS: N/A
478 ! TYPE: CRTM_AtmOptics_type
479 ! DIMENSION: Scalar
480 ! ATTRIBUTES: INTENT(IN)
481 !
482 ! SfcOptics: Structure containing the surface optical properties
483 ! data.
484 ! UNITS: N/A
485 ! TYPE: CRTM_SfcOptics_type
486 ! DIMENSION: Scalar
487 ! ATTRIBUTES: INTENT(IN)
488 !
489 ! Atmosphere_TL: Structure containing the tangent-linear atmospheric
490 ! state data.
491 ! UNITS: N/A
492 ! TYPE: CRTM_Atmosphere_type
493 ! DIMENSION: Scalar
494 ! ATTRIBUTES: INTENT(IN)
495 !
496 ! Surface_TL: Structure containing the tangent-linear surface state data.
497 ! UNITS: N/A
498 ! TYPE: CRTM_Surface_type
499 ! DIMENSION: Scalar
500 ! ATTRIBUTES: INTENT(IN)
501 !
502 ! AtmOptics_TL: Structure containing the tangent-linear atmospheric
503 ! optical properties.
504 ! UNITS: N/A
505 ! TYPE: CRTM_AtmOptics_type
506 ! DIMENSION: Scalar
507 ! ATTRIBUTES: INTENT(IN)
508 !
509 ! SfcOptics_TL: Structure containing the tangent-linear surface optical
510 ! properties. Argument is defined as INTENT (IN OUT ) as
511 ! different RT algorithms may compute the surface optics
512 ! properties before this routine is called.
513 ! UNITS: N/A
514 ! TYPE: CRTM_SfcOptics_type
515 ! DIMENSION: Scalar
516 ! ATTRIBUTES: INTENT(IN OUT)
517 !
518 ! GeometryInfo: Structure containing the view geometry data.
519 ! UNITS: N/A
520 ! TYPE: CRTM_GeometryInfo_type
521 ! DIMENSION: Scalar
522 ! ATTRIBUTES: INTENT(IN)
523 !
524 ! SensorIndex: Sensor index id. This is a unique index associated
525 ! with a (supported) sensor used to access the
526 ! shared coefficient data for a particular sensor.
527 ! See the ChannelIndex argument.
528 ! UNITS: N/A
529 ! TYPE: INTEGER
530 ! DIMENSION: Scalar
531 ! ATTRIBUTES: INTENT(IN)
532 !
533 ! ChannelIndex: Channel index id. This is a unique index associated
534 ! with a (supported) sensor channel used to access the
535 ! shared coefficient data for a particular sensor's
536 ! channel.
537 ! See the SensorIndex argument.
538 ! UNITS: N/A
539 ! TYPE: INTEGER
540 ! DIMENSION: Scalar
541 ! ATTRIBUTES: INTENT(IN)
542 !
543 ! RTV: Structure containing internal forward model variables
544 ! required for subsequent tangent-linear or adjoint model
545 ! calls. The contents of this structure are NOT accessible
546 ! outside of the CRTM_RTSolution module.
547 ! UNITS: N/A
548 ! TYPE: RTV_type
549 ! DIMENSION: Scalar
550 ! ATTRIBUTES: INTENT(IN)
551 !
552 ! OUTPUT ARGUMENTS:
553 !
554 !
555 ! RTSolution_TL: Structure containing the solution to the tangent-linear
556 ! RT equation for the given inputs.
557 ! UNITS: N/A
558 ! TYPE: CRTM_RTSolution_type
559 ! DIMENSION: Scalar
560 ! ATTRIBUTES: INTENT(IN OUT)
561 !
562 ! nz: Integer assigned from RTV%n_Angles
563 ! UNITS: N/A
564 ! TYPE: INTEGER
565 ! DIMENSION: Scalar
566 ! ATTRIBUTES: INTENT(OUT)
567 !
568 ! User_Emissivity_TL: Tangent linear of emissivity
569 ! UNITS: N/A
570 ! TYPE: REAL(fp)
571 ! DIMENSION: Scalar
572 ! ATTRIBUTES: INTENT(OUT)
573 !
574 ! Direct_Reflectivity_TL: Tangent linear of direct reflectivity
575 ! UNITS: N/A
576 ! TYPE: REAL(fp)
577 ! DIMENSION: Scalar
578 ! ATTRIBUTES: INTENT(OUT)
579 !
580 ! Planck_Surface_TL: Tangent linear of Planck Surface
581 ! UNITS: N/A
582 ! TYPE: REAL(fp)
583 ! DIMENSION: Scalar
584 ! ATTRIBUTES: INTENT(OUT)
585 !
586 ! Planck_Atmosphere_TL: Tangent linear of Planck Atmosphere
587 ! UNITS: N/A
588 ! TYPE: REAL(fp)
589 ! DIMENSION: Rank-1
590 ! ATTRIBUTES: INTENT(OUT)
591 !
592 ! Pff_TL: Tangent linear of forward scattering
593 ! phase function
594 ! UNITS: N/A
595 ! TYPE: REAL(fp)
596 ! DIMENSION: Rank-3
597 ! ATTRIBUTES: INTENT(OUT)
598 !
599 ! Pbb_TL: Tangent linear of backward scattering
600 ! phase function
601 ! UNITS: N/A
602 ! TYPE: REAL(fp)
603 ! DIMENSION: Rank-3
604 ! ATTRIBUTES: INTENT(OUT)
605 !
606 ! FUNCTION RESULT:
607 ! Error_Status: The return value is an integer defining the error status
608 ! The error codes are defined in the Message_Handler module.
609 ! If == SUCCESS the computation was sucessful
610 ! == FAILURE an unrecoverable error occurred
611 ! UNITS: N/A
612 ! TYPE: INTEGER
613 ! DIMENSION: Scalar
614 !
615 ! COMMENTS:
616 ! Note the INTENT on the output RTSolution_TL argument is IN OUT rather
617 ! than just OUT. This is necessary because the argument may be defined
618 ! upon input. To prevent memory leaks, the IN OUT INTENT is a must.
619 !
620 !--------------------------------------------------------------------------------
621  FUNCTION assign_common_input_tl( &
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
634  nz , & ! Output
635  User_Emissivity_TL , & ! Output
636  Direct_Reflectivity_TL, & ! Output
637  Planck_Surface_TL , & ! Output
638  Planck_Atmosphere_TL , & ! Output
639  Pff_TL , & ! Output
640  Pbb_TL , & ! Output
641  RTV ) & ! Internal variable input
642  result( error_status )
643  ! Arguments
644  TYPE(crtm_atmosphere_type), INTENT(IN) :: atmosphere
645  TYPE(crtm_surface_type), INTENT(IN) :: surface
646  TYPE(crtm_atmoptics_type), INTENT(IN) :: atmoptics
647  TYPE(crtm_sfcoptics_type), INTENT(IN) :: sfcoptics
648  TYPE(crtm_atmosphere_type), INTENT(IN) :: atmosphere_tl
649  TYPE(crtm_surface_type), INTENT(IN) :: surface_tl
650  TYPE(crtm_atmoptics_type), INTENT(IN) :: atmoptics_tl
651  TYPE(crtm_sfcoptics_type), INTENT(IN OUT) :: sfcoptics_tl
652  TYPE(crtm_geometryinfo_type), INTENT(IN) :: geometryinfo
653  INTEGER, INTENT(IN) :: sensorindex
654  INTEGER, INTENT(IN) :: channelindex
655  TYPE(crtm_rtsolution_type), INTENT(IN OUT) :: rtsolution_tl
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(:,:,:)
663  TYPE(rtv_type), INTENT(IN) :: rtv
664  ! Function result
665  INTEGER :: error_status
666  ! Local parameters
667  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Input_TL'
668  ! Local variables
669  CHARACTER(256) :: message
670  INTEGER :: i, k
671  INTEGER :: no, na, nt
672 
673  ! ------
674  ! Set up
675  ! ------
676  error_status = success
677 
678  nz = rtv%n_Angles
679  ! Set the sensor zenith angle index
680  sfcoptics_tl%Index_Sat_Ang = sfcoptics%Index_Sat_Ang
681  ! Save the TL optical depth if required
682  IF ( crtm_rtsolution_associated( rtsolution_tl ) ) THEN
683  ! Shorter names for indexing
684  no = rtsolution_tl%n_Layers ! Original no. of layers
685  na = rtv%n_Added_Layers ! No. of added layers
686  nt = rtv%n_Layers ! Current total no. of layers
687  ! Assign only the optical depth profile
688  ! defined by the user input layering
689  rtsolution_tl%Layer_Optical_Depth(1:no) = atmoptics_tl%Optical_Depth(na+1:nt)
690  END IF
691 
692  ! ---------------------------------------------------
693  ! Compute the tangent-linear phase matrix if required
694  ! ---------------------------------------------------
695  IF ( rtv%Scattering_RT ) THEN
696  CALL crtm_phase_matrix_tl( &
697  atmoptics, & ! Input
698  atmoptics_tl, & ! Input
699  pff_tl, & ! Output - TL forward scattering
700  pbb_tl, & ! Output - TL backward scattering
701  rtv ) ! Internal variable
702  END IF
703 
704 
705  ! ---------------------------------------------------------------------
706  ! Assign the number of Stokes parameters. Currently, this is set to 1
707  ! for decoupled polarization between surface and atmosphere.
708  ! Remember for polarised microwave instruments, each frequency's Stokes
709  ! parameter is treated as a separate channel.
710  ! ---------------------------------------------------------------------
711  sfcoptics_tl%n_Stokes = 1
712 
713 
714  ! ----------------------------------------------------
715  ! Assign the angles and weights to SfcOptics structure
716  ! ----------------------------------------------------
717  sfcoptics_tl%n_Angles = sfcoptics%n_Angles
718  sfcoptics_tl%Angle = sfcoptics%Angle
719  sfcoptics_tl%Weight = sfcoptics%Weight
720 
721 
722  ! -----------------------------------------
723  ! Populate the SfcOptics_TL structure based
724  ! on FORWARD model SfcOptics Compute_Switch
725  ! -----------------------------------------
726  IF ( sfcoptics%Compute ) THEN
727  error_status = crtm_compute_sfcoptics_tl( &
728  surface , & ! Input
729  sfcoptics , & ! Input
730  surface_tl , & ! Input
731  geometryinfo, & ! Input
732  sensorindex , & ! Input
733  channelindex, & ! Input
734  sfcoptics_tl, & ! In/Output
735  rtv%SOV ) ! Internal variable input
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)
740  CALL display_message( routine_name, trim(message), error_status )
741  RETURN
742  END IF
743  ELSE
744  IF( rtv%Scattering_RT ) THEN
745  ! Replicate the user emissivity for all angles
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
751  ! Replicate the user reflectivities for all angles
752  sfcoptics_tl%Direct_Reflectivity(1:nz,1) = direct_reflectivity_tl
753  IF( rtv%Diffuse_Surface) THEN
754  DO i = 1, nz
755  sfcoptics_tl%Reflectivity(1:nz, 1, i, 1) = -sfcoptics_tl%Emissivity(i,1)*sfcoptics%Weight(i)
756  END DO
757  ELSE ! Specular surface
758  DO i = 1, nz
759  sfcoptics_tl%Reflectivity(i, 1, i, 1) = -sfcoptics_tl%Emissivity(i,1)
760  END DO
761  END IF
762  ELSE
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
766  END IF
767  END IF
768 
769  ! ------------------------------
770  ! Save the TL surface properties
771  ! ------------------------------
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 )
775 
776  ! -------------------------------------------
777  ! Compute the tangent-linear planck radiances
778  ! -------------------------------------------
779  IF( rtv%mth_Azi == 0 ) THEN
780  ! Atmospheric layer TL radiances
781  DO k = 1, atmosphere%n_Layers
783  sensorindex , & ! Input
784  channelindex , & ! Input
785  atmosphere%Temperature(k) , & ! Input
786  atmosphere_tl%Temperature(k), & ! Input
787  planck_atmosphere_tl(k) ) ! Output
788  END DO
789  ! Surface TL radiance
791  sensorindex , & ! Input
792  channelindex , & ! Input
793  sfcoptics%Surface_Temperature , & ! Input
794  sfcoptics_tl%Surface_Temperature, & ! Input
795  planck_surface_tl ) ! Output
796  ELSE
797  planck_atmosphere_tl = zero
798  planck_surface_tl = zero
799  END IF
800 
801  END FUNCTION assign_common_input_tl
802 
803 !--------------------------------------------------------------------------------
804 !
805 ! NAME:
806 ! Assign_Common_Input_AD
807 !
808 ! PURPOSE:
809 ! Function to assign common input before calling the adjoint RT algorithms.
810 !
811 ! CALLING SEQUENCE:
812 ! Error_Status = Assign_Common_Input_AD( SfcOptics , & ! FWD Input
813 ! RTSolution , & ! FWD Input
814 ! GeometryInfo , & ! Input
815 ! SensorIndex , & ! Input
816 ! ChannelIndex , & ! Input
817 ! RTSolution_AD, & ! AD Input/Output
818 ! Atmosphere_AD, & ! AD Output
819 ! Surface_AD , & ! AD Output
820 ! AtmOptics_AD , & ! AD Output
821 ! SfcOptics_AD , & ! AD Output
822 ! RTV ) ! Internal variable input
823 !
824 ! INPUT ARGUMENTS:
825 !
826 ! SfcOptics: Structure containing the surface optical properties
827 ! data.
828 ! UNITS: N/A
829 ! TYPE: CRTM_SfcOptics_type
830 ! DIMENSION: Scalar
831 ! ATTRIBUTES: INTENT(IN)
832 !
833 ! RTSolution: Structure containing the solution to the RT equation
834 ! for the given inputs.
835 ! UNITS: N/A
836 ! TYPE: CRTM_RTSolution_type
837 ! DIMENSION: Scalar
838 ! ATTRIBUTES: INTENT(IN)
839 !
840 ! GeometryInfo: Structure containing the view geometry data.
841 ! UNITS: N/A
842 ! TYPE: CRTM_GeometryInfo_type
843 ! DIMENSION: Scalar
844 ! ATTRIBUTES: INTENT(IN)
845 !
846 ! SensorIndex: Sensor index id. This is a unique index associated
847 ! with a (supported) sensor used to access the
848 ! shared coefficient data for a particular sensor.
849 ! See the ChannelIndex argument.
850 ! UNITS: N/A
851 ! TYPE: INTEGER
852 ! DIMENSION: Scalar
853 ! ATTRIBUTES: INTENT(IN)
854 !
855 ! ChannelIndex: Channel index id. This is a unique index associated
856 ! with a (supported) sensor channel used to access the
857 ! shared coefficient data for a particular sensor's
858 ! channel.
859 ! See the SensorIndex argument.
860 ! UNITS: N/A
861 ! TYPE: INTEGER
862 ! DIMENSION: Scalar
863 ! ATTRIBUTES: INTENT(IN)
864 !
865 ! RTV: Structure containing internal forward model variables
866 ! required for subsequent tangent-linear or adjoint model
867 ! calls. The contents of this structure are NOT accessible
868 ! outside of the CRTM_RTSolution module.
869 ! UNITS: N/A
870 ! TYPE: RTV_type
871 ! DIMENSION: Scalar
872 ! ATTRIBUTES: INTENT(IN)
873 !
874 ! OUTPUT ARGUMENTS:
875 !
876 ! RTSolution_AD: Structure containing the RT solution adjoint inputs.
877 ! UNITS: N/A
878 ! TYPE: CRTM_RTSolution_type
879 ! DIMENSION: Scalar
880 ! ATTRIBUTES: INTENT(IN OUT)
881 !
882 ! SfcOptics_AD: Structure containing the adjoint surface optical
883 ! properties data.
884 ! UNITS: N/A
885 ! TYPE: CRTM_SfcOptics_type
886 ! DIMENSION: Scalar
887 ! ATTRIBUTES: INTENT(IN OUT)
888 !
889 ! Planck_Surface_AD: Adjoint of Planck surface.
890 ! UNITS: N/A
891 ! TYPE: REAL(fp)
892 ! DIMENSION: Scalar
893 ! ATTRIBUTES: INTENT(OUT)
894 !
895 ! Planck_Atmosphere_AD: Adjoint of Planck atmosphere
896 ! UNITS: N/A
897 ! TYPE: REAL(fp)
898 ! DIMENSION: Rank-1 (0:n_Layers)
899 ! ATTRIBUTES: INTENT(OUT)
900 !
901 ! Radiance_AD: Adjoint of Radiance
902 ! UNITS: N/A
903 ! TYPE: REAL(fp)
904 ! DIMENSION: Scalar
905 ! ATTRIBUTES: INTENT(OUT)
906 !
907 ! nz: Adjoint of Planck atmosphere
908 ! UNITS: N/A
909 ! TYPE: INTEGER
910 ! DIMENSION: Scalar
911 ! ATTRIBUTES: INTENT(OUT)
912 !
913 ! FUNCTION RESULT:
914 ! Error_Status: The return value is an integer defining the error status
915 ! The error codes are defined in the Message_Handler module.
916 ! If == SUCCESS the computation was sucessful
917 ! == FAILURE an unrecoverable error occurred
918 ! UNITS: N/A
919 ! TYPE: INTEGER
920 ! DIMENSION: Scalar
921 !
922 ! COMMENTS:
923 ! Note the INTENT on all of the adjoint arguments (whether input or output)
924 ! is IN OUT rather than just OUT. This is necessary because the Input
925 ! adjoint arguments are modified, and the Output adjoint arguments must
926 ! be defined prior to entry to this routine. So, anytime a structure is
927 ! to be output, to prevent memory leaks the IN OUT INTENT is a must.
928 !
929 !--------------------------------------------------------------------------------
930 
931  FUNCTION assign_common_input_ad( &
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
942  nz , & ! Output
943  RTV ) & ! Internal variable input
944  result( error_status )
945  ! Arguments
946  TYPE(crtm_sfcoptics_type), INTENT(IN) :: sfcoptics
947  TYPE(crtm_rtsolution_type), INTENT(IN) :: rtsolution
948  TYPE(crtm_geometryinfo_type), INTENT(IN) :: geometryinfo
949  INTEGER, INTENT(IN) :: sensorindex
950  INTEGER, INTENT(IN) :: channelindex
951  TYPE(crtm_rtsolution_type), INTENT(IN OUT) :: rtsolution_ad
952  TYPE(crtm_sfcoptics_type), INTENT(IN OUT) :: sfcoptics_ad
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
957  TYPE(rtv_type), INTENT(IN) :: rtv
958  ! Function result
959  INTEGER :: error_status
960  ! Local parameters
961  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Input_AD'
962 
963  ! -----
964  ! Setup
965  ! -----
966  error_status = success
967  ! Short named local copies
968  nz = rtv%n_Angles
969  ! Set the sensor zenith angle index
970  sfcoptics_ad%Index_Sat_Ang = sfcoptics%Index_Sat_Ang
971 
972  ! Initialise local adjoint variables
973  planck_surface_ad = zero
974  planck_atmosphere_ad = zero
975  radiance_ad = zero
976 
977  ! ------------------------------------------
978  ! Compute the brightness temperature adjoint
979  ! ------------------------------------------
980  IF ( spccoeff_isinfraredsensor( sc(sensorindex) ) .OR. &
981  spccoeff_ismicrowavesensor( sc(sensorindex) ) ) THEN
982  IF( rtv%mth_Azi == 0 ) THEN
984  sensorindex , & ! Input
985  channelindex , & ! Input
986  rtsolution%Radiance , & ! Input
987  rtsolution_ad%Brightness_Temperature, & ! Input
988  radiance_ad ) ! Output
989  rtsolution_ad%Brightness_Temperature = zero
990  END IF
991  END IF
992 
993  ! accumulate Fourier component
994  radiance_ad = radiance_ad + rtsolution_ad%Radiance * &
995  cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
996 
997  END FUNCTION assign_common_input_ad
998 
999 !-------------------------------------------------------------------------
1000 !
1001 ! NAME:
1002 ! Assign_Common_Output
1003 !
1004 ! PURPOSE:
1005 ! Function to assign output from CRTM_ADA, CRTM_SOI and CRTM_Emission
1006 !
1007 ! CALLING SEQUENCE:
1008 ! Error_Status = Assign_Common_Output( Atmosphere , & ! Input
1009 ! SfcOptics , & ! Input
1010 ! GeometryInfo , & ! Input
1011 ! SensorIndex , & ! Input
1012 ! ChannelIndex , & ! Input
1013 ! RTV , & ! Input
1014 ! RTSolution ) ! Output
1015 !
1016 ! INPUT ARGUMENTS:
1017 !
1018 ! Atmosphere: Structure containing the atmospheric state data.
1019 ! UNITS: N/A
1020 ! TYPE: CRTM_Atmosphere_type
1021 ! DIMENSION: Scalar
1022 ! ATTRIBUTES: INTENT(IN)
1023 !
1024 ! SfcOptics: Structure containing the surface optical properties
1025 ! data. Argument is defined as INTENT (IN OUT ) as
1026 ! different RT algorithms may compute the surface
1027 ! optics properties before this routine is called.
1028 ! UNITS: N/A
1029 ! TYPE: CRTM_SfcOptics_type
1030 ! DIMENSION: Scalar
1031 ! ATTRIBUTES: INTENT(IN)
1032 !
1033 ! GeometryInfo: Structure containing the view geometry data.
1034 ! UNITS: N/A
1035 ! TYPE: CRTM_GeometryInfo_type
1036 ! DIMENSION: Scalar
1037 ! ATTRIBUTES: INTENT(IN)
1038 !
1039 ! SensorIndex: Sensor index id. This is a unique index associated
1040 ! with a (supported) sensor used to access the
1041 ! shared coefficient data for a particular sensor.
1042 ! See the ChannelIndex argument.
1043 ! UNITS: N/A
1044 ! TYPE: INTEGER
1045 ! DIMENSION: Scalar
1046 ! ATTRIBUTES: INTENT(IN)
1047 !
1048 ! ChannelIndex: Channel index id. This is a unique index associated
1049 ! with a (supported) sensor channel used to access the
1050 ! shared coefficient data for a particular sensor's
1051 ! channel.
1052 ! See the SensorIndex argument.
1053 ! UNITS: N/A
1054 ! TYPE: INTEGER
1055 ! DIMENSION: Scalar
1056 ! ATTRIBUTES: INTENT(IN)
1057 !
1058 ! RTV: Structure containing internal variables required for
1059 ! subsequent tangent-linear or adjoint model calls.
1060 ! The contents of this structure are NOT accessible
1061 ! outside of the CRTM_RTSolution module.
1062 ! UNITS: N/A
1063 ! TYPE: RTV_type
1064 ! DIMENSION: Scalar
1065 ! ATTRIBUTES: INTENT(IN)
1066 !
1067 ! OUTPUT ARGUMENTS:
1068 ! RTSolution: Structure containing the soluition to the RT equation
1069 ! for the given inputs.
1070 ! UNITS: N/A
1071 ! TYPE: CRTM_RTSolution_type
1072 ! DIMENSION: Scalar
1073 ! ATTRIBUTES: INTENT(IN OUT)
1074 !
1075 ! FUNCTION RESULT:
1076 ! Error_Status: The return value is an integer defining the error status.
1077 ! The error codes are defined in the Message_Handler module.
1078 ! If == SUCCESS the computation was sucessful
1079 ! == FAILURE an unrecoverable error occurred
1080 ! UNITS: N/A
1081 ! TYPE: INTEGER
1082 ! DIMENSION: Scalar
1083 !
1084 ! COMMENTS:
1085 ! Note the INTENT on the output RTSolution argument is IN OUT rather than
1086 ! just OUT. This is necessary because the argument is defined upon
1087 ! input. To prevent memory leaks, the IN OUT INTENT is a must.
1088 !
1089 !--------------------------------------------------------------------------------
1090 
1091  FUNCTION assign_common_output( &
1092  Atmosphere , & ! Input
1093  SfcOptics , & ! Input
1094  GeometryInfo , & ! Input
1095  SensorIndex , & ! Input
1096  ChannelIndex , & ! Input
1097  RTV , & ! Input
1098  RTSolution ) & ! Output
1099  result( error_status )
1100  ! Arguments
1101  TYPE(crtm_atmosphere_type) , INTENT(IN) :: atmosphere
1102  TYPE(crtm_sfcoptics_type) , INTENT(IN) :: sfcoptics
1103  TYPE(crtm_geometryinfo_type), INTENT(IN) :: geometryinfo
1104  INTEGER , INTENT(IN) :: sensorindex
1105  INTEGER , INTENT(IN) :: channelindex
1106  TYPE(rtv_type) , INTENT(IN) :: rtv
1107  TYPE(crtm_rtsolution_type) , INTENT(IN OUT) :: rtsolution
1108 
1109  ! Function Result
1110  INTEGER :: error_status
1111  ! Local Parameters
1112  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Output'
1113  ! Local variables
1114  INTEGER :: no, na, nt
1115  REAL(fp) :: radiance
1116 
1117  error_status = success
1118 
1119  ! ADA and SOI specific assignments
1120  IF( rtv%Scattering_RT ) THEN
1121 
1122  ! Assign radiance output from ADA or SOI
1123  IF ( rtv%aircraft%rt ) THEN
1124  radiance = rtv%s_Level_Rad_UP(sfcoptics%Index_Sat_Ang, rtv%aircraft%idx)
1125  ELSE
1126  radiance = rtv%s_Level_Rad_UP(sfcoptics%Index_Sat_Ang, 0)
1127  END IF
1128 
1129  ! Emission specific assignments
1130  ELSE
1131 
1132  ! The output radiance at TOA or at Aircraft height
1133  IF ( rtv%aircraft%rt ) THEN
1134  radiance = rtv%e_Level_Rad_UP(rtv%aircraft%idx)
1135  ELSE
1136  radiance = rtv%e_Level_Rad_UP(0)
1137  END IF
1138 
1139  ! Other emission-only output
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
1144  IF ( crtm_rtsolution_associated( rtsolution ) ) THEN
1145  ! Shorter names for indexing
1146  no = rtsolution%n_Layers ! Original no. of layers
1147  na = rtv%n_Added_Layers ! No. of added layers
1148  nt = rtv%n_Layers ! Current total no. of layers
1149  ! Assign only the upwelling radiance profile
1150  ! defined by the user input layering
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)
1153  END IF
1154  END IF
1155 
1156  ! accumulate Fourier component
1157  rtsolution%Radiance = rtsolution%Radiance + radiance* &
1158  cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1159 
1160  ! ------------------------------------------------
1161  ! Compute the corresponding brightness temperature
1162  ! ------------------------------------------------
1163  IF( rtv%mth_Azi == 0 ) &
1164  CALL crtm_planck_temperature( &
1165  sensorindex, & ! Input
1166  channelindex, & ! Input
1167  rtsolution%Radiance, & ! Input
1168  rtsolution%Brightness_Temperature ) ! Output
1169 
1170  END FUNCTION assign_common_output
1171 !-------------------------------------------------------------------------
1172 !
1173 ! NAME:
1174 ! Assign_Common_Output_TL
1175 !
1176 ! PURPOSE:
1177 ! Function to assign output from CRTM_ADA, CRTM_SOI and CRTM_Emission
1178 !
1179 ! CALLING SEQUENCE:
1180 ! Error_Status = Assign_Common_Output_TL( SfcOptics , & ! Input
1181 ! RTSolution , & ! Input
1182 ! GeometryInfo , & ! Input
1183 ! Radiance_TL , & ! Input
1184 ! Scattering_Radiance_TL , & ! Input
1185 ! SensorIndex , & ! Input
1186 ! ChannelIndex , & ! Input
1187 ! RTV , & ! Input
1188 ! RTSolution_TL ) ! Output
1189 !
1190 ! INPUT ARGUMENTS:
1191 !
1192 ! SfcOptics: Structure containing the surface optical properties
1193 ! data.
1194 ! UNITS: N/A
1195 ! TYPE: CRTM_SfcOptics_type
1196 ! DIMENSION: Scalar
1197 ! ATTRIBUTES: INTENT(IN)
1198 !
1199 ! RTSolution: Structure containing the soluition to the RT equation
1200 ! UNITS: N/A
1201 ! TYPE: CRTM_RTSolution_type
1202 ! DIMENSION: Scalar
1203 ! ATTRIBUTES: INTENT(IN)
1204 !
1205 ! GeometryInfo: Structure containing the view geometry data.
1206 ! UNITS: N/A
1207 ! TYPE: CRTM_GeometryInfo_type
1208 ! DIMENSION: Scalar
1209 ! ATTRIBUTES: INTENT(IN)
1210 !
1211 ! Radiance_TL: TL radiance
1212 ! UNITS: N/A
1213 ! TYPE: REAL(fp)
1214 ! DIMENSION: Scalar
1215 ! ATTRIBUTES: INTENT(IN)
1216 !
1217 ! Scattering_Radiance_TL: Scattering TL radiance
1218 ! UNITS: N/A
1219 ! TYPE: REAL(fp)
1220 ! DIMENSION: Rank-1
1221 ! ATTRIBUTES: INTENT(IN)
1222 !
1223 ! SensorIndex: Sensor index id. This is a unique index associated
1224 ! with a (supported) sensor used to access the
1225 ! shared coefficient data for a particular sensor.
1226 ! See the ChannelIndex argument.
1227 ! UNITS: N/A
1228 ! TYPE: INTEGER
1229 ! DIMENSION: Scalar
1230 ! ATTRIBUTES: INTENT(IN)
1231 !
1232 ! ChannelIndex: Channel index id. This is a unique index associated
1233 ! with a (supported) sensor channel used to access the
1234 ! shared coefficient data for a particular sensor's
1235 ! channel.
1236 ! See the SensorIndex argument.
1237 ! UNITS: N/A
1238 ! TYPE: INTEGER
1239 ! DIMENSION: Scalar
1240 ! ATTRIBUTES: INTENT(IN)
1241 !
1242 ! RTV: Structure containing internal variables required for
1243 ! subsequent tangent-linear or adjoint model calls.
1244 ! The contents of this structure are NOT accessible
1245 ! outside of the CRTM_RTSolution module.
1246 ! UNITS: N/A
1247 ! TYPE: RTV_type
1248 ! DIMENSION: Scalar
1249 ! ATTRIBUTES: INTENT(IN)
1250 !
1251 ! OUTPUT ARGUMENTS:
1252 ! RTSolution_TL: Structure containing the TL soluition to the RT equation.
1253 ! UNITS: N/A
1254 ! TYPE: CRTM_RTSolution_type
1255 ! DIMENSION: Scalar
1256 ! ATTRIBUTES: INTENT(IN OUT)
1257 !
1258 ! FUNCTION RESULT:
1259 ! Error_Status: The return value is an integer defining the error status.
1260 ! The error codes are defined in the Message_Handler module.
1261 ! If == SUCCESS the computation was sucessful
1262 ! == FAILURE an unrecoverable error occurred
1263 ! UNITS: N/A
1264 ! TYPE: INTEGER
1265 ! DIMENSION: Scalar
1266 !
1267 ! COMMENTS:
1268 ! Note the INTENT on the output RTSolution argument is IN OUT rather than
1269 ! just OUT. This is necessary because the argument is defined upon
1270 ! input. To prevent memory leaks, the IN OUT INTENT is a must.
1271 !
1272 !--------------------------------------------------------------------------------
1273 
1274  FUNCTION assign_common_output_tl( &
1275  SfcOptics , & ! Input
1276  RTSolution , & ! Input
1277  GeometryInfo , & ! Input
1278  Radiance_TL , & ! Input
1279  Scattering_Radiance_TL , & ! Input
1280  SensorIndex , & ! Input
1281  ChannelIndex , & ! Input
1282  RTV , & ! Input
1283  RTSolution_TL ) & ! Output
1284  result( error_status )
1285  ! Arguments
1286  TYPE(crtm_sfcoptics_type) , INTENT(IN) :: sfcoptics
1287  TYPE(crtm_rtsolution_type) , INTENT(IN) :: rtsolution
1288  TYPE(crtm_geometryinfo_type), INTENT(IN) :: geometryinfo
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
1293  TYPE(rtv_type) , INTENT(IN) :: rtv
1294  TYPE(crtm_rtsolution_type) , INTENT(IN OUT) :: rtsolution_tl
1295 
1296  ! Function Result
1297  INTEGER :: error_status
1298  ! Local Parameters
1299  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Output_TL'
1300 
1301  error_status = success
1302 
1303  ! ADA and SOI specific assignments
1304  IF( rtv%Scattering_RT ) THEN
1305  ! accumulate Fourier component
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) )
1308  ! Emission specific assignments
1309  ELSE
1310  ! accumulate Fourier component
1311  rtsolution_tl%Radiance = rtsolution_tl%Radiance + radiance_tl* &
1312  cos( rtv%mth_Azi*(geometryinfo%Sensor_Azimuth_Radian-geometryinfo%Source_Azimuth_Radian) )
1313  END IF
1314 
1315  ! ---------------------------------------------------------------
1316  ! Compute the corresponding tangent-linear brightness temperature
1317  ! ---------------------------------------------------------------
1318  IF( rtv%mth_Azi == 0 ) &
1320  sensorindex , & ! Input
1321  channelindex , & ! Input
1322  rtsolution%Radiance , & ! Input
1323  rtsolution_tl%Radiance , & ! Input
1324  rtsolution_tl%Brightness_Temperature ) ! Output
1325 
1326  END FUNCTION assign_common_output_tl
1327 
1328 !------------------------------------------------------------------------------
1329 !
1330 ! NAME:
1331 ! Assign_Common_Output_AD
1332 !
1333 ! PURPOSE:
1334 ! Function to assign output from CRTM_ADA, CRTM_SOI and CRTM_Emission
1335 ! adjoint calls
1336 !
1337 ! CALLING SEQUENCE:
1338 ! Error_Status = Assign_Common_Output_AD( Atmosphere , & ! Input
1339 ! Surface , & ! Input
1340 ! AtmOptics , & ! Input
1341 ! SfcOptics , & ! Input
1342 ! Pff_AD , & ! Input
1343 ! Pbb_AD , & ! Input
1344 ! RTSolution , & ! Input
1345 ! GeometryInfo , & ! Input
1346 ! SensorIndex , & ! Input
1347 ! ChannelIndex , & ! Input
1348 ! nZ , & ! Input
1349 ! AtmOptics_AD , & ! Output
1350 ! SfcOptics_AD , & ! Output
1351 ! Planck_Surface_AD , & ! Output
1352 ! Planck_Atmosphere_AD , & ! Output
1353 ! User_Emissivity_AD , & ! Output
1354 ! Atmosphere_AD , & ! Output
1355 ! Surface_AD , & ! Output
1356 ! RTSolution_AD , & ! Output
1357 ! RTV ) & ! Input
1358 !
1359 ! INPUT ARGUMENTS:
1360 !
1361 ! Atmosphere: Structure containing the atmospheric state data.
1362 ! UNITS: N/A
1363 ! TYPE: CRTM_Atmosphere_type
1364 ! DIMENSION: Scalar
1365 ! ATTRIBUTES: INTENT(IN)
1366 !
1367 ! Surface: Structure containing the Surface state data.
1368 ! UNITS: N/A
1369 ! TYPE: CRTM_Surface_type
1370 ! DIMENSION: Scalar
1371 ! ATTRIBUTES: INTENT(IN)
1372 !
1373 ! AtmOptics: Structure containing the combined atmospheric
1374 ! optical properties for gaseous absorption, clouds,
1375 ! and aerosols.
1376 ! UNITS: N/A
1377 ! TYPE: CRTM_AtmOptics_type
1378 ! DIMENSION: Scalar
1379 ! ATTRIBUTES: INTENT(IN)
1380 !
1381 ! SfcOptics: Structure containing the surface
1382 ! optical properties data.
1383 ! UNITS: N/A
1384 ! TYPE: CRTM_SfcOptics_type
1385 ! DIMENSION: Scalar
1386 ! ATTRIBUTES: INTENT(IN)
1387 !
1388 ! GeometryInfo: Structure containing the view geometry data.
1389 ! UNITS: N/A
1390 ! TYPE: CRTM_GeometryInfo_type
1391 ! DIMENSION: Scalar
1392 ! ATTRIBUTES: INTENT(IN)
1393 !
1394 ! SensorIndex: Sensor index id. This is a unique index associated
1395 ! with a (supported) sensor used to access the
1396 ! shared coefficient data for a particular sensor.
1397 ! See the ChannelIndex argument.
1398 ! UNITS: N/A
1399 ! TYPE: INTEGER
1400 ! DIMENSION: Scalar
1401 ! ATTRIBUTES: INTENT(IN)
1402 !
1403 ! ChannelIndex: Channel index id. This is a unique index associated
1404 ! with a (supported) sensor channel used to access the
1405 ! shared coefficient data for a particular sensor's
1406 ! channel.
1407 ! See the SensorIndex argument.
1408 ! UNITS: N/A
1409 ! TYPE: INTEGER
1410 ! DIMENSION: Scalar
1411 ! ATTRIBUTES: INTENT(IN)
1412 !
1413 ! nz: Integer assigned from RTV%n_Angles
1414 ! UNITS: N/A
1415 ! TYPE: INTEGER
1416 ! DIMENSION: Scalar
1417 ! ATTRIBUTES: INTENT(IN)
1418 !
1419 ! RTV: Structure containing internal variables required for
1420 ! subsequent tangent-linear or adjoint model calls.
1421 ! The contents of this structure are NOT accessible
1422 ! outside of the CRTM_RTSolution module.
1423 ! UNITS: N/A
1424 ! TYPE: RTV_type
1425 ! DIMENSION: Scalar
1426 ! ATTRIBUTES: INTENT(IN)
1427 !
1428 ! OUTPUT ARGUMENTS
1429 !
1430 ! Pff_AD: Forward scattering AD phase matrix
1431 ! UNITS: N/A
1432 ! TYPE: REAL(fp)
1433 ! DIMENSION: Scalar
1434 ! ATTRIBUTES: INTENT(IN OUT)
1435 !
1436 ! Pbb_AD: Backward scattering AD phase matrix
1437 ! UNITS: N/A
1438 ! TYPE: REAL(fp)
1439 ! DIMENSION: Scalar
1440 ! ATTRIBUTES: INTENT(IN OUT)
1441 !
1442 ! AtmOptics_AD: Structure containing the adjoint combined atmospheric
1443 ! optical properties for gaseous absorption, clouds,
1444 ! and aerosols.
1445 ! UNITS: N/A
1446 ! TYPE: CRTM_AtmOptics_type
1447 ! DIMENSION: Scalar
1448 ! ATTRIBUTES: INTENT(IN OUT)
1449 !
1450 ! SfcOptics_AD: Structure containing the adjoint surface optical
1451 ! properties data.
1452 ! UNITS: N/A
1453 ! TYPE: CRTM_SfcOptics_type
1454 ! DIMENSION: Scalar
1455 ! ATTRIBUTES: INTENT(IN OUT)
1456 !
1457 ! Planck_Surface_AD: Adjoint of Planck surface.
1458 ! UNITS: N/A
1459 ! TYPE: REAL(fp)
1460 ! DIMENSION: Scalar
1461 ! ATTRIBUTES: INTENT(IN OUT)
1462 !
1463 ! Planck_Atmosphere_AD: Adjoint of Planck atmosphere
1464 ! UNITS: N/A
1465 ! TYPE: REAL(fp)
1466 ! DIMENSION: Rank-1 (0:n_Layers)
1467 ! ATTRIBUTES: INTENT(IN OUT)
1468 !
1469 ! User_Emissivity_AD: Adjoint of user emissivity
1470 ! UNITS: N/A
1471 ! TYPE: REAL(fp)
1472 ! DIMENSION: Scalar
1473 ! ATTRIBUTES: INTENT(OUT)
1474 !
1475 ! Atmosphere_AD: Adjoint of atmosphere
1476 ! UNITS: N/A
1477 ! TYPE: CRTM_Atmosphere_type
1478 ! DIMENSION: Scalar
1479 ! ATTRIBUTES: INTENT(IN OUT)
1480 !
1481 ! Surface_AD: Adjoint of surface
1482 ! UNITS: N/A
1483 ! TYPE: CRTM_Surface_type
1484 ! DIMENSION: Scalar
1485 ! ATTRIBUTES: INTENT(IN OUT)
1486 !
1487 ! RTSolution_AD: Structure containing the RT solution adjoint inputs.
1488 ! UNITS: N/A
1489 ! TYPE: CRTM_RTSolution_type
1490 ! DIMENSION: Scalar
1491 ! ATTRIBUTES: INTENT(IN OUT)
1492 !
1493 ! FUNCTION RESULT:
1494 ! Error_Status: The return value is an integer defining the error status.
1495 ! The error codes are defined in the Message_Handler module.
1496 ! If == SUCCESS the computation was sucessful
1497 ! == FAILURE an unrecoverable error occurred
1498 ! UNITS: N/A
1499 ! TYPE: INTEGER
1500 ! DIMENSION: Scalar
1501 !
1502 ! COMMENTS:
1503 ! Note the INTENT on the output RTSolution argument is IN OUT rather than
1504 ! just OUT. This is necessary because the argument is defined upon
1505 ! input. To prevent memory leaks, the IN OUT INTENT is a must.
1506 !
1507 !--------------------------------------------------------------------------------
1508 
1509  FUNCTION assign_common_output_ad( &
1510  Atmosphere , & ! Input
1511  Surface , & ! Input
1512  AtmOptics , & ! Input
1513  SfcOptics , & ! Input
1514  Pff_AD , & ! Input
1515  Pbb_AD , & ! Input
1516  GeometryInfo , & ! Input
1517  SensorIndex , & ! Input
1518  ChannelIndex , & ! Input
1519  nZ , & ! 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
1528  RTV ) & ! Input
1529  result( error_status )
1530  ! Arguments
1531  TYPE(crtm_atmosphere_type) , INTENT(IN) :: atmosphere
1532  TYPE(crtm_surface_type) , INTENT(IN) :: surface
1533  TYPE(crtm_atmoptics_type) , INTENT(IN) :: atmoptics
1534  TYPE(crtm_sfcoptics_type) , INTENT(IN) :: sfcoptics
1535  REAL(fp) , INTENT(IN OUT) :: pff_ad(:,:,:)
1536  REAL(fp) , INTENT(IN OUT) :: pbb_ad(:,:,:)
1537  TYPE(crtm_geometryinfo_type), INTENT(IN) :: geometryinfo
1538  INTEGER , INTENT(IN) :: sensorindex
1539  INTEGER , INTENT(IN) :: channelindex
1540  INTEGER , INTENT(IN) :: nz
1541  TYPE(crtm_atmoptics_type) , INTENT(IN OUT) :: atmoptics_ad
1542  TYPE(crtm_sfcoptics_type) , INTENT(IN OUT) :: sfcoptics_ad
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
1546  TYPE(crtm_atmosphere_type) , INTENT(IN OUT) :: atmosphere_ad
1547  TYPE(crtm_surface_type) , INTENT(IN OUT) :: surface_ad
1548  TYPE(crtm_rtsolution_type) , INTENT(IN OUT) :: rtsolution_ad
1549  TYPE(rtv_type) , INTENT(IN) :: rtv
1550 
1551  ! Function Result
1552  INTEGER :: error_status
1553  ! Local Parameters
1554  CHARACTER(*), PARAMETER :: routine_name = 'Assign_Common_Output_AD'
1555  ! Local variables
1556  CHARACTER(256) :: message
1557  INTEGER :: no, na, nt
1558  INTEGER :: k, i
1559 
1560  error_status = success
1561 
1562  ! ------------------------------------
1563  ! Compute the adjoint planck radiances
1564  ! ------------------------------------
1565  IF( rtv%mth_Azi == 0 ) THEN
1566  ! Surface AD radiance
1567  CALL crtm_planck_radiance_ad( &
1568  sensorindex , & ! Input
1569  channelindex , & ! Input
1570  sfcoptics%Surface_Temperature , & ! Input
1571  planck_surface_ad , & ! Input
1572  sfcoptics_ad%Surface_Temperature ) ! In/Output
1573  planck_surface_ad = zero
1574  ! Atmospheric layer AD radiances
1575  DO k = 1, atmosphere%n_Layers
1576  CALL crtm_planck_radiance_ad( &
1577  sensorindex , & ! Input
1578  channelindex , & ! Input
1579  atmosphere%Temperature(k) , & ! Input
1580  planck_atmosphere_ad(k) , & ! Input
1581  atmosphere_ad%Temperature(k) ) ! In/Output
1582  planck_atmosphere_ad(k) = zero
1583  END DO
1584  ELSE
1585  planck_surface_ad = zero
1586  planck_atmosphere_ad = zero
1587  END IF
1588 
1589 
1590  ! ---------------------------------------------------------------------
1591  ! Assign the number of Stokes parameters. Currently, this is set to 1
1592  ! for decoupled polarization between surface and atmosphere.
1593  ! Remember for polarised microwave instruments, each frequency's Stokes
1594  ! parameter is treated as a separate channel.
1595  ! ---------------------------------------------------------------------
1596  sfcoptics_ad%n_Stokes = 1
1597 
1598 
1599  ! ----------------------------------------------------
1600  ! Assign the angles and weights to SfcOptics structure
1601  ! ----------------------------------------------------
1602  sfcoptics_ad%n_Angles = sfcoptics%n_Angles
1603  sfcoptics_ad%Angle = sfcoptics%Angle
1604  sfcoptics_ad%Weight = sfcoptics%Weight
1605 
1606 
1607  ! --------------------------------------------------------------------------------
1608  ! This part is designed for specific requirement for the total sensitivity
1609  ! of emissivity. The requirement is regardless whether having user input or not.
1610  ! Therefore, this part is common part with/without optional input.
1611  ! The outputs of CRTM_Compute_RTSolution are the partial derivative
1612  ! of emissivity and reflectivity. Since SfcOptics_AD is used as input only in
1613  ! CRTM_Compute_SfcOptics_AD, the total sensitivity of the emissivity is taken into
1614  ! account for here.
1615  ! --------------------------------------------------------------------------------
1616  IF( rtv%Scattering_RT ) THEN
1617  user_emissivity_ad = zero
1618  IF( rtv%Diffuse_Surface) THEN
1619  DO i = nz, 1, -1
1620  user_emissivity_ad = user_emissivity_ad - &
1621  (sum(sfcoptics_ad%Reflectivity(1:nz,1,i,1))*sfcoptics%Weight(i))
1622  END DO
1623  ELSE ! Specular surface
1624  DO i = nz, 1, -1
1625  user_emissivity_ad = user_emissivity_ad - sfcoptics_ad%Reflectivity(i,1,i,1)
1626  END DO
1627  END IF
1628 ! Direct_Reflectivity_AD = SUM(SfcOptics_AD%Direct_Reflectivity(1:nZ,1))
1629 ! SfcOptics_AD%Direct_Reflectivity(1,1) = SfcOptics_AD%Direct_Reflectivity(1,1) +
1630 ! (Direct_Reflectivity_AD/PI)
1631  rtsolution_ad%Surface_Emissivity = user_emissivity_ad
1632  ELSE
1633  rtsolution_ad%Surface_Emissivity = sfcoptics_ad%Emissivity(sfcoptics_ad%Index_Sat_Ang,1) - &
1634  sfcoptics_ad%Reflectivity(1,1,1,1)
1635  END IF
1636 
1637  ! -----------------------------------------
1638  ! Populate the SfcOptics_AD structure based
1639  ! on FORWARD model SfcOptics Compute_Switch
1640  ! -----------------------------------------
1641  IF ( sfcoptics%Compute ) THEN
1642  error_status = crtm_compute_sfcoptics_ad( &
1643  surface , & ! Input
1644  sfcoptics , & ! Input
1645  sfcoptics_ad, & ! Input
1646  geometryinfo, & ! Input
1647  sensorindex , & ! Input
1648  channelindex, & ! Input
1649  surface_ad , & ! In/Output
1650  rtv%SOV ) ! Internal variable input
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)
1655  CALL display_message( routine_name, trim(message), error_status )
1656  RETURN
1657  END IF
1658  END IF
1659 
1660  ! --------------------------------------------
1661  ! Compute the adjoint phase matrix if required
1662  ! --------------------------------------------
1663  IF ( rtv%Scattering_RT ) THEN
1664  CALL crtm_phase_matrix_ad( &
1665  atmoptics, & ! Input - FWD
1666  pff_ad, & ! Input - AD forward scattering
1667  pbb_ad, & ! Input - AD backward scattering
1668  atmoptics_ad, & ! Output - AD
1669  rtv ) ! Internal variable
1670  END IF
1671 
1672  ! ------------------------------------------
1673  ! Save the adjoint optical depth if required
1674  ! ------------------------------------------
1675  IF ( crtm_rtsolution_associated( rtsolution_ad ) ) THEN
1676  ! Shorter names for indexing
1677  no = rtsolution_ad%n_Layers ! Original no. of layers
1678  na = rtv%n_Added_Layers ! No. of added layers
1679  nt = rtv%n_Layers ! Current total no. of layers
1680  ! Assign only the optical depth profile
1681  ! defined by the user input layering
1682  rtsolution_ad%Layer_Optical_Depth(1:no) = atmoptics_ad%Optical_Depth(na+1:nt)
1683  END IF
1684 
1685  END FUNCTION assign_common_output_ad
1686 
1687 !################################################################################
1688 !################################################################################
1689 !## ##
1690 !## ## PRIVATE MODULE ROUTINES ## ##
1691 !## ##
1692 !################################################################################
1693 !################################################################################
1694 
1695 ! -------------------------------------------------------------------------------
1696 !
1697 ! NAME:
1698 ! CRTM_Phase_Matrix
1699 !
1700 ! PURPOSE:
1701 ! Subroutine to calculate the phase function for the scattering model.
1702 !
1703 ! CALLING SEQUENCE:
1704 ! CALL CRTM_Phase_Matrix( AtmOptics, & ! Input
1705 ! RTV ) ! Internal variable
1706 !
1707 ! INPUT ARGUMENTS:
1708 ! AtmOptics: Structure containing the atmospheric optical
1709 ! parameters
1710 ! UNITS: N/A
1711 ! TYPE: CRTM_AtmOptics_type
1712 ! DIMENSION: Scalar
1713 ! ATTRIBUTES: INTENT(IN)
1714 !
1715 ! RTV: Structure containing internal forward model variables
1716 ! required for tangent-linear or adjoint model calls in
1717 ! CRTM_RTSolution module.
1718 ! UNITS: N/A
1719 ! TYPE: RTV_type
1720 ! DIMENSION: Scalar
1721 ! ATTRIBUTES: INTENT(IN OUT)
1722 !
1723 ! OUTPUT ARGUMENTS:
1724 ! RTV: Structure containing internal forward model variables
1725 ! required for tangent-linear or adjoint model
1726 ! calls in CRTM_RTSolution module.
1727 ! UNITS: N/A
1728 ! TYPE: RTV_type
1729 ! DIMENSION: Scalar
1730 ! ATTRIBUTES: INTENT(IN OUT)
1731 !
1732 ! COMMENTS:
1733 ! Note that the INTENT on the RTV argument is IN OUT as it
1734 ! contains data prior to this call and is filled with data within this
1735 ! routine.
1736 !
1737 !S-
1738 !--------------------------------------------------------------------------------
1739  SUBROUTINE crtm_phase_matrix( &
1740  AtmOptics, & ! Input
1741  RTV ) ! Internal variable
1742  ! Arguments
1743  TYPE(CRTM_AtmOptics_type), INTENT(IN) :: AtmOptics
1744  TYPE(RTV_type), INTENT(IN OUT) :: RTV
1745  ! Local variables
1746  INTEGER :: i, j, k, l, ifac, jn
1747 
1748  !#--------------------------------------------------------------------------#
1749  !# -- COMPUTING LEGENDRE FUNCTION -- #
1750  !# Pplus for positive cosine of angles #
1751  !# Pminus for negative cosine of angles #
1752  !#--------------------------------------------------------------------------#
1753 
1754  DO i = 1, rtv%n_Angles
1755  CALL legendre_m( rtv%mth_Azi , &
1756  atmoptics%n_Legendre_Terms, &
1757  rtv%COS_Angle(i) , &
1758  rtv%Pleg(0:,i) )
1759  END DO
1760 
1761  IF(rtv%Solar_Flag_true) &
1762  CALL legendre_m( rtv%mth_Azi , &
1763  atmoptics%n_Legendre_Terms , &
1764  rtv%COS_SUN , &
1765  rtv%Pleg(0:,rtv%n_Angles+1) )
1766 
1767  IF( rtv%Solar_Flag_true ) THEN
1768  jn = rtv%n_Angles + 1
1769  ELSE
1770  jn = rtv%n_Angles
1771  END IF
1772 
1773  !#--------------------------------------------------------------------------#
1774  !# -- COMPUTE THE PHASE MATRICES -- #
1775  !#--------------------------------------------------------------------------#
1776 
1777  layer_loop: DO k = 1, rtv%n_Layers
1778 
1779 
1780  ! ------------------------------
1781  ! Only proceed if the scattering
1782  ! coefficient is significant
1783  ! ------------------------------
1784 
1785  significant_scattering: IF( atmoptics%Single_Scatter_Albedo(k) > scattering_albedo_threshold) THEN
1786 
1787  DO j = 1, jn
1788  ! add solar angle
1789  DO i = 1, rtv%n_Angles
1790 
1791  rtv%Off(i,j,k)=zero
1792  rtv%Obb(i,j,k)=zero
1793 
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 )
1798  END DO
1799 
1800  rtv%Pff(i,j,k) = rtv%Off(i,j,k)
1801  rtv%Pbb(i,j,k) = rtv%Obb(i,j,k)
1802 
1803  ! For intensity, the phase matrix element must >= ZERO
1804  IF ( rtv%mth_Azi == 0 ) THEN
1805  IF(rtv%Pff(i,j,k) < zero) rtv%Pff(i,j,k) = phase_threshold
1806  IF(rtv%Pbb(i,j,k) < zero) rtv%Pbb(i,j,k) = phase_threshold
1807  END IF
1808 
1809  END DO
1810  END DO
1811 
1812  ! Normalization to ensure energy conservation
1813  IF ( rtv%mth_Azi == 0 ) CALL normalize_phase( k, rtv )
1814 
1815  END IF significant_scattering
1816 
1817  END DO layer_loop
1818 
1819  END SUBROUTINE crtm_phase_matrix
1820 
1821 !--------------------------------------------------------------------------------
1822 !
1823 ! NAME:
1824 ! CRTM_Phase_Matrix_TL
1825 !
1826 ! PURPOSE:
1827 ! Subroutine to calculate the tangent-linear phase function for the
1828 ! scattering model.
1829 !
1830 ! CALLING SEQUENCE:
1831 ! CALL CRTM_Phase_Matrix_TL( AtmOptics, & ! FWD Input
1832 ! AtmOptics_TL, & ! TL Input
1833 ! Pff_TL, & ! TL Output
1834 ! Pff_TL, & ! TL Output
1835 ! RTV ) ! Internal variable
1836 !
1837 ! INPUT ARGUMENTS:
1838 ! AtmOptics: Structure containing the atmospheric optical
1839 ! parameters
1840 ! UNITS: N/A
1841 ! TYPE: CRTM_AtmOptics_type
1842 ! DIMENSION: Scalar
1843 ! ATTRIBUTES: INTENT(IN)
1844 !
1845 ! AtmOptics_TL: Structure containing the tangent-linear atmospheric
1846 ! optical parameters
1847 ! UNITS: N/A
1848 ! TYPE: CRTM_AtmOptics_type
1849 ! DIMENSION: Scalar
1850 ! ATTRIBUTES: INTENT(IN)
1851 !
1852 ! RTV: Structure containing internal forward model variables
1853 ! required for subsequent tangent-linear or adjoint model
1854 ! calls. The contents of this structure are NOT accessible
1855 ! outside of the CRTM_RTSolution module.
1856 ! UNITS: N/A
1857 ! TYPE: RTV_type
1858 ! DIMENSION: Scalar
1859 ! ATTRIBUTES: INTENT(IN)
1860 !
1861 ! OUTPUT ARGUMENTS:
1862 ! Pff_TL: Array containing the tangent-linear of the
1863 ! forward phase matrix.
1864 ! UNITS: N/A
1865 ! TYPE: REAL
1866 ! DIMENSION: Rank-3, n_Angles x n_Angles x n_Layers
1867 ! ATTRIBUTES: INTENT(OUT)
1868 !
1869 ! Pbb_TL: Array containing the tangent-linear of the
1870 ! backward phase matrix.
1871 ! UNITS: N/A
1872 ! TYPE: REAL
1873 ! DIMENSION: Rank-3, n_Angles x n_Angles x n_Layers
1874 ! ATTRIBUTES: INTENT(OUT)
1875 !
1876 !--------------------------------------------------------------------------------
1877 
1878  SUBROUTINE crtm_phase_matrix_tl( &
1879  AtmOptics, & ! FWD Input
1880  AtmOptics_TL, & ! TL Input
1881  Pff_TL, & ! TL Output
1882  Pbb_TL, & ! TL Output
1883  RTV ) ! Internal variable
1884  ! Arguments
1885  TYPE(CRTM_AtmOptics_type) , INTENT(IN) :: AtmOptics
1886  TYPE(CRTM_AtmOptics_type) , INTENT(IN) :: AtmOptics_TL
1887  REAL(fp) , INTENT(OUT) :: Pff_TL(:,:,:) ! n_Angles x n_Angles x n_Layers
1888  REAL(fp) , INTENT(OUT) :: Pbb_TL(:,:,:) ! n_Angles x n_Angles x n_Layers
1889  TYPE(RTV_type) , INTENT(IN) :: RTV
1890  ! Local variables
1891  INTEGER :: i, j, k, l, nZ, ifac, jn
1892  REAL(fp), DIMENSION(RTV%n_Angles,RTV%n_Angles+1) :: Lff, Lbb
1893 
1894 
1895 
1896  !#--------------------------------------------------------------------------#
1897  !# -- COMPUTE THE TANGENT-LINEAR PHASE MATRICES -- #
1898  !#--------------------------------------------------------------------------#
1899 
1900  nz = rtv%n_Angles
1901  IF( rtv%Solar_Flag_true ) THEN
1902  jn = rtv%n_Angles + 1
1903  ELSE
1904  jn = rtv%n_Angles
1905  END IF
1906 
1907  layer_loop: DO k = 1, rtv%n_Layers
1908 
1909  ! ------------------------------
1910  ! Only proceed if the scattering
1911  ! coefficient is significant
1912  ! ------------------------------
1913 
1914  significant_scattering: IF( atmoptics%Single_Scatter_Albedo(k) > scattering_albedo_threshold) THEN
1915 
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)
1918 
1919 
1920  DO j = 1, jn
1921  DO i = 1, rtv%n_Angles
1922 
1923  pff_tl(i,j,k) = zero
1924  pbb_tl(i,j,k) = zero
1925 
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 )
1930  END DO
1931 
1932  ! For intensity, the FWD phase matrix element must >= ZERO
1933  ! so the TL form is always == zero.
1934  IF ( rtv%mth_Azi == 0 ) THEN
1935  IF ( rtv%Off(i,j,k) < zero ) THEN
1936  pff_tl(i,j,k) = zero
1937  lff(i,j) = phase_threshold
1938  END IF
1939 
1940  IF ( rtv%Obb(i,j,k) < zero ) THEN
1941  pbb_tl(i,j,k) = zero
1942  lbb(i,j) = phase_threshold
1943  END IF
1944  END IF
1945 
1946  END DO
1947  END DO
1948 
1949  IF ( rtv%mth_Azi == 0 ) THEN
1950  ! Normalisation for energy conservation
1951  CALL normalize_phase_tl( &
1952  k, rtv, &
1953  lff, & ! FWD Input
1954  lbb, & ! FWD Input
1955  pff_tl(:,:,k), & ! TL Output
1956  pbb_tl(:,:,k) ) ! TL Output
1957  END IF
1958  END IF significant_scattering
1959 
1960  END DO layer_loop
1961 
1962  END SUBROUTINE crtm_phase_matrix_tl
1963 
1964 !--------------------------------------------------------------------------------
1965 !
1966 ! NAME:
1967 ! CRTM_Phase_Matrix_AD
1968 !
1969 ! PURPOSE:
1970 ! Subroutine to calculate the adjoint of the phase function for the
1971 ! scattering model.
1972 !
1973 ! CALLING SEQUENCE:
1974 ! CALL CRTM_Phase_Matrix_AD( AtmOptics, & ! FWD Input
1975 ! Pff_AD, & ! AD Input
1976 ! Pbb_AD, & ! AD Input
1977 ! AtmOptics_AD, & ! AD Output
1978 ! RTV ) ! Internal variable
1979 !
1980 ! INPUT ARGUMENTS:
1981 ! AtmOptics: Structure containing the atmospheric optical
1982 ! parameters
1983 ! UNITS: N/A
1984 ! TYPE: CRTM_AtmOptics_type
1985 ! DIMENSION: Scalar
1986 ! ATTRIBUTES: INTENT(IN)
1987 !
1988 ! Pff_AD: Array containing the adjoint of the
1989 ! forward phase matrix.
1990 ! ** NOTE: This argument will be zeroed upon exit
1991 ! ** from this routine.
1992 ! UNITS: N/A
1993 ! TYPE: REAL
1994 ! DIMENSION: Rank-3, n_Angles x n_Angles x n_Layers
1995 ! ATTRIBUTES: INTENT(IN OUT)
1996 !
1997 ! Pbb_AD: Array containing the adjoint of the
1998 ! backward phase matrix.
1999 ! ** NOTE: This argument will be zeroed upon exit
2000 ! ** from this routine.
2001 ! UNITS: N/A
2002 ! TYPE: REAL
2003 ! DIMENSION: Rank-3, n_Angles x n_Angles x n_Layers
2004 ! ATTRIBUTES: INTENT(IN OUT)
2005 !
2006 ! RTV: Structure containing internal forward model variables
2007 ! required for subsequent tangent-linear or adjoint model
2008 ! calls. The contents of this structure are NOT accessible
2009 ! outside of the CRTM_RTSolution module.
2010 ! UNITS: N/A
2011 ! TYPE: RTV_type
2012 ! DIMENSION: Scalar
2013 ! ATTRIBUTES: INTENT(IN)
2014 !
2015 ! OUTPUT ARGUMENTS:
2016 ! AtmOptics_AD: Structure containing the adjoint atmospheric optical
2017 ! parameters
2018 ! UNITS: N/A
2019 ! TYPE: CRTM_AtmOptics_type
2020 ! DIMENSION: Scalar
2021 ! ATTRIBUTES: INTENT(IN OUT)
2022 !
2023 ! COMMENTS:
2024 ! Note the INTENT on the output AtmOptics argument is IN OUT rather than
2025 ! just OUT. This is necessary because the argument MUST be defined upon
2026 ! input. To prevent memory leaks, and in this case errors in accessing
2027 ! unallocated memory, the IN OUT INTENT is a must.
2028 !
2029 !S-
2030 !--------------------------------------------------------------------------------
2031 
2032  SUBROUTINE crtm_phase_matrix_ad( &
2033  AtmOptics, & ! FWD Input
2034  Pff_AD, & ! AD Input
2035  Pbb_AD, & ! AD Input
2036  AtmOptics_AD, & ! AD Output
2037  RTV ) ! Internal variable
2038  ! Arguments
2039  TYPE(CRTM_AtmOptics_type), INTENT(IN) :: AtmOptics
2040  REAL(fp) , INTENT(IN OUT) :: Pff_AD(:,:,:) ! n_Angles x n_Angles x n_Layers
2041  REAL(fp) , INTENT(IN OUT) :: Pbb_AD(:,:,:) ! n_Angles x n_Angles x n_Layers
2042  TYPE(CRTM_AtmOptics_type), INTENT(IN OUT) :: AtmOptics_AD
2043  TYPE(RTV_type) , INTENT(IN) :: RTV
2044  ! Local variables
2045  INTEGER :: i, j, k, l, nZ, ifac, jn
2046  REAL(fp), DIMENSION(RTV%n_Angles,RTV%n_Angles+1) :: Lff, Lbb
2047 
2048 
2049 
2050  !#--------------------------------------------------------------------------#
2051  !# -- COMPUTE THE ADJOINT OF THE PHASE MATRICES -- #
2052  !#--------------------------------------------------------------------------#
2053 
2054  nz = rtv%n_Angles
2055  IF( rtv%Solar_Flag_true ) THEN
2056  jn = rtv%n_Angles + 1
2057  ELSE
2058  jn = rtv%n_Angles
2059  END IF
2060 
2061  layer_loop: DO k = 1, rtv%n_Layers
2062 
2063 
2064  ! ------------------------------
2065  ! Only proceed if the scattering
2066  ! coefficient is significant
2067  ! ------------------------------
2068 
2069  significant_scattering: IF( atmoptics%Single_Scatter_Albedo(k) > scattering_albedo_threshold) THEN
2070 
2071 
2072  ! AD normalization to ensure energy conservation
2073  !! AtmOptics_AD%Phase_Coefficient(0:,1,k) = ZERO
2074 
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)
2077 
2078  IF ( rtv%mth_Azi == 0 ) THEN
2079  DO j = 1, jn
2080  DO i = 1, rtv%n_Angles
2081 
2082  ! For intensity, the FWD phase matrix element must >= ZERO
2083  ! so the TL, and thus the AD, for is always == zero.
2084 
2085  IF ( rtv%Off(i,j,k) < zero) lff(i,j) = phase_threshold
2086  IF ( rtv%Obb(i,j,k) < zero) lbb(i,j) = phase_threshold
2087 
2088  END DO
2089  END DO
2090 
2091  CALL normalize_phase_ad( &
2092  k, rtv, &
2093  lff, lbb, & ! FWD Input
2094  pff_ad(:,:,k), & ! AD Output
2095  pbb_ad(:,:,k) ) ! AD Output
2096  END IF
2097 
2098  DO j = 1, jn
2099  DO i = 1, rtv%n_Angles
2100 
2101  ! For intensity, the FWD phase matrix element must >= ZERO
2102  ! so the TL, and thus the AD, for is always == zero.
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
2106  END IF
2107 
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 )
2114  END DO
2115 
2116  pff_ad(i,j,k) = zero
2117  pbb_ad(i,j,k) = zero
2118 
2119  END DO
2120  END DO
2121 
2122  END IF significant_scattering
2123 
2124  END DO layer_loop
2125 
2126  END SUBROUTINE crtm_phase_matrix_ad
2127 
2128 !**************************************************************
2129 ! Normalize Phase Subroutines
2130 !**************************************************************
2131 
2132  SUBROUTINE normalize_phase( k, RTV )
2133  ! Arguments
2134  INTEGER, INTENT(IN) :: k
2135  TYPE(RTV_type), INTENT(IN OUT) :: RTV
2136  ! Local variables
2137  INTEGER :: i, j, nZ
2138 
2139  nz = rtv%n_Angles
2140 
2141  ! Normalisation for stream angles
2142  rtv%Sum_Fac(0,k)=zero
2143  DO i = 1, rtv%n_Streams
2144  rtv%n_Factor(i,k)=zero
2145  DO j = i,nz
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)
2147  END DO
2148  DO j=i,nz
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))
2151  END DO
2152  rtv%Sum_Fac(i,k)=zero
2153  IF( i < nz ) THEN
2154  DO j=i+1,nz
2155  rtv%Pff(j,i,k)=rtv%Pff(i,j,k)
2156  rtv%Pbb(j,i,k)=rtv%Pbb(i,j,k)
2157  END DO
2158  DO j = 1, i
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)
2160  END DO
2161  END IF
2162  END DO
2163 
2164  IF( rtv%n_Streams < nz ) THEN
2165  ! Sensor viewing angle differs from the Gaussian angles
2166  rtv%n_Factor(nz,k) = rtv%Sum_Fac(nz-1,k)
2167  DO j = 1, nz
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)
2170  ! Symmetric condition
2171  IF( j < nz ) THEN
2172  rtv%Pff(nz,j,k) = rtv%Pff(j,nz,k)
2173  rtv%Pbb(nz,j,k) = rtv%Pbb(j,nz,k)
2174  END IF
2175  END DO
2176  END IF
2177 
2178  END SUBROUTINE normalize_phase
2179 
2180  SUBROUTINE normalize_phase_tl( k, RTV, Pff, Pbb, Pff_TL, Pbb_TL )
2181  ! Arguments
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(:,:)
2188  ! Local variables
2189  REAL(fp) :: n_Factor_TL
2190  REAL(fp) :: Sum_Fac_TL(0:RTV%n_Angles)
2191  INTEGER :: i, j, nZ
2192 
2193  nz = rtv%n_Angles
2194 
2195  ! Normalisation for stream angles
2196  sum_fac_tl(0) = zero
2197  DO i = 1, rtv%n_Streams
2198  n_factor_tl = zero
2199  DO j = i,nz
2200  n_factor_tl=n_factor_tl+(pff_tl(i,j)+pbb_tl(i,j))*rtv%COS_Weight(j)
2201  END DO
2202  DO j=i,nz
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)
2206 
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)
2210  END DO
2211  sum_fac_tl(i)=zero
2212  ! Symmetric condition
2213  IF( i < nz ) THEN
2214  DO j = i+1, nz
2215  pff_tl(j,i) = pff_tl(i,j)
2216  pbb_tl(j,i) = pbb_tl(i,j)
2217  END DO
2218  DO j = 1, i
2219  sum_fac_tl(i)=sum_fac_tl(i) + (pff_tl(j,i+1)+pbb_tl(j,i+1))*rtv%COS_Weight(j)
2220  END DO
2221  END IF
2222  END DO
2223 
2224  IF( rtv%n_Streams < nz ) THEN
2225  ! Sensor viewing angle differs from the Gaussian angles
2226  n_factor_tl = sum_fac_tl(nz-1)
2227  DO j = 1, nz
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
2230 
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
2233  ! Symmetric condition
2234  IF( j < nz ) THEN
2235  pff_tl(nz,j) = pff_tl(j,nz)
2236  pbb_tl(nz,j) = pbb_tl(j,nz)
2237  END IF
2238  END DO
2239  END IF
2240 
2241  END SUBROUTINE normalize_phase_tl
2242 
2243  SUBROUTINE normalize_phase_ad( k, RTV, Pff, Pbb, Pff_AD, Pbb_AD )
2244  ! Arguments
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(:,:)
2251  ! Local variables
2252  INTEGER :: i, j, nZ
2253  REAL(fp) :: n_Factor_AD
2254  REAL(fp) :: Sum_Fac_AD(0:RTV%n_Angles)
2255 
2256 
2257  nz = rtv%n_Angles
2258  sum_fac_ad = zero
2259 
2260  n_factor_ad = zero
2261  IF( rtv%n_Streams < nz ) THEN
2262  ! Sensor viewing angle diffs from the Gaussian angles
2263  DO j = nz, 1, -1
2264  ! Symmetric condition
2265  IF( j < nz ) THEN
2266  pff_ad(j,nz) = pff_ad(j,nz) + pff_ad(nz,j)
2267  pff_ad(nz,j) = zero
2268  pbb_ad(j,nz) = pbb_ad(j,nz) + pbb_ad(nz,j)
2269  pbb_ad(nz,j) = zero
2270  END IF
2271 
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)
2274 
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)
2277  END DO
2278  sum_fac_ad(nz-1) = n_factor_ad
2279  n_factor_ad = zero
2280  END IF
2281 
2282  DO i = rtv%n_Streams, 1, -1
2283  ! Symmetric condition
2284  IF( i < nz ) THEN
2285  DO j = i, 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)
2288  END DO
2289  DO j = nz,i+1,-1
2290  pff_ad(i,j) = pff_ad(i,j) + pff_ad(j,i)
2291  pff_ad(j,i) = zero
2292  pbb_ad(i,j) = pbb_ad(i,j) + pbb_ad(j,i)
2293  pbb_ad(j,i) = zero
2294  END DO
2295  END IF
2296  sum_fac_ad(i) = zero
2297  DO j = nz, i, -1
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))
2302 
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))
2307  END DO
2308  DO j = nz, i, -1
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)
2311  END DO
2312  n_factor_ad = zero
2313  END DO
2314  sum_fac_ad(0) = zero
2315  END SUBROUTINE normalize_phase_ad
2316 
2317 END MODULE common_rtsolution
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
Definition: Type_Kinds.f90:124
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
Definition: RTV_Define.f90:72
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
Definition: RTV_Define.f90:76
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