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 
13 use qg_fields_mod
15 
16 implicit none
17 
18 private
19 public :: qg_change_var_config
20 public :: qg_change_var_registry
22 ! ------------------------------------------------------------------------------
24  character(len=1024) :: varchange !< Variable change name
25 end type qg_change_var_config
26 
27 #define LISTED_TYPE qg_change_var_config
28 
29 !> Linked list interface - defines registry_t type
30 #include "oops/util/linkedList_i.f"
31 
32 !> Global registry
33 type(registry_t) :: qg_change_var_registry
34 ! ------------------------------------------------------------------------------
35 contains
36 ! ------------------------------------------------------------------------------
37 ! Public
38 ! ------------------------------------------------------------------------------
39 !> Linked list implementation
40 #include "oops/util/linkedList_c.f"
41 ! ------------------------------------------------------------------------------
42 !> Setup change of variable
43 subroutine qg_change_var_setup(self,vars_in,vars_out)
44 
45 implicit none
46 
47 ! Passed variables
48 type(qg_change_var_config),intent(inout) :: self !< Variable change
49 type(oops_variables),intent(in) :: vars_in !< Input variables
50 type(oops_variables),intent(in) :: vars_out !< Output variables
51 
52 ! Check
53 if ((vars_in%nvars() /= 1) .or. (vars_out%nvars() /= 1)) then
54  call abor1_ftn('qg_change_var_setup: wrong change of variable')
55 endif
56 if (vars_in%variable(1) == vars_out%variable(1)) then
57  self%varchange = 'identity'
58 elseif ((vars_in%variable(1) == 'x') .and. (vars_out%variable(1) == 'q')) then
59  self%varchange = 'x_to_q'
60 elseif ((vars_in%variable(1) == 'q') .and. (vars_out%variable(1) == 'x')) then
61  self%varchange = 'q_to_x'
62 else
63  call abor1_ftn('qg_change_var_setup: wrong change of variable')
64 endif
65 
66 end subroutine qg_change_var_setup
67 ! ------------------------------------------------------------------------------
68 !> Change of variable
69 subroutine qg_change_var(conf,fld_in,fld_out)
70 
71 implicit none
72 
73 ! Passed variables
74 type(qg_change_var_config),intent(in) :: conf !< Variable change
75 type(qg_fields),intent(in) :: fld_in !< Input fields
76 type(qg_fields),intent(inout) :: fld_out !< Output fields
77 
78 ! Check fields resolution
79 call qg_fields_check_resolution(fld_in,fld_out)
80 
81 select case (trim(conf%varchange))
82 case ('identity')
83  ! Copy fields
84  call qg_fields_copy(fld_out,fld_in)
85 case ('x_to_q')
86  ! Check fields variables
87  if (fld_in%lq) call abor1_ftn('qg_change_var: wrong input fields variables for '//trim(conf%varchange))
88  if (.not.fld_out%lq) call abor1_ftn('qg_change_var: wrong output fields variables for '//trim(conf%varchange))
89 
90  ! Conversion
91  call convert_x_to_q_tl(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
92 
93  ! Copy boundary conditions
94  call qg_fields_copy(fld_out,fld_in,.true.)
95 case ('q_to_x')
96  ! Check fields variables
97  if (.not.fld_in%lq) call abor1_ftn('qg_change_var: wrong input fields variables for '//trim(conf%varchange))
98  if (fld_out%lq) call abor1_ftn('qg_change_var: wrong output fields variables for '//trim(conf%varchange))
99 
100  ! Conversion
101  call convert_q_to_x_tl(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
102 
103  ! Copy boundary conditions
104  call qg_fields_copy(fld_out,fld_in,.true.)
105 case default
106  call abor1_ftn('qg_change_var: wrong variable change')
107 end select
108 
109 end subroutine qg_change_var
110 ! ------------------------------------------------------------------------------
111 !> Change of variable - inverse
112 subroutine qg_change_var_inv(conf,fld_in,fld_out)
113 
114 implicit none
115 
116 ! Passed variables
117 type(qg_change_var_config),intent(in) :: conf !< Variable change
118 type(qg_fields),intent(in) :: fld_in !< Input fields
119 type(qg_fields),intent(inout) :: fld_out !< Output fields
120 
121 ! Check fields resolution
122 call qg_fields_check_resolution(fld_in,fld_out)
123 
124 select case (trim(conf%varchange))
125 case ('identity')
126  ! Copy fields
127  call qg_fields_copy(fld_out,fld_in)
128 case ('x_to_q')
129  ! Check fields variables
130  if (.not.fld_in%lq) call abor1_ftn('qg_change_var_inv: wrong input fields variables for '//trim(conf%varchange))
131  if (fld_out%lq) call abor1_ftn('qg_change_var_inv: wrong output fields variables for '//trim(conf%varchange))
132 
133  ! Conversion
134  call convert_q_to_x_tl(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
135 
136  ! Copy boundary conditions
137  call qg_fields_copy(fld_out,fld_in,.true.)
138 case ('q_to_x')
139  ! Check fields variables
140  if (fld_in%lq) call abor1_ftn('qg_change_var_inv: wrong input fields variables for '//trim(conf%varchange))
141  if (.not.fld_out%lq) call abor1_ftn('qg_change_var_inv: wrong output fields variables for '//trim(conf%varchange))
142 
143  ! Conversion
144  call convert_x_to_q_tl(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
145 
146  ! Copy boundary conditions
147  call qg_fields_copy(fld_out,fld_in,.true.)
148 case default
149  call abor1_ftn('qg_change_var_inv: wrong variable change')
150 end select
151 
152 end subroutine qg_change_var_inv
153 ! ------------------------------------------------------------------------------
154 !> Change of variable - adjoint
155 subroutine qg_change_var_ad(conf,fld_in,fld_out)
156 
157 implicit none
158 
159 ! Passed variables
160 type(qg_change_var_config),intent(in) :: conf !< Variable change
161 type(qg_fields),intent(in) :: fld_in !< Input fields
162 type(qg_fields),intent(inout) :: fld_out !< Output fields
163 
164 ! Check fields resolution
165 call qg_fields_check_resolution(fld_in,fld_out)
166 
167 select case (trim(conf%varchange))
168 case ('identity')
169  ! Copy fields
170  call qg_fields_copy(fld_out,fld_in)
171 case ('x_to_q')
172  ! Check fields variables
173  if (.not.fld_in%lq) call abor1_ftn('qg_change_var_ad: wrong input fields variables for '//trim(conf%varchange))
174  if (fld_out%lq) call abor1_ftn('qg_change_var_ad: wrong output fields variables for '//trim(conf%varchange))
175 
176  ! Conversion
177  call convert_x_to_q_ad(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
178 
179  ! Copy boundary conditions
180  call qg_fields_copy(fld_out,fld_in,.true.)
181 case ('q_to_x')
182  ! Check fields variables
183  if (fld_in%lq) call abor1_ftn('qg_change_var_ad: wrong input fields variables for '//trim(conf%varchange))
184  if (.not.fld_out%lq) call abor1_ftn('qg_change_var_ad: wrong output fields variables for '//trim(conf%varchange))
185 
186  ! Conversion
187  call convert_q_to_x_ad(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
188 
189  ! Copy boundary conditions
190  call qg_fields_copy(fld_out,fld_in,.true.)
191 case default
192  call abor1_ftn('qg_change_var_ad: wrong variable change')
193 end select
194 
195 end subroutine qg_change_var_ad
196 ! ------------------------------------------------------------------------------
197 !> Change of variable - inverse adjoint
198 subroutine qg_change_var_inv_ad(conf,fld_in,fld_out)
199 
200 implicit none
201 
202 ! Passed variables
203 type(qg_change_var_config),intent(in) :: conf !< Variable change
204 type(qg_fields),intent(in) :: fld_in !< Input fields
205 type(qg_fields),intent(inout) :: fld_out !< Output fields
206 
207 ! Check fields resolution
208 call qg_fields_check_resolution(fld_in,fld_out)
209 
210 select case (trim(conf%varchange))
211 case ('identity')
212  ! Copy fields
213  call qg_fields_copy(fld_out,fld_in)
214 case ('x_to_q')
215  ! Check fields variables
216  if (fld_in%lq) call abor1_ftn('qg_change_var_inv_ad: wrong input fields variables for '//trim(conf%varchange))
217  if (.not.fld_out%lq) call abor1_ftn('qg_change_var_inv_ad: wrong output fields variables for '//trim(conf%varchange))
218 
219  ! Conversion
220  call convert_q_to_x_ad(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
221 
222  ! Copy boundary conditions
223  call qg_fields_copy(fld_out,fld_in,.true.)
224 case ('q_to_x')
225  ! Check fields variables
226  if (.not.fld_in%lq) call abor1_ftn('qg_change_var_inv_ad: wrong input fields variables for '//trim(conf%varchange))
227  if (fld_out%lq) call abor1_ftn('qg_change_var_inv_ad: wrong output fields variables for '//trim(conf%varchange))
228 
229  ! Conversion
230  call convert_x_to_q_ad(fld_in%geom,fld_in%gfld3d,fld_out%gfld3d)
231 
232  ! Copy boundary conditions
233  call qg_fields_copy(fld_out,fld_in,.true.)
234 case default
235  call abor1_ftn('qg_change_var_inv_ad: wrong variable change')
236 end select
237 
238 end subroutine qg_change_var_inv_ad
239 ! ------------------------------------------------------------------------------
240 end module qg_change_var_mod
qg_fields_mod
Definition: qg_fields_mod.F90:9
qg_fields_mod::qg_fields
Definition: qg_fields_mod.F90:51
qg_convert_x_to_q_mod
Definition: qg_convert_x_to_q_mod.F90:9
qg_convert_q_to_x_mod::convert_q_to_x_tl
subroutine, public convert_q_to_x_tl(geom, q, x)
Convert potential vorticity to streamfunction - tangent linear.
Definition: qg_convert_q_to_x_mod.F90:86
qg_change_var_mod::qg_change_var_ad
subroutine, public qg_change_var_ad(conf, fld_in, fld_out)
Change of variable - adjoint.
Definition: qg_change_var_mod.F90:156
qg_change_var_mod
Definition: qg_change_var_mod.F90:9
qg_change_var_mod::qg_change_var
subroutine, public qg_change_var(conf, fld_in, fld_out)
Change of variable.
Definition: qg_change_var_mod.F90:70
qg_fields_mod::qg_fields_check_resolution
subroutine, public qg_fields_check_resolution(fld1, fld2)
Check fields resolution.
Definition: qg_fields_mod.F90:1490
qg_convert_x_to_q_mod::convert_x_to_q_tl
subroutine, public convert_x_to_q_tl(geom, x, q)
Convert streamfunction to potential vorticity - tangent linear.
Definition: qg_convert_x_to_q_mod.F90:73
oops_variables_mod
Fortran interface to Variables.
Definition: variables_mod.F90:9
qg_convert_q_to_x_mod
Definition: qg_convert_q_to_x_mod.F90:9
qg_change_var_mod::qg_change_var_inv
subroutine, public qg_change_var_inv(conf, fld_in, fld_out)
Change of variable - inverse.
Definition: qg_change_var_mod.F90:113
qg_change_var_mod::qg_change_var_config
Definition: qg_change_var_mod.F90:23
qg_change_var_mod::qg_change_var_setup
subroutine, public qg_change_var_setup(self, vars_in, vars_out)
Linked list implementation.
Definition: qg_change_var_mod.F90:44
qg_change_var_mod::qg_change_var_registry
type(registry_t), public qg_change_var_registry
Linked list interface - defines registry_t type.
Definition: qg_change_var_mod.F90:33
qg_change_var_mod::qg_change_var_inv_ad
subroutine, public qg_change_var_inv_ad(conf, fld_in, fld_out)
Change of variable - inverse adjoint.
Definition: qg_change_var_mod.F90:199
qg_convert_x_to_q_mod::convert_x_to_q_ad
subroutine, public convert_x_to_q_ad(geom, q, x)
Convert streamfunction to potential vorticity - adjoint.
Definition: qg_convert_x_to_q_mod.F90:107
oops_variables_mod::oops_variables
Definition: variables_mod.F90:16
qg_fields_mod::qg_fields_copy
subroutine, public qg_fields_copy(self, other, bconly)
Copy fields.
Definition: qg_fields_mod.F90:321
qg_convert_q_to_x_mod::convert_q_to_x_ad
subroutine, public convert_q_to_x_ad(geom, x, q)
Convert potential vorticity to streamfunction - adjoint.
Definition: qg_convert_q_to_x_mod.F90:133