UFO
ufo_luts_utils_mod.F90
Go to the documentation of this file.
1 ! (c) copyright 2018 ucar
2 !
3 ! this software is licensed under the terms of the apache licence version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/license-2.0.
5 
6 !> fortran module to provide code shared between nonlinear and tlm/adm radiance calculations
7 
9 
10  USE fckit_configuration_module, ONLY: fckit_configuration
11  USE iso_c_binding
12  USE kinds
13 
14  USE crtm_module
15 
16  USE ufo_vars_mod
18  USE obsspace_mod
19 
21  USE ufo_crtm_utils_mod, ONLY: grav,rv_rd,&
24 
25  IMPLICIT NONE
26  PRIVATE
27 
28  PUBLIC luts_conf
29  PUBLIC luts_conf_setup
30  PUBLIC luts_conf_delete
31  PUBLIC :: calculate_aero_layers
32 
33  INTEGER, PARAMETER, PUBLIC :: max_string=800
34 
35 
36 !type for general config
37  TYPE luts_conf
38  CHARACTER(len=maxvarlen) :: aerosol_option
39  CHARACTER(len=max_string) :: rcfile
40  CHARACTER(len=255), ALLOCATABLE :: sensor_id(:)
41  INTEGER :: n_sensors
42  CHARACTER(len=255) :: endian_type
43  CHARACTER(len=255) :: coefficient_path
44  END TYPE luts_conf
45 
46 CONTAINS
47 
48 ! ------------------------------------------------------------------------------
49 
50  SUBROUTINE luts_conf_setup(conf, f_confopts, f_confoper)
51 
52  IMPLICIT NONE
53  TYPE(luts_conf), INTENT(inout) :: conf
54  TYPE(fckit_configuration), INTENT(in) :: f_confopts
55  TYPE(fckit_configuration), INTENT(in) :: f_confoper
56 
57  CHARACTER(*), PARAMETER :: routine_name = 'luts_conf_setup'
58  CHARACTER(len=:), ALLOCATABLE :: str
59 
60  CALL f_confopts%get_or_die("AerosolOption",str)
61  conf%aerosol_option = upper2lower(str)
62  CALL f_confopts%get_or_die("RCFile",str)
63  conf%rcfile = str
64  conf%n_sensors = 1
65  ALLOCATE(conf%sensor_id(conf%n_sensors))
66  CALL f_confopts%get_or_die("Sensor_ID",str)
67  conf%sensor_id(conf%n_sensors) = str
68  CALL f_confopts%get_or_die("EndianType",str)
69  conf%endian_type = str
70  CALL f_confopts%get_or_die("CoefficientPath",str)
71  conf%coefficient_path = str
72 
73  DEALLOCATE(str)
74 
75  END SUBROUTINE luts_conf_setup
76 
77 ! -----------------------------------------------------------------------------
78 
79  SUBROUTINE luts_conf_delete(conf)
80 
81  IMPLICIT NONE
82  TYPE(luts_conf), INTENT(inout) :: conf
83 
84  IF (ALLOCATED(conf%sensor_id)) DEALLOCATE(conf%sensor_id)
85 
86  END SUBROUTINE luts_conf_delete
87 
88 ! ------------------------------------------------------------------------------
89 
90  SUBROUTINE calculate_aero_layers(aerosol_option,&
91  &n_aerosols, n_profiles, n_layers, &
92  &geovals, aero_layers, rh, layer_factors)
93 
94  IMPLICIT NONE
95 
96  CHARACTER(*), INTENT(in) :: aerosol_option
97  INTEGER, INTENT(in) :: n_aerosols, n_profiles, n_layers
98  TYPE(ufo_geovals), INTENT(in) :: geovals
99  REAL(kind_real), OPTIONAL, INTENT(out) :: &
100  &aero_layers(n_aerosols,n_layers,n_profiles)
101  REAL(kind_real), OPTIONAL, INTENT(out) :: rh(n_layers,n_profiles)
102  REAL(kind_real), OPTIONAL, INTENT(out) :: layer_factors(n_layers,n_profiles)
103 
104  TYPE(luts_conf) :: conf
105 
106 ! local variables
107  INTEGER :: m, ivar,k
108  TYPE(ufo_geoval), POINTER :: geoval
109  CHARACTER(max_string) :: err_msg
110  CHARACTER(len=maxvarlen) :: varname
111  REAL(kind_real), DIMENSION(n_layers,n_profiles) :: t,sphum,pmid
112  REAL(kind_real), DIMENSION(n_layers+1,n_profiles) :: pint
113  CHARACTER(len=maxvarlen), ALLOCATABLE :: var_aerosols(:)
114  REAL(kind_real) :: factors(n_layers,n_profiles)
115 
116  CALL ufo_geovals_get_var(geovals, var_ts, geoval)
117  IF (geoval%nval /= n_layers) THEN
118  WRITE(err_msg,*) 'get_atm_aero_data error: layers inconsistent!'
119  CALL abor1_ftn(err_msg)
120  ENDIF
121  t=geoval%vals
122 
123  CALL ufo_geovals_get_var(geovals, var_prs, geoval)
124  pmid=geoval%vals
125 
126  CALL ufo_geovals_get_var(geovals, var_prsi, geoval)
127  pint=geoval%vals
128 
129  CALL ufo_geovals_get_var(geovals, var_q, geoval)
130  sphum=geoval%vals
131 
132  IF (trim(aerosol_option) /= "aerosols_gocart_merra_2") THEN
133  WRITE(err_msg,*) 'this aerosol not implemented - check later'
134  CALL abor1_ftn(err_msg)
135  ENDIF
136 
137  CALL assign_aerosol_names(aerosol_option,var_aerosols)
138 
139  DO k=1,n_layers
140  DO m = 1, n_profiles
141 !correct for mixing ratio factor ugkg_kgm2
142 !being calculated from dry pressure, cotton eq. (2.4)
143 !p_dry=p_total/(1+r_v/r_d*mixing_ratio)
144  factors(k,m)=1e-9_kind_real*(pint(k+1,m)-pint(k,m))/grav/&
145  &(1_kind_real+rv_rd*sphum(k,m)/(1_kind_real-sphum(k,m)))
146  ENDDO
147  ENDDO
148 
149  IF ( PRESENT(aero_layers) ) THEN
150  DO ivar=1,n_aerosols
151  CALL ufo_geovals_get_var(geovals, var_aerosols(ivar), geoval)
152  aero_layers(ivar,:,:)=geoval%vals*factors
153  ENDDO
154  ENDIF
155 
156  IF ( PRESENT(rh) ) THEN
157  CALL qsmith(t,sphum,pmid,rh)
158  WHERE (rh > 1_kind_real) rh=1_kind_real
159  ENDIF
160 
161  IF ( PRESENT(layer_factors) ) layer_factors=factors
162 
163  END SUBROUTINE calculate_aero_layers
164 
165 ! ------------------------------------------------------------------------------
166 
167 END MODULE ufo_luts_utils_mod
168 
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
real(kind_real), parameter, public grav
real(kind_real), parameter, public aerosol_concentration_minvalue_layer
character(len(str)) function, public upper2lower(str)
real(kind_real), parameter, public aerosol_concentration_minvalue
subroutine, public assign_aerosol_names(aerosol_option, var_aerosols)
real(kind_real), parameter, public rv_rd
subroutine, public ufo_geovals_get_var(self, varname, geoval)
fortran module to provide code shared between nonlinear and tlm/adm radiance calculations
subroutine, public luts_conf_delete(conf)
integer, parameter, public max_string
subroutine, public calculate_aero_layers(aerosol_option, n_aerosols, n_profiles, n_layers, geovals, aero_layers, rh, layer_factors)
subroutine, public luts_conf_setup(conf, f_confopts, f_confoper)
character(len=maxvarlen), parameter, public var_prsi
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_ts
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators