UFO
locations_test.F90
Go to the documentation of this file.
1 !
2 ! (C) Copyright 2020 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 !
8 
9 use iso_c_binding
10 
11 implicit none
12 private
13 
14 contains
15 
16 ! ------------------------------------------------------------------------------
17 !> Tests whether latitudes and longitudes from ufo::Locations object
18 !! are matching references in yaml file
19 integer function test_locations_lonslats_c(c_locs, c_conf) bind(c,name='test_locations_lonslats_f90')
20 use fckit_configuration_module, only: fckit_configuration
21 use fckit_log_module, only: fckit_log
23 implicit none
24 type(c_ptr), value, intent(in) :: c_locs !< ufo::Locations object
25 type(c_ptr), value, intent(in) :: c_conf !< test configuration (described in LocationsTestParameters
26 
27 !> local variables
28 type(ufo_locations) :: locs
29 type(fckit_configuration) :: f_conf
30 integer :: nlocs, iloc
31 real :: tolerance = 1.e-8
32 real(c_double), dimension(:), allocatable :: lons_ref, lats_ref, lons, lats
33 character(len=200) :: logmessage
34 
35 !> default value: test passed
37 
38 !> read reference longitudes and latitudes
39 f_conf = fckit_configuration(c_conf)
40 nlocs = f_conf%get_size("reference lons")
41 call f_conf%get_or_die("reference lons", lons_ref)
42 call f_conf%get_or_die("reference lats", lats_ref)
43 
44 !> get longitudes and latitudes from locations object
45 locs = ufo_locations(c_locs)
46 allocate(lons(nlocs), lats(nlocs))
47 call locs%get_lons(lons)
48 call locs%get_lats(lats)
49 
50 do iloc = 1, nlocs
51  write(logmessage, *) "lons: loc and ref: ", lons(iloc), ", ", lons_ref(iloc)
52  call fckit_log%debug(logmessage)
53  write(logmessage, *) "lats: loc and ref: ", lats(iloc), ", ", lats_ref(iloc)
54  call fckit_log%debug(logmessage)
55 enddo
56 
57 !> compare to reference
58 if(any(abs(lons-lons_ref) > tolerance)) test_locations_lonslats_c = 0
59 if(any(abs(lats-lats_ref) > tolerance)) test_locations_lonslats_c = 0
60 
61 deallocate(lons, lats)
62 
63 end function test_locations_lonslats_c
64 
65 ! ------------------------------------------------------------------------------
66 
67 ! ------------------------------------------------------------------------------
68 !> Tests whether ufo::Locations time mask method results
69 !! are matching references in yaml file
70 integer function test_locations_timemask_c(c_locs, c_conf) bind(c,name='test_locations_timemask_f90')
71 use fckit_configuration_module, only: fckit_configuration
72 use fckit_log_module, only: fckit_log
73 use datetime_mod
75 implicit none
76 type(c_ptr), value, intent(in) :: c_locs !< ufo::Locations object
77 type(c_ptr), value, intent(in) :: c_conf !< test configuration (described in LocationsTestParameters
78 
79 !> local variables
80 type(ufo_locations) :: locs
81 type(fckit_configuration) :: f_conf
82 type(fckit_configuration), allocatable :: testconfigs(:)
83 character(kind=c_char,len=:), allocatable :: t1str, t2str
84 logical(kind=c_bool), allocatable :: mask(:)
85 logical, allocatable :: mask_ref(:)
86 type(datetime) :: t1, t2
87 integer :: nlocs, iloc, itest
88 character(len=200) :: logmessage
89 
90 !> default value: test passed
92 
93 locs = ufo_locations(c_locs)
94 f_conf = fckit_configuration(c_conf)
95 call f_conf%get_or_die("time mask tests", testconfigs)
96 !> loop over all the tests for time masks
97 do itest = 1, size(testconfigs)
98  call testconfigs(itest)%get_or_die("t1", t1str)
99  call testconfigs(itest)%get_or_die("t2", t2str)
100  call datetime_create(t1str, t1)
101  call datetime_create(t2str, t2)
102  call testconfigs(itest)%get_or_die("reference mask", mask_ref)
103 
104  !> call ufo::Locations method to compute time mask
105  nlocs = locs%nlocs()
106  allocate(mask(nlocs))
107  call locs%get_timemask(t1, t2, mask)
108 
109  write(logmessage, *) "Test ", itest
110  call fckit_log%debug(logmessage)
111  do iloc = 1, nlocs
112  write(logmessage, *) "mask: loc and ref: ", mask(iloc), ", ", mask_ref(iloc)
113  call fckit_log%debug(logmessage)
114  enddo
115 
116  !> compare to reference
117  if(any(mask .neqv. mask_ref)) test_locations_timemask_c = 0
118 
119  deallocate(mask)
120 
121 enddo
122 
123 end function test_locations_timemask_c
124 
125 end module test_locations
integer function test_locations_timemask_c(c_locs, c_conf)
Tests whether ufo::Locations time mask method results are matching references in yaml file.
integer function test_locations_lonslats_c(c_locs, c_conf)
Tests whether latitudes and longitudes from ufo::Locations object are matching references in yaml fil...
Fortran interface to ufo::Locations.
integer function nlocs(this)
Return the number of observational locations in this Locations object.