UFO
ufo_locs_mod.F90
Go to the documentation of this file.
1 !
2 ! (C) Copyright 2017 UCAR
3 !
4 ! This software is licensed under the terms of the Apache Licence Version 2.0
5 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6 
7 !> Fortran module handling observation locations
8 
10 
11 use datetime_mod
12 use iso_c_binding
13 use kinds
14 use fckit_log_module, only : fckit_log
15 use obsspace_mod
16 
17 implicit none
18 private
21 
22 ! --------------------------------------------------------------------------------------------------
23 
24 !> Fortran derived type to hold observation locations
25 type :: ufo_locs
26  integer :: nlocs
27  integer :: max_indx
28  real(kind_real), allocatable, dimension(:) :: lat !< latitude
29  real(kind_real), allocatable, dimension(:) :: lon !< longitude
30  type(datetime), allocatable, dimension(:) :: time !< obs-time
31  integer, allocatable, dimension(:) :: indx !< indices of locations in the full [geovals] array
32 end type ufo_locs
33 
34 ! --------------------------------------------------------------------------------------------------
35 contains
36 ! --------------------------------------------------------------------------------------------------
37 
38 subroutine ufo_locs_create(self, obss, nlocs, lats, lons)
39 implicit none
40 type(ufo_locs), intent(inout) :: self
41 type(c_ptr), value, intent(in) :: obss
42 integer, intent(in) :: nlocs
43 real(kind_real), intent(in) :: lats(nlocs)
44 real(kind_real), intent(in) :: lons(nlocs)
45 
46 integer :: n
47 character(len=20), allocatable :: fstring(:)
48 type(datetime), dimension(:), allocatable :: date_time
49 self%nlocs = nlocs
50 self%max_indx = nlocs
51 allocate(self%lat(nlocs), self%lon(nlocs), self%time(nlocs))
52 allocate(self%indx(nlocs))
53 self%lat(:) = lats(:)
54 self%lon(:) = lons(:)
55 
56 allocate(fstring(nlocs))
57 
58 if (obsspace_has(obss,"MetaData", "datetime")) then
59  allocate(date_time(nlocs))
60  call obsspace_get_db(obss, "MetaData", "datetime", date_time)
61  do n = 1, self%nlocs
62  call datetime_to_string(date_time(n), fstring(n))
63  enddo
64  deallocate(date_time)
65 else
66  fstring(:) = "9999-09-09T09:09:09Z"
67 endif
68 
69 do n = 1, self%nlocs
70  call datetime_create(fstring(n), self%time(n))
71 enddo
72 do n = 1, self%nlocs
73  self%indx(n) = n
74 enddo
75 
76 end subroutine ufo_locs_create
77 
78 ! --------------------------------------------------------------------------------------------------
79 
80 subroutine ufo_locs_setup(self, nlocs)
81 implicit none
82 type(ufo_locs), intent(inout) :: self
83 integer, intent(in) :: nlocs
84 
85 character(len=20) :: fstring
86 integer :: n
87 
88 call ufo_locs_delete(self)
89 
90 self%max_indx = nlocs
91 self%nlocs = nlocs
92 allocate(self%lat(nlocs), self%lon(nlocs), self%time(nlocs), self%indx(nlocs))
93 self%lat(:) = 0.0
94 self%lon(:) = 0.0
95 fstring="9999-09-09T09:09:09Z"
96 do n = 1, self%nlocs
97  call datetime_create(fstring, self%time(n))
98 enddo
99 self%indx(:) = 0
100 
101 end subroutine ufo_locs_setup
102 
103 ! --------------------------------------------------------------------------------------------------
104 
105 subroutine ufo_locs_copy(self, other)
106 
107 implicit none
108 type(ufo_locs), intent(inout) :: self
109 type(ufo_locs), intent(in) :: other
110 
111 self%nlocs = other%nlocs
112 self%max_indx = other%max_indx
113 
114 allocate(self%lat (self%nlocs))
115 allocate(self%lon (self%nlocs))
116 allocate(self%time(self%nlocs))
117 allocate(self%indx(self%nlocs))
118 
119 self%lat = other%lat
120 self%lon = other%lon
121 self%time = other%time
122 self%indx = other%indx
123 
124 end subroutine ufo_locs_copy
125 
126 ! --------------------------------------------------------------------------------------------------
127 
128 subroutine ufo_locs_delete(self)
129 implicit none
130 type(ufo_locs), intent(inout) :: self
131 
132 if (allocated(self%lat)) deallocate(self%lat)
133 if (allocated(self%lon)) deallocate(self%lon)
134 if (allocated(self%time)) deallocate(self%time)
135 if (allocated(self%indx)) deallocate(self%indx)
136 self%nlocs = 0
137 self%max_indx = -1 ! not set
138 
139 end subroutine ufo_locs_delete
140 
141 ! --------------------------------------------------------------------------------------------------
142 
143 subroutine ufo_locs_concatenate(self, other)
144 implicit none
145 type(ufo_locs), intent(inout) :: self
146 type(ufo_locs), intent(in) :: other
147 
148 type(ufo_locs) :: temp_self, temp_other
149 character(255) :: message
150 integer :: n
151 
152 if ((self%max_indx < 0) .OR. (other%max_indx < 0)) then
153  write(message,'(A, A, I6, A, I6)') &
154  'ufo_locs_concatenate: either self or other needs to be constructed valid indices ', &
155  ' self%max_indx =', self%max_indx, ' other%max_indx = ', other%max_indx
156  call fckit_log%info(message)
157  stop
158 end if
159 
160 ! make a temporary copy of self
161 temp_self%nlocs = self%nlocs
162 temp_self%max_indx = self%max_indx
163 allocate(temp_self%lat(self%nlocs), temp_self%lon(self%nlocs), &
164  temp_self%time(self%nlocs), temp_self%indx(self%nlocs))
165 
166 temp_self%lat(:) = self%lat(:)
167 temp_self%lon(:) = self%lon(:)
168 temp_self%indx(:) = self%indx(:)
169 do n = 1, self%nlocs
170  temp_self%time(n) = self%time(n)
171 end do
172 
173 ! make a temporary copy of other
174 temp_other%nlocs = other%nlocs
175 temp_other%max_indx = other%max_indx
176 allocate(temp_other%lat(other%nlocs), temp_other%lon(other%nlocs), &
177  temp_other%time(other%nlocs), temp_other%indx(other%nlocs))
178 
179 temp_other%lat(:) = other%lat(:)
180 temp_other%lon(:) = other%lon(:)
181 temp_other%indx(:) = other%indx(:)
182 do n = 1, other%nlocs
183  temp_other%time(n) = other%time(n)
184 end do
185 
186 ! deallocate self
187 call ufo_locs_delete(self)
188 
189 ! reallocate self with combined concatenation
190 self%nlocs = temp_self%nlocs + temp_other%nlocs
191 allocate(self%lat(self%nlocs), self%lon(self%nlocs), &
192  self%time(self%nlocs), self%indx(self%nlocs))
193 
194 self%lat(1:temp_self%nlocs) = temp_self%lat(1:temp_self%nlocs)
195 self%lat(temp_self%nlocs+1:) = temp_other%lat(1:temp_other%nlocs)
196 
197 self%lon(1:temp_self%nlocs) = temp_self%lon(1:temp_self%nlocs)
198 self%lon(temp_self%nlocs+1:) = temp_other%lon(1:temp_other%nlocs)
199 
200 do n = 1, temp_self%nlocs
201  self%time(n) = temp_self%time(n)
202 end do
203 do n = 1, temp_other%nlocs
204  self%time(temp_self%nlocs + n) = temp_other%time(n)
205 end do
206 
207 self%indx(1:temp_self%nlocs) = temp_self%indx(1:temp_self%nlocs)
208 do n = 1, temp_other%nlocs
209  self%indx(temp_self%nlocs + n) = temp_other%indx(n) + &
210  temp_self%max_indx
211 end do
212 self%max_indx = temp_self%max_indx + temp_other%max_indx
213 
214 call ufo_locs_delete(temp_other)
215 call ufo_locs_delete(temp_self)
216 
217 end subroutine ufo_locs_concatenate
218 
219 
220 ! --------------------------------------------------------------------------------------------------
221 
222 subroutine ufo_locs_init(self, obss, t1, t2)
223 
224  implicit none
225 
226  type(ufo_locs), intent(inout) :: self
227  type(c_ptr), value, intent(in) :: obss
228  type(datetime), intent(in) :: t1, t2
229 
230  integer :: nlocs
231 
232  character(len=*),parameter:: &
233  myname = "ufo_locs_init"
234  integer :: i
235  integer :: tw_nlocs
236  integer, dimension(:), allocatable :: tw_indx
237  real(kind_real), dimension(:), allocatable :: lon, lat
238  type(datetime), dimension(:), allocatable :: date_time
239 
240  ! Local copies pre binning
241  nlocs = obsspace_get_nlocs(obss)
242 
243  allocate(date_time(nlocs), lon(nlocs), lat(nlocs))
244 
245  call obsspace_get_db(obss, "MetaData", "datetime", date_time)
246 
247  ! Generate the timing window indices
248  allocate(tw_indx(nlocs))
249  tw_nlocs = 0
250  do i = 1, nlocs
251  if (date_time(i) > t1 .and. date_time(i) <= t2) then
252  tw_nlocs = tw_nlocs + 1
253  tw_indx(tw_nlocs) = i
254  endif
255  enddo
256 
257  call obsspace_get_db(obss, "MetaData", "longitude", lon)
258  call obsspace_get_db(obss, "MetaData", "latitude", lat)
259 
260  !Setup ufo locations
261  call ufo_locs_setup(self, tw_nlocs)
262  do i = 1, tw_nlocs
263  self%lon(i) = lon(tw_indx(i))
264  self%lat(i) = lat(tw_indx(i))
265  self%time(i) = date_time(tw_indx(i))
266  enddo
267  self%indx = tw_indx(1:tw_nlocs)
268 
269  do i = 1, nlocs
270  call datetime_delete(date_time(i))
271  enddo
272  deallocate(date_time, lon, lat, tw_indx)
273 
274  self%max_indx = obsspace_get_gnlocs(obss)
275 
276 
277 end subroutine ufo_locs_init
278 
279 ! --------------------------------------------------------------------------------------------------
280 
281 subroutine ufo_locs_time_mask(self, t1, t2, time_mask)
282 
283 type(ufo_locs), intent(in) :: self
284 type(datetime), intent(in) :: t1
285 type(datetime), intent(in) :: t2
286 logical, allocatable, intent(inout) :: time_mask(:)
287 
288 ! Locals
289 integer :: n
290 
291 ! Return a mask that is true where the location times are between t1 and t2
292 
293 ! Check for sensible inputs
294 if (t1>t2) call abor1_ftn("ufo_locs_mod.ufo_locs_time_mask t2 is not greater than or equal to t1")
295 
296 ! Allocate the array to output
297 if (.not.allocated(time_mask)) allocate(time_mask(self%nlocs))
298 time_mask = .false.
299 
300 ! Loop over times and check if between two times
301 do n = 1, self%nlocs
302 
303  ! Check if in the time range
304  if (self%time(n) > t1 .and. self%time(n) <= t2 ) then
305  time_mask(n) = .true.
306  endif
307 
308 enddo
309 
310 end subroutine ufo_locs_time_mask
311 
312 ! --------------------------------------------------------------------------------------------------
313 
314 end module ufo_locs_mod
ufo_locs_mod::ufo_locs_create
subroutine, public ufo_locs_create(self, obss, nlocs, lats, lons)
Definition: ufo_locs_mod.F90:39
ufo_locs_mod::ufo_locs_copy
subroutine, public ufo_locs_copy(self, other)
Definition: ufo_locs_mod.F90:106
ufo_locs_mod::ufo_locs_time_mask
subroutine, public ufo_locs_time_mask(self, t1, t2, time_mask)
Definition: ufo_locs_mod.F90:282
ufo_locs_mod::ufo_locs_delete
subroutine, public ufo_locs_delete(self)
Definition: ufo_locs_mod.F90:129
ufo_locs_mod::ufo_locs_concatenate
subroutine, public ufo_locs_concatenate(self, other)
Definition: ufo_locs_mod.F90:144
ufo_locs_mod::ufo_locs_setup
subroutine, public ufo_locs_setup(self, nlocs)
Definition: ufo_locs_mod.F90:81
ufo_locs_mod
Fortran module handling observation locations.
Definition: ufo_locs_mod.F90:9
ufo_locs_mod::ufo_locs
Fortran derived type to hold observation locations.
Definition: ufo_locs_mod.F90:25
ufo_locs_mod::ufo_locs_init
subroutine, public ufo_locs_init(self, obss, t1, t2)
Definition: ufo_locs_mod.F90:223