UFO
ufo_aodgeos_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2019 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 for aodgeos tl/ad observation operator
7 
9 
10  use iso_c_binding
11 
12  use kinds
13 
15  use ufo_vars_mod
16  use ufo_constants_mod, only: grav
17 
18  use geos_mieobs_mod
19  use oops_variables_mod
20 
21 
22  implicit none
23  private
24  integer, parameter :: max_string=800
25 
26  !> Fortran derived type for the tl/ad observation operator
27  type, public :: ufo_aodgeos_tlad
28  integer :: nlocs, nlayers, ntracers, nvars
29  type(oops_variables), public :: obsvars
30  type(oops_variables), public :: geovars
31  real(kind=kind_real), allocatable :: bext(:,:,:,:)
32  real(kind=kind_real), public, allocatable :: wavelength(:)
33  character(len=maxvarlen),public:: rcfile
34  real(kind=kind_real), dimension(:,:), allocatable :: delp(:,:)
35  contains
36  procedure :: setup => ufo_aodgeos_tlad_setup
37  procedure :: delete => ufo_aodgeos_tlad_delete
38  procedure :: settraj => ufo_aodgeos_tlad_settraj
39  procedure :: simobs_tl => ufo_aodgeos_simobs_tl
40  procedure :: simobs_ad => ufo_aodgeos_simobs_ad
41  end type ufo_aodgeos_tlad
42 
43 contains
44 
45 ! ------------------------------------------------------------------------------
46 
47 subroutine ufo_aodgeos_tlad_setup(self, f_conf)
48 use fckit_configuration_module, only: fckit_configuration
49 implicit none
50 class(ufo_aodgeos_tlad), intent(inout) :: self
51 type(fckit_configuration), intent(in) :: f_conf
52 
53 !Locals
54 integer :: iq
55 character(len=:), allocatable :: tracer_variables(:)
56 character(len=:), allocatable :: str
57 
58  ! Let user choose specific aerosols needed.
59  call f_conf%get_or_die("tracer_geovals",tracer_variables)
60  self%ntracers = f_conf%get_size("tracer_geovals")
61  do iq = 1, self%ntracers
62  call self%geovars%push_back(tracer_variables(iq)) ! aer MR
63  enddo
64  deallocate(tracer_variables)
65 
66  ! Size of variables (number of obs type (wavelength for AOD))
67  self%nvars = self%obsvars%nvars()
68 
69  allocate(self%wavelength(self%nvars))
70  call f_conf%get_or_die("wavelengths", self%wavelength)
71 
72 ! RC File for ChemBase
73  call f_conf%get_or_die("RCFile",str)
74  self%rcfile = str
75  deallocate(str)
76 
77 end subroutine ufo_aodgeos_tlad_setup
78 
79 ! ------------------------------------------------------------------------------
80 
81 subroutine ufo_aodgeos_tlad_delete(self)
82 
83 implicit none
84 class(ufo_aodgeos_tlad), intent(inout) :: self
85  if (allocated(self%bext)) deallocate(self%bext)
86  if (allocated(self%delp)) deallocate(self%delp)
87  if (allocated(self%wavelength)) deallocate(self%wavelength)
88 end subroutine ufo_aodgeos_tlad_delete
89 
90 ! ------------------------------------------------------------------------------
91 
92 subroutine ufo_aodgeos_tlad_settraj(self, geovals, obss)
93 use obsspace_mod
94 implicit none
95 class(ufo_aodgeos_tlad), intent(inout) :: self
96 
97 type(ufo_geovals), intent(in) :: geovals
98 type(c_ptr), value, intent(in) :: obss
99 
100 ! Local variables
101 type(ufo_geoval), pointer :: aer_profile
102 type(ufo_geoval), pointer :: delp_profile
103 type(ufo_geoval), pointer :: rh_profile
104 integer :: rc, iq
105 character(len=MAXVARLEN) :: geovar
106 character(len=MAXVARLEN), dimension(:), allocatable:: tracer_name
107 
108 real(kind=kind_real), dimension(:,:), allocatable :: rh
109 real(kind=kind_real), dimension(:,:,:), allocatable :: qm ! aer concentration (kg/kg *delp/g) profiles at obs loc
110 
111  ! Get number of locations
112  self%nlocs = obsspace_get_nlocs(obss)
113 
114  ! Get delp and rh from model interp at obs loc (from geovals)
115  call ufo_geovals_get_var(geovals, var_delp, delp_profile)
116  self%nlayers = delp_profile%nval ! number of model layers
117  allocate(self%delp(self%nlayers,self%nlocs))
118  self%delp = delp_profile%vals
119 
120  ! Get RH from geovals
121  allocate(rh(self%nlayers,self%nlocs))
122  call ufo_geovals_get_var(geovals, var_rh, rh_profile)
123  rh = rh_profile%vals
124 
125  ! Get Aer profiles interpolated at obs loc
126  allocate(qm(self%ntracers, self%nlayers, self%nlocs))
127  allocate(tracer_name(self%ntracers))
128  do iq = 1, self%ntracers
129  geovar = self%geovars%variable(iq) !self%geovars contains tracers
130  tracer_name(iq) = geovar
131  call ufo_geovals_get_var(geovals, geovar, aer_profile)
132  qm(iq,:,:) = aer_profile%vals
133  qm(iq,:,:) = qm(iq,:,:) * self%delp / grav
134  enddo
135 
136  allocate(self%bext(self%nlayers, self%nvars, self%ntracers, self%nlocs)) !mass extinction efficiency
137  call get_geos_aod(self%nlayers, self%nlocs, self%nvars, self%ntracers, self%rcfile, &
138  self%wavelength, tracer_name, qm, rh, ext=self%bext, rc = rc)
139 
140  deallocate(rh)
141  deallocate(qm)
142  deallocate(self%wavelength)
143  deallocate(tracer_name)
144 
145 end subroutine ufo_aodgeos_tlad_settraj
146 
147 ! ------------------------------------------------------------------------------
148 
149 subroutine ufo_aodgeos_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
150 
151 implicit none
152 class(ufo_aodgeos_tlad), intent(in) :: self
153 type(ufo_geovals), intent(in) :: geovals
154 type(c_ptr), value, intent(in) :: obss
155 integer, intent(in) :: nvars, nlocs
156 real(c_double), intent(inout) :: hofx(nvars, nlocs)
157 
158 integer :: iq
159 real(kind_real), dimension(:,:,:), allocatable :: qm_tl
160 type(ufo_geoval), pointer :: aer_profile
161 
162 character(len=MAXVARLEN) :: geovar
163 
164  ! Get Aer profiles interpolated at obs loc
165  allocate(qm_tl(self%ntracers, self%nlayers, nlocs))
166 
167  do iq = 1, self%ntracers
168  geovar = self%geovars%variable(iq) !self%geovars contains tracers
169  call ufo_geovals_get_var(geovals, geovar, aer_profile)
170  qm_tl(iq,:,:) = aer_profile%vals ! aer mass mixing ratio
171  qm_tl(iq,:,:) = qm_tl(iq,:,:) * self%delp / grav ! aer concentration
172  enddo
173 
174  call get_geos_aod_tl(self%nlayers,self%nlocs, nvars, self%ntracers, self%bext, qm_tl, aod_tot_tl=hofx)
175  deallocate(qm_tl)
176 
177 
178 end subroutine ufo_aodgeos_simobs_tl
179 
180 ! ------------------------------------------------------------------------------
181 
182 subroutine ufo_aodgeos_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
183 
184 implicit none
185 class(ufo_aodgeos_tlad), intent(in) :: self
186 type(ufo_geovals), intent(inout) :: geovals
187 type(c_ptr), value, intent(in) :: obss
188 integer, intent(in) :: nvars, nlocs
189 real(c_double), intent(in) :: hofx(nvars, nlocs)
190 
191 integer :: iq
192 real(kind_real), dimension(:,:,:), allocatable :: qm_ad
193 type(ufo_geoval), pointer :: aer_profile
194 character(len=MAXVARLEN) :: geovar
195 
196  allocate(qm_ad(self%ntracers, self%nlayers, nlocs))
197 
198  call get_geos_aod_ad(self%nlayers, nlocs, nvars, self%ntracers, self%bext, hofx, qm_ad)
199 
200  do iq = self%ntracers,1,-1
201 
202  geovar = self%geovars%variable(iq) !self%geovars contains tracers
203  call ufo_geovals_get_var(geovals, geovar, aer_profile)
204  qm_ad(iq,:,:) = qm_ad(iq,:,:) * self%delp / grav
205  aer_profile%vals = qm_ad(iq,:,:)
206 
207  enddo
208 
209  deallocate(qm_ad)
210 
211 end subroutine ufo_aodgeos_simobs_ad
212 
213 ! ------------------------------------------------------------------------------
214 
215 
216 end module ufo_aodgeos_tlad_mod
Fortran module for aodgeos tl/ad observation operator.
subroutine ufo_aodgeos_tlad_setup(self, f_conf)
subroutine ufo_aodgeos_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
subroutine ufo_aodgeos_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
integer, parameter max_string
subroutine ufo_aodgeos_tlad_settraj(self, geovals, obss)
subroutine ufo_aodgeos_tlad_delete(self)
real(kind_real), parameter, public grav
subroutine, public ufo_geovals_get_var(self, varname, geoval)
character(len=maxvarlen), parameter, public var_rh
character(len=maxvarlen), parameter, public var_delp
Fortran derived type for the tl/ad observation operator.
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators