FV3 Bundle
CRTM_AerosolScatter.f90
Go to the documentation of this file.
1 !
2 ! CRTM_AerosolScatter
3 !
4 !
5 ! Module to compute the aerosol absorption and scattering properties
6 ! required for radiative transfer in an atmosphere with aerosols.
7 !
8 !
9 ! CREATION HISTORY:
10 ! Written by: Paul van Delst, 15-Feb-2005
11 ! paul.vandelst@noaa.gov
12 ! Modified by Quanhua Liu, 03-Oct-2006
13 ! Quanhua.Liu@noaa.gov
14 !
15 
17 
18  ! -----------------
19  ! Environment setup
20  ! -----------------
21  ! Module use
22  USE type_kinds, ONLY: fp
25  max_n_layers, &
28  bs_threshold, &
31  hgphase ! <<< NEED TO REMOVE THIS IN FUTURE
32  USE crtm_spccoeff, ONLY: sc, &
33  spccoeff_ismicrowavesensor , &
34  spccoeff_isinfraredsensor , &
35  spccoeff_isvisiblesensor , &
36  spccoeff_isultravioletsensor
37  USE crtm_aerosolcoeff, ONLY: aeroc
39  dust_aerosol , &
40  seasalt_ssam_aerosol , &
41  seasalt_sscm1_aerosol , &
42  seasalt_sscm2_aerosol , &
43  seasalt_sscm3_aerosol , &
44  organic_carbon_aerosol , &
45  black_carbon_aerosol , &
46  sulfate_aerosol
48  USE crtm_interpolation, ONLY: npts , &
49  lpoly_type , &
50  find_index , &
51  interp_1d , &
52  interp_2d , &
53  interp_3d , &
54  interp_2d_tl, &
55  interp_3d_tl, &
56  interp_2d_ad, &
57  interp_3d_ad, &
58  clear_lpoly , &
59  lpoly , &
60  lpoly_tl , &
61  lpoly_ad
63 
64  ! Internal variable definition module
65  USE asvar_define, ONLY: asvar_type, &
66  asinterp_type, &
68  asvar_destroy , &
70 
71 
72  ! Disable implicit typing
73  IMPLICIT NONE
74 
75 
76  ! ------------
77  ! Visibilities
78  ! ------------
79  ! Everything private by default
80  PRIVATE
81  ! Procedures
85 
86 
87  ! -----------------
88  ! Module parameters
89  ! -----------------
90  ! Version Id for the module
91  CHARACTER(*), PARAMETER :: module_version_id = &
92  '$Id: CRTM_AerosolScatter.f90 60152 2015-08-13 19:19:13Z paul.vandelst@noaa.gov $'
93  ! Message string length
94  INTEGER, PARAMETER :: ml = 256
95  ! Number of stream angle definitions
96  INTEGER, PARAMETER :: two_streams = 2
97  INTEGER, PARAMETER :: four_streams = 4
98  INTEGER, PARAMETER :: six_streams = 6
99  INTEGER, PARAMETER :: eight_streams = 8
100  INTEGER, PARAMETER :: sixteen_streams = 16
101  INTEGER, PARAMETER :: thirtytwo_streams = 32
102 
103 
104 CONTAINS
105 
106 
107 !------------------------------------------------------------------------------
108 !:sdoc+:
109 !
110 ! NAME:
111 ! CRTM_Compute_AerosolScatter
112 !
113 ! PURPOSE:
114 ! Function to compute the aerosol absorption and scattering properties
115 ! and populate the output AerosolScatter structure for a single channel.
116 !
117 ! CALLING SEQUENCE:
118 ! Error_Status = CRTM_Compute_AerosolScatter( Atmosphere , &
119 ! SensorIndex , &
120 ! ChannelIndex , &
121 ! AerosolScatter, &
122 ! ASvar )
123 !
124 ! INPUT ARGUMENTS:
125 ! Atmosphere: CRTM_Atmosphere structure containing the atmospheric
126 ! profile data.
127 ! UNITS: N/A
128 ! TYPE: CRTM_Atmosphere_type
129 ! DIMENSION: Scalar
130 ! ATTRIBUTES: INTENT(IN)
131 !
132 ! ChannelIndex: Channel index id. This is a unique index associated
133 ! with a (supported) sensor channel used to access the
134 ! shared coefficient data.
135 ! UNITS: N/A
136 ! TYPE: INTEGER
137 ! DIMENSION: Scalar
138 ! ATTRIBUTES: INTENT(IN)
139 !
140 ! SensorIndex: Sensor index id. This ia a unique index associated
141 ! with a sensor.
142 ! UNITS: N/A
143 ! TYPE: INTEGER
144 ! DIMENSION: Scalar
145 ! ATTRIBUTES: INTENT(IN)
146 !
147 ! OUTPUT ARGUMENTS:
148 ! AerosolScatter: CRTM_AtmOptics structure containing the aerosol
149 ! absorption and scattering properties required by
150 ! the radiative transfer.
151 ! UNITS: N/A
152 ! TYPE: CRTM_AtmOptics_type
153 ! DIMENSION: Scalar
154 ! ATTRIBUTES: INTENT(IN OUT)
155 !
156 ! ASvar: Structure containing internal variables required for
157 ! subsequent tangent-linear or adjoint model calls.
158 ! UNITS: N/A
159 ! TYPE: ASvar_type
160 ! DIMENSION: Scalar
161 ! ATTRIBUTES: INTENT(OUT)
162 !
163 ! FUNCTION RESULT:
164 ! Error_Status: The return value is an integer defining the error status.
165 ! The error codes are defined in the ERROR_HANDLER module.
166 ! If == SUCCESS the computation was sucessful
167 ! == FAILURE an unrecoverable error occurred
168 ! UNITS: N/A
169 ! TYPE: INTEGER
170 ! DIMENSION: Scalar
171 !
172 !:sdoc-:
173 !------------------------------------------------------------------------------
174 
175  FUNCTION crtm_compute_aerosolscatter( &
176  Atm , & ! Input
177  SensorIndex , & ! Input
178  ChannelIndex, & ! Input
179  AScat , & ! Output
180  ASV ) & ! Internal variable output
181  result( error_status )
182  !Arguments
183  TYPE(crtm_atmosphere_type), INTENT(IN) :: atm
184  INTEGER , INTENT(IN) :: channelindex
185  INTEGER , INTENT(IN) :: sensorindex
186  TYPE(crtm_atmoptics_type) , INTENT(IN OUT) :: ascat
187  TYPE(asvar_type) , INTENT(IN OUT) :: asv
188  ! Function result
189  INTEGER :: error_status
190  ! Local parameters
191  CHARACTER(*), PARAMETER :: routine_name = 'CRTM_Compute_AerosolScatter'
192  ! Local Variables
193  CHARACTER(ML) :: message
194  INTEGER :: k, ka, l, m, n
195  REAL(fp) :: frequency
196  LOGICAL :: layer_mask(atm%n_layers)
197  INTEGER :: layer_index(atm%n_layers)
198  INTEGER :: naerosol_layers
199  REAL(fp) :: bs
200 
201  ! ------
202  ! Set up
203  ! ------
204  error_status = success
205  IF ( spccoeff_ismicrowavesensor(sc(sensorindex)) ) RETURN
206  IF ( atm%n_Aerosols == 0 ) RETURN
207  asv%Total_bs = zero
208  ! Frequency in inverse centimetres
209  frequency = sc(sensorindex)%Wavenumber(channelindex)
210  ! Determine offset for Legendre coefficients in the
211  ! AeroC lookup table corresponding to the
212  ! number of streams
213  SELECT CASE(ascat%n_Legendre_Terms)
214  CASE (two_streams) ; ascat%lOffset = 0
215  CASE (four_streams) ; ascat%lOffset = 0
216  CASE (six_streams) ; ascat%lOffset = 5
217  CASE (eight_streams) ; ascat%lOffset = 12
218  CASE (sixteen_streams); ascat%lOffset = 21
219  CASE DEFAULT
220  ascat%lOffset = 0
221  ! Use two-stream model or HG and RAYLEIGH Phase function
222  IF( hgphase ) THEN
223  ascat%n_Legendre_Terms = 0
224  ELSE
225  error_status = failure
226  WRITE(message,'("The n_Legendre_Terms in AerosolScatter, ",i0,", do not fit model")') &
227  ascat%n_Legendre_Terms
228  CALL display_message(routine_name, message, error_status )
229  RETURN
230  END IF
231  END SELECT
232 
233 
234  ! -----------------------------------------------
235  ! Loop over the different Aerosols in the profile
236  ! -----------------------------------------------
237  aerosol_loop: DO n = 1, atm%n_Aerosols
238 
239  ! Only process aerosols with more
240  ! than the threshold Aerosol amount
241  layer_mask = atm%Aerosol(n)%Concentration > aerosol_content_threshold
242  naerosol_layers = count(layer_mask)
243  IF ( naerosol_layers == 0 ) cycle aerosol_loop
244 
245 
246  ! --------------------------------------
247  ! Loop over the current Aerosol's Layers
248  ! --------------------------------------
249  layer_index(1:naerosol_layers) = pack((/(k,k=1,atm%Aerosol(n)%n_Layers)/), layer_mask)
250  aerosol_layer_loop: DO k = 1, naerosol_layers
251  ka = layer_index(k)
252 
253  ! Get the aerosol optical properties
254  CALL get_aerosol_opt( ascat , & ! Input
255  frequency , & ! Input
256  atm%Aerosol(n)%Type , & ! Input
257  atm%Aerosol(n)%Effective_Radius(ka), & ! Input
258  asv%ke(ka,n) , & ! Output
259  asv%w(ka,n) , & ! Output
260  asv%pcoeff(:,:,ka,n) , & ! Output
261  asv%asi(ka,n) ) ! Interpolation
262 
263  ! interpolation quality control
264  IF( asv%ke(ka,n) <= zero ) THEN
265  asv%ke(ka,n) = zero
266  asv%w(ka,n) = zero
267  END IF
268  IF( asv%w(ka,n) <= zero ) THEN
269  asv%w(ka,n) = zero
270  asv%pcoeff(:,:,ka,n) = zero
271  END IF
272  IF( asv%w(ka,n) >= one ) THEN
273  asv%w(ka,n) = one
274  END IF
275 
276  ! Compute the optical depth (absorption + scattering)
277  ! tau = rho.ke
278  ! where
279  ! rho = Integrated Aerosol Concentration for a layer(kg/m^2) [M.L^-2]
280  ! ke = mass extintion coefficient (m^2/kg) [L^2.M^-1]
281  ! Note that since all these computations are done for a given
282  ! layer, the optical depth is the same as the volume extinction
283  ! coefficient, be. Usually,
284  ! tau = be.d(z)
285  ! but we are working with height/thickness independent quantities
286  ! so that
287  ! tau = be
288  ! This is why the optical depth is used in the denominator to
289  ! compute the single scatter albedo in the Layer_loop below.
290  ascat%Optical_Depth(ka) = ascat%Optical_Depth(ka) + &
291  (asv%ke(ka,n)*atm%Aerosol(n)%Concentration(ka))
292 
293  ! Compute the phase matrix coefficients
294  ! p = p + p(LUT)*bs
295  ! where
296  ! p(LUT) = the phase coefficient from the LUT
297  IF( ascat%n_Phase_Elements > 0 .and. ascat%Include_Scattering ) THEN
298  ! Compute the volume scattering coefficient for the current
299  ! aerosol layer and accumulate it for the layer total for the
300  ! profile (i.e. all aerosols)
301  ! bs = rho.w.ke
302  ! where
303  ! bs = volume scattering coefficient for a layer [dimensionless]
304  ! rho = integrated aerosol concentration for a layer (kg/m^2) [M.L^-2]
305  ! w = single scatter albedo [dimensionless]
306  ! ke = mass extintion coefficient (m^2/kg) [L^2.M^-1]
307  ! qliu
308  bs = atm%Aerosol(n)%Concentration(ka) * asv%ke(ka,n) * asv%w(ka,n)
309  asv%Total_bs(ka) = asv%Total_bs(ka) + bs
310  ascat%Single_Scatter_Albedo(ka) = ascat%Single_Scatter_Albedo(ka) + bs
311 
312  DO m = 1, ascat%n_Phase_Elements
313  DO l = 0, ascat%n_Legendre_Terms
314  ascat%Phase_Coefficient(l,m,ka) = ascat%Phase_Coefficient(l,m,ka) + &
315  (asv%pcoeff(l,m,ka,n) * bs)
316  END DO
317  END DO
318  END IF
319  END DO aerosol_layer_loop
320  END DO aerosol_loop
321 
322 
323  END FUNCTION crtm_compute_aerosolscatter
324 
325 
326 !------------------------------------------------------------------------------
327 !:sdoc+:
328 !
329 ! NAME:
330 ! CRTM_Compute_AerosolScatter_TL
331 !
332 ! PURPOSE:
333 ! Function to compute the tangent-linear aerosol absorption and
334 ! scattering properties and populate the output AerosolScatter_TL
335 ! structure for a single channel.
336 !
337 ! CALLING SEQUENCE:
338 ! Error_Status = CRTM_Compute_AerosolScatter_TL( Atmosphere , &
339 ! AerosolScatter , &
340 ! Atmosphere_TL , &
341 ! SensorIndex , &
342 ! ChannelIndex , &
343 ! AerosolScatter_TL, &
344 ! ASvar )
345 !
346 ! INPUT ARGUMENTS:
347 ! Atmosphere: CRTM_Atmosphere structure containing the atmospheric
348 ! profile data.
349 ! UNITS: N/A
350 ! TYPE: CRTM_Atmosphere_type
351 ! DIMENSION: Scalar
352 ! ATTRIBUTES: INTENT(IN)
353 !
354 ! AerosolScatter: CRTM_AtmOptics structure containing the forward model
355 ! aerosol absorption and scattering properties required
356 ! for radiative transfer.
357 ! UNITS: N/A
358 ! TYPE: CRTM_AtmOptics_type
359 ! DIMENSION: Scalar
360 ! ATTRIBUTES: INTENT(IN)
361 !
362 ! Atmosphere_TL: CRTM Atmosphere structure containing the tangent-linear
363 ! atmospheric state data.
364 ! UNITS: N/A
365 ! TYPE: CRTM_Atmosphere_type
366 ! DIMENSION: Scalar
367 ! ATTRIBUTES: INTENT(IN)
368 !
369 ! SensorIndex: Sensor index id. This ia a unique index associated
370 ! with a sensor.
371 ! UNITS: N/A
372 ! TYPE: INTEGER
373 ! DIMENSION: Scalar
374 ! ATTRIBUTES: INTENT(IN)
375 !
376 !
377 ! Channel_Index: Channel index id. This is a unique index associated
378 ! with a (supported) sensor channel used to access the
379 ! shared coefficient data.
380 ! UNITS: N/A
381 ! TYPE: INTEGER
382 ! DIMENSION: Scalar
383 ! ATTRIBUTES: INTENT(IN)
384 !
385 ! ASvar: Structure containing internal variables required for
386 ! subsequent tangent-linear or adjoint model calls.
387 ! UNITS: N/A
388 ! TYPE: ASvar_type
389 ! DIMENSION: Scalar
390 ! ATTRIBUTES: INTENT(IN)
391 !
392 ! OUTPUT ARGUMENTS:
393 ! AerosolScatter_TL: CRTM_AtmOptics structure containing the tangent-linear
394 ! aerosol absorption and scattering properties required
395 ! for radiative transfer.
396 ! UNITS: N/A
397 ! TYPE: CRTM_AtmOptics_type
398 ! DIMENSION: Scalar
399 ! ATTRIBUTES: INTENT(IN OUT)
400 !
401 ! FUNCTION RESULT:
402 ! Error_Status: The return value is an integer defining the error status.
403 ! The error codes are defined in the ERROR_HANDLER module.
404 ! If == SUCCESS the computation was sucessful
405 ! == FAILURE an unrecoverable error occurred
406 ! UNITS: N/A
407 ! TYPE: INTEGER
408 ! DIMENSION: Scalar
409 !
410 !:sdoc-:
411 !------------------------------------------------------------------------------
412 
414  Atm , & ! FWD Input
415  AScat , & ! FWD Input
416  Atm_TL , & ! TL Input
417  SensorIndex , & ! Input
418  ChannelIndex, & ! Input
419  AScat_TL , & ! TL Input
420  ASV ) & ! Internal variable input
421  result( error_status )
422  ! Arguments
423  TYPE(crtm_atmosphere_type) , INTENT(IN) :: atm
424  TYPE(crtm_atmoptics_type) , INTENT(IN) :: ascat
425  TYPE(crtm_atmosphere_type) , INTENT(IN) :: atm_tl
426  INTEGER , INTENT(IN) :: sensorindex
427  INTEGER , INTENT(IN) :: channelindex
428  TYPE(crtm_atmoptics_type) , INTENT(IN OUT) :: ascat_tl
429  TYPE(asvar_type) , INTENT(IN) :: asv
430  ! Function result
431  INTEGER :: error_status
432  ! Local parameters
433  CHARACTER(*), PARAMETER :: routine_name = 'CRTM_Compute_AerosolScatter_TL'
434  ! Local variables
435  INTEGER :: k, ka, l, m, n
436  INTEGER :: n_legendre_terms, n_phase_elements
437  REAL(fp) :: frequency
438  LOGICAL :: layer_mask(atm%n_layers)
439  INTEGER :: layer_index(atm%n_layers)
440  INTEGER :: naerosol_layers
441  REAL(fp) :: ke_tl, w_tl
442  REAL(fp) :: pcoeff_tl(0:ascat%n_legendre_terms, ascat%n_phase_elements)
443  REAL(fp) :: bs, bs_tl
444 
445  ! ------
446  ! Set up
447  ! ------
448  error_status = success
449  IF ( spccoeff_ismicrowavesensor(sc(sensorindex)) ) RETURN
450  IF ( atm%n_Aerosols == 0 ) RETURN
451  ! Frequency
452  frequency = sc(sensorindex)%Wavenumber(channelindex)
453  ! Phase matrix dimensions
454  n_legendre_terms = ascat_tl%n_Legendre_Terms
455  n_phase_elements = ascat_tl%n_Phase_Elements
456  ascat_tl%lOffset = ascat%lOffset
457 
458 
459  ! -----------------------------------------------
460  ! Loop over the different Aerosols in the profile
461  ! -----------------------------------------------
462  aerosol_loop: DO n = 1, atm%n_Aerosols
463 
464  ! Only process aerosols with more
465  ! than the threshold aerosol amount
466  layer_mask = atm%Aerosol(n)%Concentration > aerosol_content_threshold
467  naerosol_layers = count(layer_mask)
468  IF (naerosol_layers == 0) cycle aerosol_loop
469 
470 
471  ! --------------------------------------
472  ! Loop over the current aerosol's layers
473  ! --------------------------------------
474  layer_index(1:naerosol_layers) = pack((/(k,k=1,atm%Aerosol(n)%n_Layers)/), layer_mask)
475  aerosol_layer_loop: DO k = 1, naerosol_layers
476  ka = layer_index(k)
477 
478  ! Obtain bulk aerosol optical properties
479  CALL get_aerosol_opt_tl(ascat_tl , & ! Input
480  atm%Aerosol(n)%Type , & ! Input
481  asv%ke(ka,n) , & ! Input
482  asv%w(ka,n) , & ! Input
483  atm_tl%Aerosol(n)%Effective_Radius(ka), & ! TL Input
484  ke_tl , & ! TL Output
485  w_tl , & ! TL Output
486  pcoeff_tl , & ! TL Output
487  asv%asi(ka,n) ) ! Interpolation
488 
489  ! interpolation quality control
490  IF( asv%ke(ka,n) <= zero ) THEN
491  ke_tl = zero
492  w_tl = zero
493  END IF
494  IF( asv%w(ka,n) <= zero ) THEN
495  w_tl = zero
496  pcoeff_tl = zero
497  END IF
498  IF( asv%w(ka,n) >= one ) THEN
499  w_tl = zero
500  END IF
501 
502  ! Compute the optical depth (absorption + scattering)
503  ascat_tl%Optical_Depth(ka) = ascat_tl%Optical_Depth(ka) + &
504  (ke_tl * atm%Aerosol(n)%Concentration(ka)) + &
505  (asv%ke(ka,n) * atm_tl%Aerosol(n)%Concentration(ka))
506 
507  ! Compute the phase matrix coefficients
508  IF( n_phase_elements > 0 .and. ascat%Include_Scattering ) THEN
509  ! Compute the volume scattering coefficient
510  bs = atm%Aerosol(n)%Concentration(ka) * asv%ke(ka,n) * asv%w(ka,n)
511  bs_tl = (atm_tl%Aerosol(n)%Concentration(ka) * asv%ke(ka,n) * asv%w(ka,n) ) + &
512  (atm%Aerosol(n)%Concentration(ka) * ke_tl * asv%w(ka,n) ) + &
513  (atm%Aerosol(n)%Concentration(ka) * asv%ke(ka,n) * w_tl )
514  ascat_tl%Single_Scatter_Albedo(ka) = ascat_tl%Single_Scatter_Albedo(ka) + bs_tl
515 
516  DO m = 1, n_phase_elements
517  DO l = 0, n_legendre_terms
518  ascat_tl%Phase_Coefficient(l,m,ka) = ascat_tl%Phase_Coefficient(l,m,ka) + &
519  (pcoeff_tl(l,m) * bs ) + &
520  (asv%pcoeff(l,m,ka,n) * bs_tl)
521  END DO
522  END DO
523  END IF
524  END DO aerosol_layer_loop
525  END DO aerosol_loop
526 
527  END FUNCTION crtm_compute_aerosolscatter_tl
528 
529 
530 !------------------------------------------------------------------------------
531 !:sdoc+:
532 !
533 ! NAME:
534 ! CRTM_Compute_AerosolScatter_AD
535 !
536 ! PURPOSE:
537 ! Function to compute the adjoint aerosol absorption and scattering
538 ! properties for a single channel.
539 !
540 ! CALLING SEQUENCE:
541 ! Error_Status = CRTM_Compute_AerosolScatter_AD( Atmosphere , &
542 ! AerosolScatter , &
543 ! AerosolScatter_AD, &
544 ! SensorIndex , &
545 ! ChannelIndex , &
546 ! Atmosphere_AD , &
547 ! ASVariables )
548 !
549 ! INPUT ARGUMENTS:
550 ! Atmosphere: CRTM_Atmosphere structure containing the atmospheric
551 ! profile data.
552 ! UNITS: N/A
553 ! TYPE: CRTM_Atmosphere_type
554 ! DIMENSION: Scalar
555 ! ATTRIBUTES: INTENT(IN)
556 !
557 ! AerosolScatter: CRTM_AtmOptics structure containing the forward model
558 ! aerosol absorption and scattering properties required
559 ! for radiative transfer.
560 ! UNITS: N/A
561 ! TYPE: CRTM_AtmOptics_type
562 ! DIMENSION: Scalar
563 ! ATTRIBUTES: INTENT(IN)
564 !
565 ! AerosolScatter_AD: CRTM_AtmOptics structure containing the adjoint
566 ! aerosol absorption and scattering properties.
567 ! **NOTE: On EXIT from this function, the contents of
568 ! this structure may be modified (e.g. set to
569 ! zero.)
570 ! UNITS: N/A
571 ! TYPE: CRTM_AtmOptics_type
572 ! DIMENSION: Scalar
573 ! ATTRIBUTES: INTENT(IN OUT)
574 !
575 ! SensorIndex: Sensor index id. This ia a unique index associated
576 ! with a sensor.
577 ! UNITS: N/A
578 ! TYPE: INTEGER
579 ! DIMENSION: Scalar
580 ! ATTRIBUTES: INTENT(IN)
581 !
582 ! Channel_Index: Channel index id. This is a unique index associated
583 ! with a (supported) sensor channel used to access the
584 ! shared coefficient data.
585 ! UNITS: N/A
586 ! TYPE: INTEGER
587 ! DIMENSION: Scalar
588 ! ATTRIBUTES: INTENT(IN)
589 !
590 ! ASvar: Structure containing internal variables required for
591 ! subsequent tangent-linear or adjoint model calls.
592 ! UNITS: N/A
593 ! TYPE: ASvar_type
594 ! DIMENSION: Scalar
595 ! ATTRIBUTES: INTENT(IN)
596 !
597 ! OUTPUT ARGUMENTS:
598 ! Atmosphere_AD: CRTM Atmosphere structure containing the adjoint
599 ! atmospheric state data.
600 ! UNITS: N/A
601 ! TYPE: CRTM_Atmosphere_type
602 ! DIMENSION: Scalar
603 ! ATTRIBUTES: INTENT(IN OUT)
604 !
605 !
606 ! FUNCTION RESULT:
607 ! Error_Status: The return value is an integer defining the error
608 ! status. The error codes are defined in the
609 ! ERROR_HANDLER module.
610 ! If == SUCCESS the computation was sucessful
611 ! == FAILURE an unrecoverable error occurred
612 ! UNITS: N/A
613 ! TYPE: INTEGER
614 ! DIMENSION: Scalar
615 !
616 !:sdoc-:
617 !------------------------------------------------------------------------------
618 
620  Atm , & ! FWD Input
621  AScat , & ! FWD Input
622  AScat_AD , & ! AD Input
623  SensorIndex , & ! Input
624  ChannelIndex, & ! Input
625  Atm_AD , & ! AD Output
626  ASV ) & ! Internal Variable input
627  result( error_status )
628  ! Arguments
629  TYPE(crtm_atmosphere_type), INTENT(IN) :: atm
630  TYPE(crtm_atmoptics_type) , INTENT(IN) :: ascat
631  TYPE(crtm_atmoptics_type) , INTENT(IN OUT) :: ascat_ad
632  INTEGER , INTENT(IN) :: sensorindex
633  INTEGER , INTENT(IN) :: channelindex
634  TYPE(crtm_atmosphere_type), INTENT(IN OUT) :: atm_ad
635  TYPE(asvar_type) , INTENT(IN) :: asv
636  ! Function result
637  INTEGER :: error_status
638  ! Local parameters
639  CHARACTER(*), PARAMETER :: routine_name = 'CRTM_Compute_AerosolScatter_AD'
640  ! Local variables
641  INTEGER :: k, ka, l, m, n
642  INTEGER :: n_legendre_terms, n_phase_elements
643  REAL(fp) :: frequency
644  LOGICAL :: layer_mask(atm%n_layers)
645  INTEGER :: layer_index(atm%n_layers)
646  INTEGER :: naerosol_layers
647  REAL(fp) :: ke_ad, w_ad
648  REAL(fp) :: pcoeff_ad(0:ascat%n_legendre_terms, ascat%n_phase_elements)
649  REAL(fp) :: bs, bs_ad
650 
651  ! ------
652  ! Set up
653  ! ------
654  error_status = success
655  IF ( spccoeff_ismicrowavesensor(sc(sensorindex)) ) RETURN
656  IF ( atm%n_Aerosols == 0 ) RETURN
657  ! Frequency
658  frequency = sc(sensorindex)%Wavenumber(channelindex)
659  ! Phase matrix dimensions
660  n_legendre_terms = ascat_ad%n_Legendre_Terms
661  n_phase_elements = ascat_ad%n_Phase_Elements
662  ascat_ad%lOffset = ascat%lOffset
663 
664 
665  ! ----------------------------------------------------------
666  ! Adjoint of accumulated optical properties for all aerosols
667  ! ----------------------------------------------------------
668  ! Shorten variable name
669  l = n_legendre_terms
670 
671  ! ------------------------------------------------
672  ! Loop over different types of aerosols in profile
673  ! ------------------------------------------------
674  aerosol_loop: DO n = 1, atm%n_Aerosols
675 
676  ! Only process aerosols with more than
677  ! the threshold aerosol concentration
678  layer_mask = atm%Aerosol(n)%Concentration > aerosol_content_threshold
679  naerosol_layers = count(layer_mask)
680  IF ( naerosol_layers == 0 ) cycle aerosol_loop
681 
682  ! --------------------------------------
683  ! Loop over the current aerosol's layers
684  ! --------------------------------------
685  layer_index(1:naerosol_layers) = pack((/(k,k=1,atm%Aerosol(n)%n_Layers)/), layer_mask)
686  aerosol_layer_loop: DO k = 1, naerosol_layers
687  ka = layer_index(k)
688 
689  ! Initialize the individual
690  ! Aerosol adjoint variables
691  bs_ad = zero
692  pcoeff_ad = zero
693  ke_ad = zero
694  w_ad = zero
695 
696  ! Compute the adjoint of the
697  ! phase matrix coefficients
698  IF( n_phase_elements > 0 .and. ascat%Include_Scattering ) THEN
699  ! Recompute the forward model volume scattering
700  ! coefficient for the current aerosol type ONLY
701  bs = atm%Aerosol(n)%Concentration(ka) * asv%ke(ka,n) * asv%w(ka,n)
702  DO m = 1, n_phase_elements
703  DO l = 0, n_legendre_terms
704  bs_ad = bs_ad + (asv%pcoeff(l,m,ka,n) * ascat_ad%Phase_Coefficient(l,m,ka))
705  pcoeff_ad(l,m) = pcoeff_ad(l,m) + (bs * ascat_ad%Phase_Coefficient(l,m,ka))
706  END DO
707  END DO
708  ! NOTE: bs_AD is not reinitialized after this
709  ! point since it is reinitialized at the
710  ! start of the Aerosol_Layer_loop
711  bs_ad = bs_ad + ascat_ad%Single_Scatter_Albedo(ka)
712  w_ad = w_ad + (atm%Aerosol(n)%Concentration(ka) * asv%ke(ka,n)* bs_ad )
713  END IF
714 
715  ! Compute the adjoint of the optical
716  ! depth (absorption + scattering)
717  atm_ad%Aerosol(n)%Concentration(ka) = atm_ad%Aerosol(n)%Concentration(ka) + &
718  (asv%ke(ka,n) * ascat_ad%Optical_Depth(ka))
719  ke_ad = ke_ad + (atm%Aerosol(n)%Concentration(ka) * ascat_ad%Optical_Depth(ka))
720 
721  ! Compute the adjoint of the volume
722  ! scattering coefficient.
723 
724  ke_ad = ke_ad + (atm%Aerosol(n)%Concentration(ka) * bs_ad * asv%w(ka,n) )
725  atm_ad%Aerosol(n)%Concentration(ka) = atm_ad%Aerosol(n)%Concentration(ka) + &
726  ( bs_ad * asv%ke(ka,n) * asv%w(ka,n) )
727 
728  ! interpolation quality control
729  IF( asv%w(ka,n) >= one ) THEN
730  w_ad = zero
731  END IF
732  IF( asv%ke(ka,n) <= zero ) THEN
733  ke_ad = zero
734  w_ad = zero
735  END IF
736  IF( asv%w(ka,n) <= zero ) THEN
737  w_ad = zero
738  pcoeff_ad = zero
739  END IF
740  ! interpolation quality control
741  IF( asv%w(ka,n) >= one ) THEN
742  w_ad = zero
743  END IF
744  IF( asv%ke(ka,n) <= zero ) THEN
745  ke_ad = zero
746  w_ad = zero
747  END IF
748  IF( asv%w(ka,n) <= zero ) THEN
749  w_ad = zero
750  pcoeff_ad = zero
751  END IF
752 
753  ! Adjoint AScat interpolation routine
754  CALL get_aerosol_opt_ad(ascat_ad , & ! Input
755  atm%Aerosol(n)%Type , & ! Input
756  asv%ke(ka,n) , & ! Output
757  asv%w(ka,n) , & ! Output
758  ke_ad , & ! AD Input
759  w_ad , & ! AD Input
760  pcoeff_ad , & ! AD Input
761  atm_ad%Aerosol(n)%Effective_Radius(ka), & ! AD Output
762  asv%asi(ka,n) ) ! Interpolation
763 
764  END DO aerosol_layer_loop
765  END DO aerosol_loop
766 
767  END FUNCTION crtm_compute_aerosolscatter_ad
768 
769 
770 
771 !################################################################################
772 !################################################################################
773 !## ##
774 !## ## PRIVATE MODULE ROUTINES ## ##
775 !## ##
776 !################################################################################
777 !################################################################################
778 
779  ! --------------------------------------------
780  ! Determine the aerosol type index, k, for the
781  ! aerosols based on AeroC LUT organisation
782  ! --------------------------------------------
783  FUNCTION aerosol_type_index( Aerosol_Type ) RESULT( k )
784  INTEGER, INTENT(IN) :: aerosol_type
785  INTEGER :: k
786  SELECT CASE (aerosol_type)
787  CASE(dust_aerosol) ; k=1
788  CASE(seasalt_ssam_aerosol) ; k=2
789  CASE(seasalt_sscm1_aerosol) ; k=3
790  CASE(seasalt_sscm2_aerosol) ; k=4
791  CASE(seasalt_sscm3_aerosol) ; k=5
792  CASE(organic_carbon_aerosol) ; k=6
793  CASE(black_carbon_aerosol) ; k=7
794  CASE(sulfate_aerosol) ; k=8
795  END SELECT
796  END FUNCTION aerosol_type_index
797 
798 
799  ! ----------------------------------------
800  ! Subroutine to obtain the bulk
801  ! optical properties of aerosol
802  ! extinction coefficient (ke),
803  ! scattering coefficient (w)
804  ! asymmetry factor (g), and
805  ! spherical Legendre coefficients (pcoeff)
806  ! ----------------------------------------
807  SUBROUTINE get_aerosol_opt( AerosolScatter, & ! Input AerosolScatter structure
808  Frequency , & ! Input in cm^-1
809  Aerosol_Type , & ! Input see CRTM_Aerosol_Define.f90
810  Reff , & ! Input effective radius (mm)
811  ke , & ! Output extinction coefficient (=~ optical depth)
812  w , & ! Output single scattering albedo
813  pcoeff , & ! Output spherical Legendre coefficients
814  asi ) ! Output interpolation data
815  ! Arguments
816  TYPE(CRTM_AtmOptics_type), INTENT(IN) :: AerosolScatter
817  REAL(fp) , INTENT(IN) :: Frequency
818  INTEGER , INTENT(IN) :: Aerosol_Type
819  REAL(fp) , INTENT(IN) :: Reff
820  REAL(fp) , INTENT(OUT) :: ke
821  REAL(fp) , INTENT(OUT) :: w
822  REAL(fp) , INTENT(IN OUT) :: pcoeff(0:,:)
823  TYPE(ASinterp_type) , INTENT(IN OUT) :: asi
824  ! Local variables
825  INTEGER :: k, l, m
826 
827 
828  ! Get the aerosol type LUT index
829  ! ------------------------------
830  k = aerosol_type_index( aerosol_type )
831 
832 
833  ! Find the Frequency and effective
834  ! radius indices for interpolation
835  ! --------------------------------
836  asi%f_int = max(min(aeroc%Frequency(aeroc%n_Wavelengths),frequency),aeroc%Frequency(1))
837  CALL find_index(aeroc%Frequency(:),asi%f_int,asi%i1,asi%i2, asi%f_outbound)
838  asi%f = aeroc%Frequency(asi%i1:asi%i2)
839 
840  asi%r_int = max(min(aeroc%Reff(aeroc%n_Radii,k),reff),aeroc%Reff(1,k))
841  CALL find_index(aeroc%Reff(:,k), asi%r_int, asi%j1,asi%j2, asi%r_outbound)
842  asi%r = aeroc%Reff(asi%j1:asi%j2,k)
843 
844 
845  ! Calculate the interpolating polynomials
846  ! ---------------------------------------
847  ! Frequency term
848  CALL lpoly( asi%f, asi%f_int, & ! Input
849  asi%wlp ) ! Output
850  ! Effective radius term
851  CALL lpoly( asi%r, asi%r_int, & ! Input
852  asi%xlp ) ! Output
853 
854 
855  ! Perform Interpolation
856  ! ---------------------
857  CALL interp_2d( aeroc%ke(asi%i1:asi%i2,asi%j1:asi%j2,k), asi%wlp, asi%xlp, ke )
858  CALL interp_2d( aeroc%w(asi%i1:asi%i2,asi%j1:asi%j2,k) , asi%wlp, asi%xlp, w )
859  IF (aerosolscatter%n_Phase_Elements > 0 .and. aerosolscatter%Include_Scattering ) THEN
860  pcoeff(0,1) = point_5
861  DO m = 1, aerosolscatter%n_Phase_Elements
862  DO l = 1, aerosolscatter%n_Legendre_Terms
863  CALL interp_2d( aeroc%pcoeff(asi%i1:asi%i2,asi%j1:asi%j2,k,l+aerosolscatter%lOffset,m), &
864  asi%wlp, asi%xlp, pcoeff(l,m) )
865  END DO
866  END DO
867  ELSE
868  ! Absorption coefficient
869  ke = ke * (one- w)
870  END IF
871 
872  END SUBROUTINE get_aerosol_opt
873 
874 
875  ! ---------------------------------------------
876  ! Subroutine to obtain the tangent-linear
877  ! bulk optical properties of a aerosol:
878  ! extinction coefficient (ke_TL),
879  ! scattereing coefficient (w_TL)
880  ! asymmetry factor (g_TL), and
881  ! spherical Legendre coefficients (pcoeff_TL)
882  ! ---------------------------------------------
883  SUBROUTINE get_aerosol_opt_tl(AerosolScatter_TL, & ! Input AerosolScatterTL structure
884  Aerosol_Type , & ! Input see CRTM_Aerosol_Define.f90
885  ke , & ! Input
886  w , & ! Input
887  Reff_TL , & ! Input TL effective radius (mm)
888  ke_tl , & ! Output TL extinction coefficient (=~ optical depth)
889  w_tl , & ! Output TL single scattering albedo
890  pcoeff_tl , & ! TL Output spherical Legendre coefficients
891  asi ) ! Input interpolation data
893  ! Arguments
894  TYPE(CRTM_AtmOptics_type), INTENT(IN) :: AerosolScatter_TL
895  INTEGER , INTENT(IN) :: Aerosol_Type
896  REAL(fp), INTENT(IN) :: w, ke, Reff_TL
897  REAL(fp), INTENT(OUT) :: ke_TL
898  REAL(fp), INTENT(OUT) :: w_TL
899  REAL(fp), INTENT(IN OUT) :: pcoeff_TL(0:,:)
900  TYPE(ASinterp_type), INTENT(IN) :: asi
901  ! Local variables
902  INTEGER :: k, l, m
903  REAL(fp) :: f_int_TL, r_int_TL
904  REAL(fp) :: f_TL(NPTS), r_TL(NPTS)
905  REAL(fp) :: z_TL(NPTS,NPTS)
906  TYPE(LPoly_type) :: wlp_TL, xlp_TL
907  REAL(fp), POINTER :: z(:,:) => null()
908 
909  ! Setup
910  ! -----
911  ! No TL output when all dimensions
912  ! are outside LUT bounds
913  IF ( asi%f_outbound .AND. asi%r_outbound ) THEN
914  ke_tl = zero
915  w_tl = zero
916  pcoeff_tl = zero
917  RETURN
918  END IF
919  ! The TL inputs
920  f_int_tl = zero
921  f_tl = zero
922  r_int_tl = reff_tl
923  r_tl = zero
924  z_tl = zero
925 
926 
927  ! Calculate the TL interpolating polynomials
928  ! ------------------------------------------
929  ! Frequency term (always zero. This is a placeholder for testing)
930  CALL lpoly_tl( asi%f, asi%f_int, & ! FWD Input
931  asi%wlp, & ! FWD Input
932  f_tl, f_int_tl, & ! TL Input
933  wlp_tl ) ! TL Output
934  ! Effective radius term
935  CALL lpoly_tl( asi%r, asi%r_int, & ! FWD Input
936  asi%xlp, & ! FWD Input
937  r_tl, r_int_tl, & ! TL Input
938  xlp_tl ) ! TL Output
939 
940 
941  ! Get the aerosol type LUT index
942  ! ------------------------------
943  k = aerosol_type_index( aerosol_type )
944 
945 
946  ! Perform Interpolation
947  ! ---------------------
948  ! Extinction coefficient
949  z => aeroc%ke(asi%i1:asi%i2,asi%j1:asi%j2,k)
950  CALL interp_2d_tl( z , asi%wlp, asi%xlp, & ! FWD Input
951  z_tl, wlp_tl , xlp_tl , & ! TL Input
952  ke_tl ) ! TL Output
953  ! Single scatter albedo
954  z => aeroc%w(asi%i1:asi%i2,asi%j1:asi%j2,k)
955  CALL interp_2d_tl( z , asi%wlp, asi%xlp, & ! FWD Input
956  z_tl, wlp_tl , xlp_tl , & ! TL Input
957  w_tl ) ! TL Output
958  ! Phase matrix coefficients
959  IF (aerosolscatter_tl%n_Phase_Elements > 0 .and. aerosolscatter_tl%Include_Scattering ) THEN
960  pcoeff_tl(0,1) = zero
961  DO m = 1, aerosolscatter_tl%n_Phase_Elements
962  DO l = 1, aerosolscatter_tl%n_Legendre_Terms
963  z => aeroc%pcoeff(asi%i1:asi%i2,asi%j1:asi%j2,k,l+aerosolscatter_tl%lOffset,m)
964  CALL interp_2d_tl( z , asi%wlp, asi%xlp, & ! FWD Input
965  z_tl, wlp_tl , xlp_tl , & ! TL Input
966  pcoeff_tl(l,m) ) ! TL Output
967  END DO
968  END DO
969  ELSE
970  ! Absorption coefficient
971  IF( w < one ) THEN
972  ke_tl = ke_tl * (one - w) - ke/(one -w) * w_tl
973  ELSE
974  ke_tl = zero
975  END IF
976  END IF
977  NULLIFY(z)
978 
979  END SUBROUTINE get_aerosol_opt_tl
980 
981 
982  ! ---------------------------------------------
983  ! Subroutine to obtain the adjoint
984  ! bulk optical properties of a aerosol:
985  ! extinction coefficient (ke_AD),
986  ! scattereing coefficient (w_AD)
987  ! asymmetry factor (g_AD), and
988  ! spherical Legendre coefficients (pcoeff_AD)
989  ! ---------------------------------------------
990  SUBROUTINE get_aerosol_opt_ad( AerosolScatter_AD, & ! Input AerosolScatter AD structure
991  Aerosol_Type , & ! Input see CRTM_Aerosol_Define.f90
992  ke , & ! Input
993  w , & ! Input
994  ke_AD , & ! AD Input extinction cross section
995  w_AD , & ! AD Input single scatter albedo
996  pcoeff_AD , & ! AD Input spherical Legendre coefficients
997  Reff_AD , & ! AD Output effective radius
998  asi ) ! Input interpolation data
999  ! Arguments
1000  TYPE(CRTM_AtmOptics_type), INTENT(IN) :: AerosolScatter_AD
1001  INTEGER , INTENT(IN) :: Aerosol_Type
1002  REAL(fp), INTENT(IN) :: ke, w
1003  REAL(fp), INTENT(IN OUT) :: ke_AD ! AD Input
1004  REAL(fp), INTENT(IN OUT) :: w_AD ! AD Input
1005  REAL(fp), INTENT(IN OUT) :: pcoeff_AD(0:,:) ! AD Input
1006  REAL(fp), INTENT(IN OUT) :: Reff_AD ! AD Output
1007  TYPE(ASinterp_type), INTENT(IN) :: asi
1008  ! Local variables
1009  INTEGER :: k, l, m
1010  REAL(fp) :: f_int_AD, r_int_AD
1011  REAL(fp) :: f_AD(NPTS), r_AD(NPTS)
1012  REAL(fp) :: z_AD(NPTS,NPTS)
1013  TYPE(LPoly_type) :: wlp_AD, xlp_AD
1014  REAL(fp), POINTER :: z(:,:) => null()
1015 
1016 
1017  ! Setup
1018  ! -----
1019  ! No AD output all dimensions
1020  ! are outside LUT bounds
1021  IF ( asi%f_outbound .AND. asi%r_outbound ) THEN
1022  reff_ad = zero
1023  ke_ad = zero
1024  w_ad = zero
1025  pcoeff_ad = zero
1026  RETURN
1027  END IF
1028  ! Initialise local adjoint variables
1029  f_int_ad = zero
1030  r_int_ad = zero
1031  f_ad = zero
1032  r_ad = zero
1033  z_ad = zero
1034  CALL clear_lpoly(wlp_ad)
1035  CALL clear_lpoly(xlp_ad)
1036 
1037 
1038  ! Get the aerosol type LUT index
1039  ! ------------------------------
1040  k = aerosol_type_index( aerosol_type )
1041 
1042 
1043  ! Perform interpolation
1044  ! ---------------------
1045  ! Phase matrix coefficients
1046  IF (aerosolscatter_ad%n_Phase_Elements > 0 .and. aerosolscatter_ad%Include_Scattering ) THEN
1047  DO m = 1, aerosolscatter_ad%n_Phase_Elements
1048  DO l = 1, aerosolscatter_ad%n_Legendre_Terms
1049  z => aeroc%pcoeff(asi%i1:asi%i2,asi%j1:asi%j2,k,l+aerosolscatter_ad%lOffset,m)
1050  CALL interp_2d_ad( z , asi%wlp, asi%xlp, & ! FWD Input
1051  pcoeff_ad(l,m) , & ! AD Input
1052  z_ad, wlp_ad, xlp_ad ) ! AD Output
1053  END DO
1054  END DO
1055  pcoeff_ad(0,1) = zero
1056  ELSE
1057  ! Absorption coefficient
1058  IF( w < one ) THEN
1059  w_ad = w_ad - ke/(one -w) * ke_ad
1060  ke_ad = ke_ad * (one - w)
1061  ELSE
1062  ke_ad = zero
1063  END IF
1064 
1065  END IF
1066 
1067  ! Single scatter albedo
1068  z => aeroc%w(asi%i1:asi%i2,asi%j1:asi%j2,k)
1069  CALL interp_2d_ad( z , asi%wlp, asi%xlp, & ! FWD Input
1070  w_ad , & ! AD Input
1071  z_ad, wlp_ad , xlp_ad ) ! AD Output
1072  ! Extinction coefficient
1073  z => aeroc%ke(asi%i1:asi%i2,asi%j1:asi%j2,k)
1074  CALL interp_2d_ad( z , asi%wlp, asi%xlp, & ! FWD Input
1075  ke_ad , & ! AD Input
1076  z_ad, wlp_ad , xlp_ad ) ! AD Output
1077  NULLIFY(z)
1078 
1079  ! Compute the AD of the interpolating polynomials
1080  ! -----------------------------------------------
1081  ! Efective radius term
1082  CALL lpoly_ad( asi%r, asi%r_int, & ! FWD Input
1083  asi%xlp, & ! FWD Input
1084  xlp_ad, & ! AD Input
1085  r_ad, r_int_ad ) ! AD Output
1086  ! Frequency term (always zero. This is a placeholder for testing)
1087  CALL lpoly_ad( asi%f, asi%f_int, & ! FWD Input
1088  asi%wlp, & ! FWD Input
1089  wlp_ad, & ! AD Input
1090  f_ad, f_int_ad ) ! AD Output
1091 
1092  ! The AD outputs
1093  ! --------------
1094  reff_ad = reff_ad + r_int_ad
1095 
1096  END SUBROUTINE get_aerosol_opt_ad
1097 
1098 END MODULE crtm_aerosolscatter
integer function, public crtm_compute_aerosolscatter_ad(Atm, AScat, AScat_AD, SensorIndex, ChannelIndex, Atm_AD, ASV)
subroutine get_aerosol_opt_tl(AerosolScatter_TL, Aerosol_Type, ke, w, Reff_TL, ke_TL, w_TL, pcoeff_TL, asi)
subroutine, public interp_3d_ad(z, ulp, vlp, wlp, z_int_AD, z_AD, ulp_AD, vlp_AD, wlp_AD)
subroutine, public interp_3d_tl(z, ulp, vlp, wlp, z_TL, ulp_TL, vlp_TL, wlp_TL, z_int_TL)
integer, parameter, public failure
real(fp), parameter, public onepointfive
real(fp), parameter, public zero
integer, parameter, public max_n_phase_elements
integer function, public crtm_compute_aerosolscatter_tl(Atm, AScat, Atm_TL, SensorIndex, ChannelIndex, AScat_TL, ASV)
integer, parameter, public fp
Definition: Type_Kinds.f90:124
character(*), parameter module_version_id
integer, parameter thirtytwo_streams
integer, parameter sixteen_streams
integer, parameter, public max_n_legendre_terms
logical, parameter, public hgphase
integer, parameter two_streams
subroutine, public clear_lpoly(p)
subroutine, public lpoly_ad(x, x_int, p, p_AD, x_AD, x_int_AD)
real(fp), parameter, public bs_threshold
subroutine get_aerosol_opt_ad(AerosolScatter_AD, Aerosol_Type, ke, w, ke_AD, w_AD, pcoeff_AD, Reff_AD, asi)
integer function, public crtm_compute_aerosolscatter(Atm, SensorIndex, ChannelIndex, AScat, ASV)
integer, parameter, public max_n_aerosols
type(aerosolcoeff_type), target, save, public aeroc
real(fp), parameter, public aerosol_content_threshold
real(fp), parameter, public one
integer function aerosol_type_index(Aerosol_Type)
recursive subroutine, public display_message(Routine_Name, Message, Error_State, Message_Log)
elemental subroutine, public asvar_destroy(self)
subroutine get_aerosol_opt(AerosolScatter, Frequency, Aerosol_Type, Reff, ke, w, pcoeff, asi)
subroutine, public interp_2d_ad(z, ulp, vlp, z_int_AD, z_AD, ulp_AD, vlp_AD)
subroutine, public lpoly(x, x_int, p)
integer, parameter, public npts
type(spccoeff_type), dimension(:), allocatable, save, public sc
elemental logical function, public asvar_associated(self)
integer, parameter, public max_n_layers
real(fp), parameter, public point_5
subroutine, public interp_3d(z, ulp, vlp, wlp, z_int)
#define max(a, b)
Definition: mosaic_util.h:33
integer, parameter four_streams
#define min(a, b)
Definition: mosaic_util.h:32
elemental subroutine, public asvar_create(self, n_Legendre_Terms, n_Phase_Elements, n_Layers, n_Aerosols)
integer, parameter six_streams
integer, parameter ml
integer, parameter eight_streams
subroutine, public lpoly_tl(x, x_int, p, x_TL, x_int_TL, p_TL)
subroutine, public interp_1d(z, ulp, z_int)
integer, parameter, public success
subroutine, public interp_2d_tl(z, ulp, vlp, z_TL, ulp_TL, vlp_TL, z_int_TL)
subroutine, public interp_2d(z, ulp, vlp, z_int)