12 use fckit_configuration_module,
only: fckit_configuration
13 use fckit_log_module,
only: fckit_log
36 #define LISTED_TYPE qg_model_config
39 #include "oops/util/linkedList_i.f"
49 #include "oops/util/linkedList_c.f"
58 type(fckit_configuration),
intent(in) :: f_conf
61 type(duration) :: dtstep
62 character(len=20) :: ststep
63 character(len=160) :: record
64 character(len=:),
allocatable :: str
67 call f_conf%get_or_die(
"tstep",str)
70 self%dt = duration_seconds(dtstep)
71 write(record,*)
'qg_model_setup: dt = ',self%dt
72 call fckit_log%info(record)
86 real(kind_real) :: q(fld%geom%nx,fld%geom%ny,fld%geom%nz),qnew(fld%geom%nx,fld%geom%ny,fld%geom%nz)
87 real(kind_real) :: u(fld%geom%nx,fld%geom%ny,fld%geom%nz),v(fld%geom%nx,fld%geom%ny,fld%geom%nz)
90 if (.not.
allocated(fld%x))
call abor1_ftn(
'qg_model_propagate: x required')
91 if (.not.
allocated(fld%x_north))
call abor1_ftn(
'qg_model_propagate: x_north required')
92 if (.not.
allocated(fld%x_south))
call abor1_ftn(
'qg_model_propagate: x_south required')
93 if (.not.
allocated(fld%q_north))
call abor1_ftn(
'qg_model_propagate: q_north required')
94 if (.not.
allocated(fld%q_south))
call abor1_ftn(
'qg_model_propagate: q_south required')
104 call advect_q(fld%geom,conf%dt,u,v,q,fld%q_north,fld%q_south,qnew)
125 real(kind_real) :: q_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz)
126 real(kind_real) :: u_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz),v_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz)
127 real(kind_real) :: q(fld%geom%nx,fld%geom%ny,fld%geom%nz),qnew(fld%geom%nx,fld%geom%ny,fld%geom%nz)
128 real(kind_real) :: u(fld%geom%nx,fld%geom%ny,fld%geom%nz),v(fld%geom%nx,fld%geom%ny,fld%geom%nz)
133 if (.not.
allocated(traj%x))
call abor1_ftn(
'qg_model_propagate_tl: x trajectory required')
134 if (.not.
allocated(traj%x_north))
call abor1_ftn(
'qg_model_propagate_tl: x_north required')
135 if (.not.
allocated(traj%x_south))
call abor1_ftn(
'qg_model_propagate_tl: x_south required')
136 if (.not.
allocated(traj%q_north))
call abor1_ftn(
'qg_model_propagate_tl: q_north required')
137 if (.not.
allocated(traj%q_south))
call abor1_ftn(
'qg_model_propagate_tl: q_south required')
140 call convert_x_to_q(traj%geom,traj%x,traj%x_north,traj%x_south,q_traj)
143 call convert_x_to_u(traj%geom,traj%x,traj%x_north,traj%x_south,u_traj)
149 if (.not.
allocated(fld%x))
call abor1_ftn(
'qg_model_propagate_tl: x perturbation required')
159 call advect_q_tl(fld%geom,conf%dt,u_traj,v_traj,q_traj,traj%q_north,traj%q_south,u,v,q,qnew)
180 real(kind_real) :: q_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz)
181 real(kind_real) :: u_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz),v_traj(fld%geom%nx,fld%geom%ny,fld%geom%nz)
182 real(kind_real) :: q(fld%geom%nx,fld%geom%ny,fld%geom%nz),qnew(fld%geom%nx,fld%geom%ny,fld%geom%nz)
183 real(kind_real) :: u(fld%geom%nx,fld%geom%ny,fld%geom%nz),v(fld%geom%nx,fld%geom%ny,fld%geom%nz)
188 if (.not.
allocated(traj%x))
call abor1_ftn(
'qg_model_propagate_tl: x trajectory required')
189 if (.not.
allocated(traj%x_north))
call abor1_ftn(
'qg_model_propagate_tl: x_north required')
190 if (.not.
allocated(traj%x_south))
call abor1_ftn(
'qg_model_propagate_tl: x_south required')
191 if (.not.
allocated(traj%q_north))
call abor1_ftn(
'qg_model_propagate_tl: q_north required')
192 if (.not.
allocated(traj%q_south))
call abor1_ftn(
'qg_model_propagate_tl: q_south required')
195 call convert_x_to_q(traj%geom,traj%x,traj%x_north,traj%x_south,q_traj)
198 call convert_x_to_u(traj%geom,traj%x,traj%x_north,traj%x_south,u_traj)
204 if (.not.
allocated(fld%x))
call abor1_ftn(
'qg_model_propagate_tl: x perturbation required')
216 call advect_q_ad(fld%geom,conf%dt,u_traj,v_traj,q_traj,traj%q_north,traj%q_south,qnew,u,v,q)
219 fld%x = 0.0_kind_real
subroutine, public advect_q_tl(geom, dt, u_traj, v_traj, q_traj, q_traj_north, q_traj_south, u, v, q, qnew)
Advect potential vorticity - tangent linear.
subroutine, public advect_q_ad(geom, dt, u_traj, v_traj, q_traj, q_traj_north, q_traj_south, qnew, u, v, q)
Advect potential vorticity - adjoint.
subroutine, public advect_q(geom, dt, u, v, q, q_north, q_south, qnew)
Advect potential vorticity.
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_complete(self, var)
Complete missing fields consistently.
subroutine, public qg_model_propagate(conf, fld)
Perform a timestep of the QG model.
type(registry_t), public qg_model_registry
Linked list interface - defines registry_t type.
subroutine, public qg_model_propagate_tl(conf, traj, fld)
Perform a timestep of the QG model - tangent linear.
subroutine, public qg_model_setup(self, f_conf)
Linked list implementation.
subroutine, public qg_model_propagate_ad(conf, traj, fld)
Perform a timestep of the QG model - adjoint.