Go to the documentation of this file.
24 character(len=1024) :: varchange
27 #define LISTED_TYPE qg_change_var_config
30 #include "oops/util/linkedList_i.f"
40 #include "oops/util/linkedList_c.f"
53 if ((vars_in%nvars() /= 1) .or. (vars_out%nvars() /= 1))
then
54 call abor1_ftn(
'qg_change_var_setup: wrong change of variable')
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'
63 call abor1_ftn(
'qg_change_var_setup: wrong change of variable')
81 select case (trim(conf%varchange))
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))
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))
106 call abor1_ftn(
'qg_change_var: wrong variable change')
124 select case (trim(conf%varchange))
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))
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))
149 call abor1_ftn(
'qg_change_var_inv: wrong variable change')
167 select case (trim(conf%varchange))
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))
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))
192 call abor1_ftn(
'qg_change_var_ad: wrong variable change')
210 select case (trim(conf%varchange))
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))
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))
235 call abor1_ftn(
'qg_change_var_inv_ad: wrong variable change')
subroutine, public convert_q_to_x_tl(geom, q, x)
Convert potential vorticity to streamfunction - tangent linear.
subroutine, public qg_change_var_ad(conf, fld_in, fld_out)
Change of variable - adjoint.
subroutine, public qg_change_var(conf, fld_in, fld_out)
Change of variable.
subroutine, public qg_fields_check_resolution(fld1, fld2)
Check fields resolution.
subroutine, public convert_x_to_q_tl(geom, x, q)
Convert streamfunction to potential vorticity - tangent linear.
Fortran interface to Variables.
subroutine, public qg_change_var_inv(conf, fld_in, fld_out)
Change of variable - inverse.
subroutine, public qg_change_var_setup(self, vars_in, vars_out)
Linked list implementation.
type(registry_t), public qg_change_var_registry
Linked list interface - defines registry_t type.
subroutine, public qg_change_var_inv_ad(conf, fld_in, fld_out)
Change of variable - inverse adjoint.
subroutine, public convert_x_to_q_ad(geom, q, x)
Convert streamfunction to potential vorticity - adjoint.
subroutine, public qg_fields_copy(self, other, bconly)
Copy fields.
subroutine, public convert_q_to_x_ad(geom, x, q)
Convert potential vorticity to streamfunction - adjoint.