FV3 Bundle
ioda_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 kinds
13 
14 implicit none
15 private
17 
18 ! ------------------------------------------------------------------------------
19 
20 !> Fortran derived type to hold observation locations
21 type :: ioda_locs
22  integer :: nlocs
23  real(kind_real), allocatable, dimension(:) :: lat !< latitude
24  real(kind_real), allocatable, dimension(:) :: lon !< longitude
25  real(kind_real), allocatable, dimension(:) :: time !< obs-time
26  integer, allocatable, dimension(:) :: indx !< indices of locations in the full [geovals] array
27 end type ioda_locs
28 
29 ! ------------------------------------------------------------------------------
30 contains
31 ! ------------------------------------------------------------------------------
32 
33 subroutine ioda_locs_create(self, nlocs, lats, lons, rdist)
34 implicit none
35 type(ioda_locs), intent(inout) :: self
36 integer, intent(in) :: nlocs
37 real(kind_real), intent(in) :: lats(nlocs)
38 real(kind_real), intent(in) :: lons(nlocs)
39 integer, intent(in) :: rdist
40 
41 type(random_distribution) :: ran_dist
42 integer, allocatable :: dist_indx(:)
43 real(kind_real), allocatable, dimension(:) :: latsr, lonsr, timer
44 integer :: n
45 
46 self%nlocs = nlocs
47 allocate(self%lat(nlocs), self%lon(nlocs), self%time(nlocs))
48 allocate(self%indx(nlocs))
49 self%lat(:) = lats(:)
50 self%lon(:) = lons(:)
51 self%time(:) = 0.0
52 do n = 1, self%nlocs
53  self%indx(n) = n
54 enddo
55 
56 ran_dist = random_distribution(self%nlocs)
57 dist_indx = ran_dist%indx
58 
59 if (rdist == 1) then
60 
61  !Redistribute randomly
62  self%nlocs = ran_dist%nobs_pe()
63  allocate(latsr(self%nlocs), lonsr(self%nlocs), timer(self%nlocs))
64  do n = 1,self%nlocs
65  latsr(n) = self%lat(dist_indx(n))
66  lonsr(n) = self%lon(dist_indx(n))
67  timer(n) = self%time(dist_indx(n))
68  enddo
69  deallocate(self%lat, self%lon, self%time, self%indx)
70 
71  allocate(self%lat(self%nlocs), self%lon(self%nlocs), self%time(self%nlocs), self%indx(self%nlocs))
72  self%lat(:) = latsr(:)
73  self%lon(:) = lonsr(:)
74  self%time(:) = 0.0
75 
76  do n = 1, self%nlocs
77  self%indx(n) = n
78  enddo
79 
80  deallocate(latsr, lonsr, timer)
81 
82 endif
83 
84 end subroutine ioda_locs_create
85 
86 ! ------------------------------------------------------------------------------
87 
88 subroutine ioda_locs_setup(self, nlocs)
89 implicit none
90 type(ioda_locs), intent(inout) :: self
91 integer, intent(in) :: nlocs
92 
93 call ioda_locs_delete(self)
94 
95 self%nlocs = nlocs
96 allocate(self%lat(nlocs), self%lon(nlocs), self%time(nlocs), self%indx(nlocs))
97 self%lat(:) = 0.0
98 self%lon(:) = 0.0
99 self%time(:) = 0.0
100 self%indx(:) = 0
101 
102 end subroutine ioda_locs_setup
103 
104 ! ------------------------------------------------------------------------------
105 
106 subroutine ioda_locs_delete(self)
107 implicit none
108 type(ioda_locs), intent(inout) :: self
109 
110 self%nlocs = 0
111 if (allocated(self%lat)) deallocate(self%lat)
112 if (allocated(self%lon)) deallocate(self%lon)
113 if (allocated(self%time)) deallocate(self%time)
114 if (allocated(self%indx)) deallocate(self%indx)
115 
116 end subroutine ioda_locs_delete
117 
118 ! ------------------------------------------------------------------------------
119 
120 end module ioda_locs_mod
Fortran derived type to hold observation locations.
subroutine, public ioda_locs_create(self, nlocs, lats, lons, rdist)
subroutine, public ioda_locs_delete(self)
Fortran module handling observation locations.
subroutine, public ioda_locs_setup(self, nlocs)