OOPS
qg_getvalues_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2020 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 
7 
8 use atlas_module, only: atlas_field
9 use datetime_mod
10 use fckit_log_module,only: fckit_log
11 use iso_c_binding
12 use kinds
13 !$ use omp_lib
16 use qg_fields_mod
17 use qg_geom_mod
18 use qg_gom_mod
19 use qg_interp_mod
20 use qg_locs_mod
21 
22 implicit none
23 
24 private
26 
27 ! ------------------------------------------------------------------------------
28 contains
29 ! ------------------------------------------------------------------------------
30 ! Public
31 ! ------------------------------------------------------------------------------
32 !> Interpolation from fields
33 subroutine qg_getvalues_interp(locs,fld,t1,t2,gom)
34 
35 implicit none
36 
37 ! Passed variables
38 type(qg_locs), intent(in) :: locs !< Locations
39 type(qg_fields),intent(in) :: fld !< Fields
40 type(datetime),intent(in) :: t1, t2 !< times
41 type(qg_gom),intent(inout) :: gom !< Interpolated values
42 
43 ! Local variables
44 integer :: jloc
45 real(kind_real), pointer :: lonlat(:,:), z(:)
46 type(atlas_field) :: lonlat_field, z_field
47 type(qg_fields) :: fld_gom
48 
49 ! Get locations
50 lonlat_field = locs%lonlat()
51 call lonlat_field%data(lonlat)
52 z_field = locs%altitude()
53 call z_field%data(z)
54 
55 ! Check field
56 call qg_fields_check(fld)
57 
58 ! Create field with GOM variables
59 call qg_fields_create(fld_gom,fld%geom,gom%vars,.true.)
60 call qg_fields_copy_lbc(fld_gom,fld)
61 
62 ! Apply change of variables
63 call qg_change_var(fld,fld_gom)
64 
65 !$omp parallel do schedule(static) private(jloc)
66 do jloc=1,locs%nlocs()
67  ! Check if current obs is in this time frame (t1,t2]
68  if (t1 < locs%times(jloc) .and. locs%times(jloc) <= t2) then
69  ! Interpolate variables
70  if (gom%vars%has('x')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
71  & z(jloc),fld_gom%x,gom%x(jloc))
72  if (gom%vars%has('q')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
73  & z(jloc),fld_gom%q,gom%q(jloc))
74  if (gom%vars%has('u')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
75  & z(jloc),fld_gom%u,gom%u(jloc))
76  if (gom%vars%has('v')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
77  & z(jloc),fld_gom%v,gom%v(jloc))
78  endif
79 enddo
80 !$omp end parallel do
81 
82 ! Release memory
83 call lonlat_field%final()
84 call z_field%final()
85 
86 end subroutine qg_getvalues_interp
87 ! ------------------------------------------------------------------------------
88 !> Interpolation from fields - tangent linear
89 subroutine qg_getvalues_interp_tl(locs,fld,t1,t2,gom)
90 
91 implicit none
92 
93 ! Passed variables
94 type(qg_locs), intent(in) :: locs !< Locations
95 type(qg_fields),intent(in) :: fld !< Fields
96 type(datetime),intent(in) :: t1, t2 !< times
97 type(qg_gom),intent(inout) :: gom !< Interpolated values
98 
99 ! Local variables
100 integer :: jloc
101 real(kind_real), pointer :: lonlat(:,:), z(:)
102 type(atlas_field) :: lonlat_field, z_field
103 type(qg_fields) :: fld_gom
104 
105 ! Get locations
106 lonlat_field = locs%lonlat()
107 call lonlat_field%data(lonlat)
108 z_field = locs%altitude()
109 call z_field%data(z)
110 
111 ! Check field
112 call qg_fields_check(fld)
113 
114 ! Create field with GOM variables
115 call qg_fields_create(fld_gom,fld%geom,gom%vars,.false.)
116 
117 ! Apply change of variables
118 call qg_change_var_tl(fld,fld_gom)
119 
120 !$omp parallel do schedule(static) private(jloc)
121 do jloc=1,locs%nlocs()
122  ! Check if current obs is in this time frame (t1,t2]
123  if (t1 < locs%times(jloc) .and. locs%times(jloc) <= t2) then
124  ! Interpolate variables
125  if (gom%vars%has('x')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
126  & z(jloc),fld_gom%x,gom%x(jloc))
127  if (gom%vars%has('q')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
128  & z(jloc),fld_gom%q,gom%q(jloc))
129  if (gom%vars%has('u')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
130  & z(jloc),fld_gom%u,gom%u(jloc))
131  if (gom%vars%has('v')) call qg_interp_trilinear(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
132  & z(jloc),fld_gom%v,gom%v(jloc))
133  endif
134 enddo
135 !$omp end parallel do
136 
137 ! Release memory
138 call lonlat_field%final()
139 call z_field%final()
140 
141 end subroutine qg_getvalues_interp_tl
142 ! ------------------------------------------------------------------------------
143 !> Interpolation from fields - adjoint
144 subroutine qg_getvalues_interp_ad(locs,fld,t1,t2,gom)
145 
146 implicit none
147 
148 ! Passed variables
149 type(qg_locs), intent(in) :: locs !< Locations
150 type(qg_fields),intent(inout) :: fld !< Fields
151 type(datetime),intent(in) :: t1, t2 !< times
152 type(qg_gom),intent(in) :: gom !< Interpolated values
153 
154 ! Local variables
155 integer :: jloc
156 real(kind_real),allocatable :: x(:,:,:),q(:,:,:),u(:,:,:),v(:,:,:)
157 real(kind_real), pointer :: lonlat(:,:), z(:)
158 type(atlas_field) :: lonlat_field, z_field
159 type(qg_fields) :: fld_gom,fld_tmp
160 
161 ! Get locations
162 lonlat_field = locs%lonlat()
163 call lonlat_field%data(lonlat)
164 z_field = locs%altitude()
165 call z_field%data(z)
166 
167 ! Check field
168 call qg_fields_check(fld)
169 
170 ! Create field with GOM variables
171 call qg_fields_create(fld_gom,fld%geom,gom%vars,.false.)
172 
173 ! Initialization
174 call qg_fields_zero(fld_gom)
175 
176 do jloc=locs%nlocs(),1,-1
177  ! Check if current obs is in this time frame (t1,t2]
178  if (t1 < locs%times(jloc) .and. locs%times(jloc) <= t2) then
179  ! Interpolate variables
180  if (gom%vars%has('x')) call qg_interp_trilinear_ad(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
181  & z(jloc),gom%x(jloc),fld_gom%x)
182  if (gom%vars%has('q')) call qg_interp_trilinear_ad(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
183  & z(jloc),gom%q(jloc),fld_gom%q)
184  if (gom%vars%has('u')) call qg_interp_trilinear_ad(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
185  & z(jloc),gom%u(jloc),fld_gom%u)
186  if (gom%vars%has('v')) call qg_interp_trilinear_ad(fld%geom,lonlat(1,jloc),lonlat(2,jloc), &
187  & z(jloc),gom%v(jloc),fld_gom%v)
188  endif
189 enddo
190 
191 ! Apply change of variables
192 call qg_fields_create_from_other(fld_tmp,fld,fld%geom)
193 call qg_change_var_ad(fld_gom,fld_tmp)
194 call qg_fields_self_add(fld,fld_tmp)
195 
196 ! Release memory
197 call lonlat_field%final()
198 call z_field%final()
199 
200 end subroutine qg_getvalues_interp_ad
201 ! ------------------------------------------------------------------------------
202 end module qg_getvalues_mod
Fortran interface to Variables.
subroutine, public qg_change_var(fld_in, fld_out)
Change of variable.
subroutine, public qg_change_var_tl(fld_in, fld_out)
Change of variable.
subroutine, public qg_change_var_ad(fld_in, fld_out)
Change of variable - adjoint.
subroutine, public qg_fields_copy_lbc(self, other)
Copy fields LBC.
subroutine, public qg_fields_zero(self)
Set fields to zero.
subroutine, public qg_fields_create(self, geom, vars, lbc)
Linked list implementation.
subroutine, public qg_fields_check(self)
Check fields.
subroutine, public qg_fields_self_add(self, rhs)
Add fields.
subroutine, public qg_fields_create_from_other(self, other, geom)
Create fields from another one.
subroutine, public qg_getvalues_interp_tl(locs, fld, t1, t2, gom)
Interpolation from fields - tangent linear.
subroutine, public qg_getvalues_interp(locs, fld, t1, t2, gom)
Interpolation from fields.
subroutine, public qg_getvalues_interp_ad(locs, fld, t1, t2, gom)
Interpolation from fields - adjoint.
subroutine, public qg_interp_trilinear_ad(geom, lon, lat, z, val, field)
Interpolation - adjoint.
subroutine, public qg_interp_trilinear(geom, lon, lat, z, field, val)
Trilinear interpolation.