OOPS
qg_change_var_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
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 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
10 
11 use kinds
16 use qg_fields_mod
18 
19 implicit none
20 
21 private
22 public :: qg_change_var_registry
23 public :: qg_change_var
25 ! ------------------------------------------------------------------------------
27  ! Input variables flags
28  logical :: x_in
29  logical :: q_in
30  logical :: u_in
31  logical :: v_in
32 
33  ! Output variables flags
34  logical :: x_out
35  logical :: q_out
36  logical :: u_out
37  logical :: v_out
38 
39  ! Conversions
40  logical :: x_to_q
41  logical :: q_to_x
42  logical :: x_to_u
43  logical :: x_to_v
44 end type qg_change_var_config
45 
46 #define LISTED_TYPE qg_change_var_config
47 
48 !> Linked list interface - defines registry_t type
49 #include "oops/util/linkedList_i.f"
50 
51 !> Global registry
52 type(registry_t) :: qg_change_var_registry
53 ! ------------------------------------------------------------------------------
54 contains
55 ! ------------------------------------------------------------------------------
56 ! Public
57 ! ------------------------------------------------------------------------------
58 !> Linked list implementation
59 #include "oops/util/linkedList_c.f"
60 ! ------------------------------------------------------------------------------
61 !> Setup change of variables setup
62 subroutine qg_change_var_setup(conf,fld_in,fld_out,ad)
63 
64 implicit none
65 
66 ! Passed variables
67 type(qg_change_var_config),intent(inout) :: conf !< Variable change
68 type(qg_fields),intent(in) :: fld_in !< Input field
69 type(qg_fields),intent(in) :: fld_out !< Output field
70 logical,intent(in),optional :: ad !< Adjoint flag
71 
72 ! Local variables
73 logical :: lad
74 
75 ! Local flag
76 lad = .false.
77 if (present(ad)) lad = ad
78 
79 ! Check what is allocated in fields
80 if (lad) then
81  ! Adjoint case
82  conf%x_in = allocated(fld_out%x)
83  conf%q_in = allocated(fld_out%q)
84  conf%u_in = allocated(fld_out%u)
85  conf%v_in = allocated(fld_out%v)
86  conf%x_out = allocated(fld_in%x)
87  conf%q_out = allocated(fld_in%q)
88  conf%u_out = allocated(fld_in%u)
89  conf%v_out = allocated(fld_in%v)
90 else
91  ! Normal case
92  conf%x_in = allocated(fld_in%x)
93  conf%q_in = allocated(fld_in%q)
94  conf%u_in = allocated(fld_in%u)
95  conf%v_in = allocated(fld_in%v)
96  conf%x_out = allocated(fld_out%x)
97  conf%q_out = allocated(fld_out%q)
98  conf%u_out = allocated(fld_out%u)
99  conf%v_out = allocated(fld_out%v)
100 endif
101 
102 ! Check input/output consistency
103 if (conf%x_out.and.(.not.(conf%x_in.or.conf%q_in))) call abor1_ftn('qg_change_var_setup: x or q required to compute x')
104 if (conf%q_out.and.(.not.(conf%x_in.or.conf%q_in))) call abor1_ftn('qg_change_var_setup: x or q required to compute q')
105 if (conf%u_out.and.(.not.(conf%x_in.or.conf%q_in.or.conf%u_in))) &
106  & call abor1_ftn('qg_change_var_setup: x, q or u required to compute u')
107 if (conf%v_out.and.(.not.(conf%x_in.or.conf%q_in.or.conf%v_in))) &
108  & call abor1_ftn('qg_change_var_setup: x, q or v required to compute v')
109 
110 ! Initialize required conversions
111 conf%q_to_x = .false.
112 conf%x_to_q = .false.
113 conf%x_to_u = .false.
114 conf%x_to_v = .false.
115 
116 ! Define required conversions
117 if (conf%x_out) conf%q_to_x = (.not.conf%x_in)
118 if (conf%q_out) conf%x_to_q = (.not.conf%q_in)
119 if (conf%u_out) then
120  if (.not.conf%u_in) then
121  conf%x_to_u = .true.
122  conf%q_to_x = (.not.conf%x_in)
123  endif
124 endif
125 if (conf%v_out) then
126  if (.not.conf%v_in) then
127  conf%x_to_v = .true.
128  conf%q_to_x = (.not.conf%x_in)
129  endif
130 endif
131 
132 end subroutine qg_change_var_setup
133 ! ------------------------------------------------------------------------------
134 !> Get variables
135 subroutine qg_change_var_get(conf,fld,x,q,u,v)
136 
137 implicit none
138 
139 ! Passed variables
140 type(qg_change_var_config),intent(in) :: conf !< Variable change
141 type(qg_fields),intent(in) :: fld !< Fields
142 real(kind_real),intent(out) :: x(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Streamfunction
143 real(kind_real),intent(out) :: q(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Potential vorticity
144 real(kind_real),intent(out) :: u(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Zonal wind
145 real(kind_real),intent(out) :: v(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Meridional wind
146 
147 if (allocated(fld%x)) then
148  x = fld%x
149 else
150  x = 0.0_kind_real
151 endif
152 if (allocated(fld%q)) then
153  q = fld%q
154 else
155  q = 0.0_kind_real
156 endif
157 if (allocated(fld%u)) then
158  u = fld%u
159 else
160  u = 0.0_kind_real
161 endif
162 if (allocated(fld%v)) then
163  v = fld%v
164 else
165  v = 0.0_kind_real
166 endif
167 
168 end subroutine qg_change_var_get
169 ! ------------------------------------------------------------------------------
170 !> Set variables
171 subroutine qg_change_var_set(conf,fld,x,q,u,v)
172 
173 implicit none
174 
175 ! Passed variables
176 type(qg_change_var_config),intent(in) :: conf !< Variable change
177 type(qg_fields),intent(inout) :: fld !< Fields
178 real(kind_real),intent(in) :: x(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Streamfunction
179 real(kind_real),intent(in) :: q(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Potential vorticity
180 real(kind_real),intent(in) :: u(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Zonal wind
181 real(kind_real),intent(in) :: v(fld%geom%nx,fld%geom%ny,fld%geom%nz) !< Meridional wind
182 
183 if (allocated(fld%x)) fld%x = x
184 if (allocated(fld%q)) fld%q = q
185 if (allocated(fld%u)) fld%u = u
186 if (allocated(fld%v)) fld%v = v
187 
188 end subroutine qg_change_var_set
189 ! ------------------------------------------------------------------------------
190 !> Change of variable
191 subroutine qg_change_var(fld_in,fld_out)
192 
193 implicit none
194 
195 ! Passed variables
196 type(qg_fields),intent(in) :: fld_in !< Input fields
197 type(qg_fields),intent(inout) :: fld_out !< Output fields
198 
199 ! Local variables
200 real(kind_real) :: x(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),q(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
201 real(kind_real) :: u(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),v(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
202 type(qg_change_var_config) :: conf
203 
204 ! Check resolution
205 call qg_fields_check_resolution(fld_in,fld_out)
206 
207 ! Copy boundary conditions
208 call qg_fields_copy_lbc(fld_out,fld_in)
209 
210 ! Define change of variable configuration
211 call qg_change_var_setup(conf,fld_in,fld_out)
212 
213 ! Get x, q, u and v
214 call qg_change_var_get(conf,fld_in,x,q,u,v)
215 
216 ! Conversions
217 if (conf%x_to_q) call convert_x_to_q(fld_in%geom,x,fld_in%x_north,fld_in%x_south,q)
218 if (conf%q_to_x) call convert_q_to_x(fld_in%geom,q,fld_in%x_north,fld_in%x_south,x)
219 if (conf%x_to_u) call convert_x_to_u(fld_in%geom,x,fld_in%x_north,fld_in%x_south,u)
220 if (conf%x_to_v) call convert_x_to_v(fld_in%geom,x,v)
221 
222 ! Set x, q, u and v
223 call qg_change_var_set(conf,fld_out,x,q,u,v)
224 
225 end subroutine qg_change_var
226 ! ------------------------------------------------------------------------------
227 !> Change of variable
228 subroutine qg_change_var_tl(fld_in,fld_out)
229 
230 implicit none
231 
232 ! Passed variables
233 type(qg_fields),intent(in) :: fld_in !< Input fields
234 type(qg_fields),intent(inout) :: fld_out !< Output fields
235 
236 ! Local variables
237 real(kind_real) :: x(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),q(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
238 real(kind_real) :: u(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),v(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
239 type(qg_change_var_config) :: conf
240 
241 ! Check resolution
242 call qg_fields_check_resolution(fld_in,fld_out)
243 
244 ! Copy boundary conditions
245 call qg_fields_copy_lbc(fld_out,fld_in)
246 
247 ! Define change of variable configuration
248 call qg_change_var_setup(conf,fld_in,fld_out)
249 
250 ! Get x, q, u and v
251 call qg_change_var_get(conf,fld_in,x,q,u,v)
252 
253 ! Conversions
254 if (conf%x_to_q) call convert_x_to_q_tl(fld_in%geom,x,q)
255 if (conf%q_to_x) call convert_q_to_x_tl(fld_in%geom,q,x)
256 if (conf%x_to_u) call convert_x_to_u_tl(fld_in%geom,x,u)
257 if (conf%x_to_v) call convert_x_to_v_tl(fld_in%geom,x,v)
258 
259 ! Set x, q, u and v
260 call qg_change_var_set(conf,fld_out,x,q,u,v)
261 
262 end subroutine qg_change_var_tl
263 ! ------------------------------------------------------------------------------
264 !> Change of variable - adjoint
265 subroutine qg_change_var_ad(fld_in,fld_out)
266 
267 implicit none
268 
269 ! Passed variables
270 type(qg_fields),intent(in) :: fld_in !< Input fields
271 type(qg_fields),intent(inout) :: fld_out !< Output fields
272 
273 ! Local variables
274 real(kind_real) :: x(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),q(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
275 real(kind_real) :: u(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz),v(fld_in%geom%nx,fld_in%geom%ny,fld_in%geom%nz)
276 type(qg_change_var_config) :: conf
277 
278 ! Checks
279 call qg_fields_check_resolution(fld_in,fld_out)
280 
281 ! Copy boundary conditions
282 call qg_fields_copy_lbc(fld_out,fld_in)
283 
284 ! Define change of variable configuration
285 call qg_change_var_setup(conf,fld_in,fld_out,.true.)
286 
287 ! Get x, q, u and v
288 call qg_change_var_get(conf,fld_in,x,q,u,v)
289 
290 ! Conversions
291 if (conf%x_to_v) call convert_x_to_v_ad(fld_in%geom,v,x)
292 if (conf%x_to_u) call convert_x_to_u_ad(fld_in%geom,u,x)
293 if (conf%q_to_x) call convert_q_to_x_ad(fld_in%geom,x,q)
294 if (conf%x_to_q) call convert_x_to_q_ad(fld_in%geom,q,x)
295 
296 ! Set x, q, u and v
297 call qg_change_var_set(conf,fld_out,x,q,u,v)
298 
299 end subroutine qg_change_var_ad
300 ! ------------------------------------------------------------------------------
301 end module qg_change_var_mod
Fortran interface to Variables.
subroutine qg_change_var_set(conf, fld, x, q, u, v)
Set variables.
subroutine qg_change_var_setup(conf, fld_in, fld_out, ad)
Linked list implementation.
type(registry_t), public qg_change_var_registry
Linked list interface - defines registry_t type.
subroutine, public qg_change_var(fld_in, fld_out)
Change of variable.
subroutine qg_change_var_get(conf, fld, x, q, u, v)
Get variables.
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 convert_q_to_x_tl(geom, q, x)
Convert potential vorticity to streamfunction - tangent linear.
subroutine, public convert_q_to_x(geom, q, x_north, x_south, x)
Convert potential vorticity to streamfunction.
subroutine, public convert_q_to_x_ad(geom, x, q)
Convert potential vorticity to streamfunction - adjoint.
subroutine, public convert_x_to_q_tl(geom, x, q)
Convert streamfunction to potential vorticity - tangent linear.
subroutine, public convert_x_to_q(geom, x, x_north, x_south, q)
Convert streamfunction to potential vorticity.
subroutine, public convert_x_to_q_ad(geom, q, x)
Convert streamfunction to potential vorticity - adjoint.
subroutine, public convert_x_to_u(geom, x, x_north, x_south, u)
Convert streafunction to zonal wind.
subroutine, public convert_x_to_u_tl(geom, x, u)
Convert streafunction to zonal wind - tangent Linear.
subroutine, public convert_x_to_u_ad(geom, u, x)
Convert streafunction to zonal wind - adjoint.
subroutine, public convert_x_to_v(geom, x, v)
Convert streafunction to meridional wind.
subroutine, public convert_x_to_v_tl(geom, x, v)
Convert streafunction to meridional wind - tangent Linear.
subroutine, public convert_x_to_v_ad(geom, v, x)
Convert streafunction to meridional wind - adjoint.
subroutine, public qg_fields_copy_lbc(self, other)
Copy fields LBC.
subroutine, public qg_fields_check_resolution(fld1, fld2)
Check fields resolution.