OOPS
qg_model_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 duration_mod
12 use fckit_configuration_module, only: fckit_configuration
13 use fckit_log_module,only: fckit_log
14 use iso_c_binding
15 use kinds
22 use qg_fields_mod
23 use random_mod
24 
25 implicit none
26 
27 private
28 public :: qg_model_config
29 public :: qg_model_registry
31 ! ------------------------------------------------------------------------------
33  real(kind_real) :: dt !< Time step (seconds)
34 end type qg_model_config
35 
36 #define LISTED_TYPE qg_model_config
37 
38 !> Linked list interface - defines registry_t type
39 #include "oops/util/linkedList_i.f"
40 
41 !> Global registry
42 type(registry_t) :: qg_model_registry
43 ! ------------------------------------------------------------------------------
44 contains
45 ! ------------------------------------------------------------------------------
46 ! Public
47 ! ------------------------------------------------------------------------------
48 !> Linked list implementation
49 #include "oops/util/linkedList_c.f"
50 ! ------------------------------------------------------------------------------
51 !> Setup model
52 subroutine qg_model_setup(self,f_conf)
53 
54 implicit none
55 
56 ! Passed variables
57 type(qg_model_config),intent(inout) :: self !< Model configuration
58 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
59 
60 ! Local variables
61 type(duration) :: dtstep
62 character(len=20) :: ststep
63 character(len=160) :: record
64 character(len=:),allocatable :: str
65 
66 ! Define time step
67 call f_conf%get_or_die("tstep",str)
68 ststep = str
69 dtstep = trim(ststep)
70 self%dt = duration_seconds(dtstep)
71 write(record,*) 'qg_model_setup: dt = ',self%dt
72 call fckit_log%info(record)
73 
74 end subroutine qg_model_setup
75 ! ------------------------------------------------------------------------------
76 !> Perform a timestep of the QG model
77 subroutine qg_model_propagate(conf,fld)
78 
79 implicit none
80 
81 ! Passed variables
82 type(qg_model_config),intent(in) :: conf !< Model configuration
83 type(qg_fields),intent(inout) :: fld !< State fields
84 
85 ! Local variables
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)
88 
89 ! Check input
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')
95 
96 ! Compute potential vorticity
97 call convert_x_to_q(fld%geom,fld%x,fld%x_north,fld%x_south,q)
98 
99 ! Compute wind
100 call convert_x_to_u(fld%geom,fld%x,fld%x_north,fld%x_south,u)
101 call convert_x_to_v(fld%geom,fld%x,v)
102 
103 ! Advect potential vorticity
104 call advect_q(fld%geom,conf%dt,u,v,q,fld%q_north,fld%q_south,qnew)
105 
106 ! Compute streamfunction
107 call convert_q_to_x(fld%geom,qnew,fld%x_north,fld%x_south,fld%x)
108 
109 ! Complete other fields
110 call qg_fields_complete(fld,'x')
111 
112 end subroutine qg_model_propagate
113 ! ------------------------------------------------------------------------------
114 !> Perform a timestep of the QG model - tangent linear
115 subroutine qg_model_propagate_tl(conf,traj,fld)
116 
117 implicit none
118 
119 ! Passed variables
120 type(qg_model_config),intent(in) :: conf !< Model configuration
121 type(qg_fields),intent(in) :: traj !< Trajectory fields
122 type(qg_fields),intent(inout) :: fld !< Increment fields
123 
124 ! Local variables
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)
129 
130 ! Trajectory
131 
132 ! Check input
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')
138 
139 ! Compute potential vorticity
140 call convert_x_to_q(traj%geom,traj%x,traj%x_north,traj%x_south,q_traj)
141 
142 ! Compute wind
143 call convert_x_to_u(traj%geom,traj%x,traj%x_north,traj%x_south,u_traj)
144 call convert_x_to_v(traj%geom,traj%x,v_traj)
145 
146 ! Perturbation
147 
148 ! Check input
149 if (.not.allocated(fld%x)) call abor1_ftn('qg_model_propagate_tl: x perturbation required')
150 
151 ! Compute potential vorticity
152 call convert_x_to_q_tl(fld%geom,fld%x,q)
153 
154 ! Compute wind
155 call convert_x_to_u_tl(fld%geom,fld%x,u)
156 call convert_x_to_v_tl(fld%geom,fld%x,v)
157 
158 ! Advect potential vorticity
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)
160 
161 ! Compute streamfunction
162 call convert_q_to_x_tl(fld%geom,qnew,fld%x)
163 
164 ! Complete other fields
165 call qg_fields_complete(fld,'x')
166 
167 end subroutine qg_model_propagate_tl
168 ! ------------------------------------------------------------------------------
169 !> Perform a timestep of the QG model - adjoint
170 subroutine qg_model_propagate_ad(conf,traj,fld)
171 
172 implicit none
173 
174 ! Passed variables
175 type(qg_model_config),intent(in) :: conf !< Model configuration
176 type(qg_fields),intent(in) :: traj !< Trajectory fields
177 type(qg_fields),intent(inout) :: fld !< Increment fields
178 
179 ! Local variables
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)
184 
185 ! Trajectory
186 
187 ! Check input
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')
193 
194 ! Compute potential vorticity
195 call convert_x_to_q(traj%geom,traj%x,traj%x_north,traj%x_south,q_traj)
196 
197 ! Compute wind
198 call convert_x_to_u(traj%geom,traj%x,traj%x_north,traj%x_south,u_traj)
199 call convert_x_to_v(traj%geom,traj%x,v_traj)
200 
201 ! Perturbation
202 
203 ! Check input
204 if (.not.allocated(fld%x)) call abor1_ftn('qg_model_propagate_tl: x perturbation required')
205 
206 ! Initialization
207 u = 0.0_kind_real
208 v = 0.0_kind_real
209 q = 0.0_kind_real
210 qnew = 0.0_kind_real
211 
212 ! Compute streamfunction
213 call convert_q_to_x_ad(fld%geom,fld%x,qnew)
214 
215 ! Advect potential vorticity
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)
217 
218 ! Initialize x
219 fld%x = 0.0_kind_real
220 
221 ! Compute wind
222 call convert_x_to_v_ad(fld%geom,v,fld%x)
223 call convert_x_to_u_ad(fld%geom,u,fld%x)
224 
225 ! Compute potential vorticity
226 call convert_x_to_q_ad(fld%geom,q,fld%x)
227 
228 ! Complete other fields
229 call qg_fields_complete(fld,'x')
230 
231 end subroutine qg_model_propagate_ad
232 ! ------------------------------------------------------------------------------
233 end module qg_model_mod
234 
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.