UFO
ufo_aodgeos_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 observation operator
7 
9 
10  use iso_c_binding
11  use kinds
12  use geos_mieobs_mod
13  use oops_variables_mod
14  use ufo_vars_mod
15  use ufo_basis_mod, only: ufo_basis
16 
17  implicit none
18  private
19  integer, parameter :: max_string=800
20 
21 !> Fortran derived type for the observation type
22  type, public :: ufo_aodgeos
23  type(oops_variables), public :: geovars
24  type(oops_variables), public :: obsvars
25  integer, public :: ntracers
26  real(kind_real), public, allocatable :: wavelength(:)
27  character(len=maxvarlen),public :: rcfile
28  contains
29  procedure :: setup => ufo_aodgeos_setup
30  procedure :: simobs => ufo_aodgeos_simobs
31  final :: destructor
32  end type ufo_aodgeos
33 
34 !> Default variables required from model
35  character(len=maxvarlen), dimension(2), parameter :: varindefault = (/var_delp, var_rh/)
36 
37 contains
38 
39 ! ------------------------------------------------------------------------------
40 
41 subroutine ufo_aodgeos_setup(self, f_conf)
42 use fckit_configuration_module, only: fckit_configuration
43 implicit none
44 
45 class(ufo_aodgeos), intent(inout) :: self
46 type(fckit_configuration), intent(in) :: f_conf
47 
48 !Locals
49 integer :: iq, nvars
50 character(kind=c_char,len=:), allocatable :: tracer_variables(:)
51 character(len=:), allocatable :: str
52 
53  ! Fill in geovars: variables we need from the model
54  ! Need slots for RH and delp
55  ! Let user choose specific aerosols needed in aod calculation.
56 
57  call f_conf%get_or_die("tracer_geovals",tracer_variables)
58  self%ntracers = f_conf%get_size("tracer_geovals")
59 
60  do iq = 1, self%ntracers
61  call self%geovars%push_back(tracer_variables(iq)) ! aer MR
62  enddo
63  call self%geovars%push_back(varindefault) ! delp and rh (for concentration)
64  deallocate(tracer_variables)
65 
66  ! Size of variables (number of obs type (wavelength for AOD))
67  nvars = self%obsvars%nvars()
68 
69  ! List of wavelengths for the calculation of aod
70  allocate(self%wavelength(nvars))
71  call f_conf%get_or_die("wavelengths", self%wavelength)
72 
73  ! RC File needed for ChemBase
74  call f_conf%get_or_die("RCFile",str)
75  self%rcfile = str
76  deallocate(str)
77 
78 end subroutine ufo_aodgeos_setup
79 
80 ! ------------------------------------------------------------------------------
81 
82 subroutine destructor(self)
83 implicit none
84 type(ufo_aodgeos), intent(inout) :: self
85 
86  if (allocated(self%wavelength)) deallocate(self%wavelength)
87 
88 end subroutine destructor
89 
90 ! ------------------------------------------------------------------------------
91 
92 subroutine ufo_aodgeos_simobs(self, geovals, obss, nvars, nlocs, hofx)
93 use kinds
94 use ufo_constants_mod, only: grav
96 use iso_c_binding
97 
98 implicit none
99 class(ufo_aodgeos), intent(in) :: self
100 integer, intent(in) :: nvars, nlocs
101 type(ufo_geovals), intent(in) :: geovals
102 real(c_double), intent(inout) :: hofx(nvars, nlocs)
103 type(c_ptr), value, intent(in) :: obss
104 
105 ! Local variables
106 type(ufo_geoval), pointer :: aer_profile
107 type(ufo_geoval), pointer :: delp_profile
108 type(ufo_geoval), pointer :: rh_profile
109 integer :: nlayers, rc, iq
110 
111 real(kind_real), dimension(:,:,:), allocatable :: qm ! aer mass mix ratio(kg/kg) that becomes concentration (*delp/g) profiles at obs loc
112 real(kind_real), dimension(:,:), allocatable :: rh ! relative humidity profile interp at obs loc
113 real(kind_real), dimension(:,:), allocatable :: delp ! air pressure thickness profiles at obs loc
114 
115 character(len=MAXVARLEN) :: geovar
116 character(len=MAXVARLEN), dimension(:), allocatable:: tracer_name
117 
118 
119  ! Get delp and rh from model interp at obs loc (from geovals)
120  call ufo_geovals_get_var(geovals, var_delp, delp_profile)
121  nlayers = delp_profile%nval ! number of model layers
122 
123  allocate(delp(nlayers,nlocs))
124  delp = delp_profile%vals
125 
126  ! Get RH from geovals
127  allocate(rh(nlayers,nlocs))
128  call ufo_geovals_get_var(geovals, var_rh, rh_profile)
129  rh = rh_profile%vals
130 
131  ! Get Aer profiles interpolated at obs loc
132  allocate(qm(self%ntracers, nlayers, nlocs))
133  allocate(tracer_name(self%ntracers))
134  do iq = 1, self%ntracers
135  geovar = self%geovars%variable(iq) !self%geovars contains tracers
136  tracer_name(iq) = geovar
137  call ufo_geovals_get_var(geovals, geovar, aer_profile)
138  qm(iq,:,:) = aer_profile%vals ! mass mixing ratio
139  qm(iq,:,:) = qm(iq,:,:) * delp / grav ! aer concentration (kg/m2)
140  enddo
141 
142  ! Call observation operator code
143  ! -----------------------------
144  hofx(:,:) = 0.0
145  call get_geos_aod(nlayers, nlocs, nvars, self%ntracers, self%rcfile, &
146  self%wavelength, tracer_name, qm, rh, &
147  aod_tot = hofx, rc = rc) !self%varin includes rh and delp
148 
149  ! cleanup memory
150  ! --------
151  deallocate(qm, rh, delp, tracer_name)
152 
153 end subroutine ufo_aodgeos_simobs
154 
155 ! ------------------------------------------------------------------------------
156 
157 end module ufo_aodgeos_mod
Fortran module for aodgeos observation operator.
subroutine ufo_aodgeos_simobs(self, geovals, obss, nvars, nlocs, hofx)
subroutine destructor(self)
subroutine ufo_aodgeos_setup(self, f_conf)
integer, parameter max_string
character(len=maxvarlen), dimension(2), parameter varindefault
Default variables required from model.
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 observation type.
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators