FV3 Bundle
ufo_seaicefrac_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-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 handle ice concentration observations
7 
9 
10 use iso_c_binding
11 use ufo_vars_mod
13 use kinds
14 
15 implicit none
16 public :: ufo_seaicefrac_tlad
20 private
21 integer, parameter :: max_string=800
22 
23 !> Fortran derived type for sea ice fraction observation operator
25  integer :: ncat = -1 !< number of ice categories
26 end type ufo_seaicefrac_tlad
27 
28 
29 ! ------------------------------------------------------------------------------
30 
31 contains
32 
33 
34 ! ------------------------------------------------------------------------------
35 
36 subroutine ufo_seaicefrac_tlad_settraj(self, geovals)
37 implicit none
38 type(ufo_seaicefrac_tlad), intent(inout) :: self
39 type(ufo_geovals), intent(in) :: geovals
40 
41 character(len=*), parameter :: myname_="ufo_seaicefrac_tlad_settraj"
42 character(max_string) :: err_msg
43 
44 type(ufo_geoval), pointer :: geoval
45 
46 ! since observation operator is linear, only need to save the number
47 ! of ice categories here, don't care about trajectory itself
48 
49 ! check if sea ice fraction variables is in geovals and get it
50 call ufo_geovals_get_var(geovals, var_seaicefrac, geoval)
51 
52 self%ncat = geoval%nval
53 
54 end subroutine ufo_seaicefrac_tlad_settraj
55 
56 
57 ! ------------------------------------------------------------------------------
58 
59 subroutine ufo_seaicefrac_simobs_tl(self, geovals, hofx)
60 implicit none
61 type(ufo_seaicefrac_tlad), intent(in) :: self
62 type(ufo_geovals), intent(in) :: geovals
63 real(c_double), intent(inout) :: hofx(:)
64 
65 character(len=*), parameter :: myname_="ufo_seaicefrac_simobs_tl"
66 character(max_string) :: err_msg
67 
68 integer :: iobs
69 type(ufo_geoval), pointer :: geoval
70 
71 print *, myname_, ' nobs: ', geovals%nobs, size(hofx,1)
72 
73 ! check if nobs is consistent in geovals & hofx
74 if (geovals%nobs /= size(hofx,1)) then
75  write(err_msg,*) myname_, ' error: nobs inconsistent!'
76  call abor1_ftn(err_msg)
77 endif
78 
79 ! check if sea ice fraction variables is in geovals and get it
80 call ufo_geovals_get_var(geovals, var_seaicefrac, geoval)
81 
82 ! total sea ice fraction obs operator
83 do iobs = 1, size(hofx,1)
84  hofx(iobs) = sum(geoval%vals(:,iobs))
85 enddo
86 
87 end subroutine ufo_seaicefrac_simobs_tl
88 
89 ! ------------------------------------------------------------------------------
90 
91 subroutine ufo_seaicefrac_simobs_ad(self, geovals, hofx)
92 implicit none
93 type(ufo_seaicefrac_tlad), intent(in) :: self
94 type(ufo_geovals), intent(inout) :: geovals
95 real(c_double), intent(inout) :: hofx(:)
96 
97 character(len=*), parameter :: myname_="ufo_seaicefrac_simobs_ad"
98 character(max_string) :: err_msg
99 
100 integer :: iobs
101 type(ufo_geoval), pointer :: geoval
102 
103 ! check if nobs is consistent in geovals & hofx
104 if (geovals%nobs /= size(hofx,1)) then
105  write(err_msg,*) myname_, ' error: nobs inconsistent!'
106  call abor1_ftn(err_msg)
107 endif
108 
109 if (.not. geovals%linit ) geovals%linit=.true.
110 
111 ! check if sea ice fraction variables is in geovals and get it
112 call ufo_geovals_get_var(geovals, var_seaicefrac, geoval)
113 
114 if (.not.(allocated(geoval%vals))) then
115  if (self%ncat < 1) then
116  write(err_msg,*) myname_, ' unknown number of categories'
117  call abor1_ftn(err_msg)
118  endif
119  allocate(geoval%vals(self%ncat,size(hofx,1)))
120 end if
121 
122 ! backward sea ice fraction obs operator
123 geoval%vals=0.0
124 do iobs = 1, size(hofx,1)
125  geoval%vals(:,iobs) = geoval%vals(:,iobs) + hofx(iobs)
126 enddo
127 
128 end subroutine ufo_seaicefrac_simobs_ad
129 
130 end module ufo_seaicefrac_tlad_mod
subroutine, public ufo_geovals_get_var(self, varname, geoval, status)
Fortran module to handle ice concentration observations.
integer, parameter max_string
subroutine, public ufo_seaicefrac_simobs_tl(self, geovals, hofx)
type to hold interpolated fields required by the obs operators
Fortran derived type for sea ice fraction observation operator.
subroutine, public ufo_seaicefrac_tlad_settraj(self, geovals)
character(len=maxvarlen), public var_seaicefrac
type to hold interpolated field for one variable, one observation
subroutine, public ufo_seaicefrac_simobs_ad(self, geovals, hofx)