OOPS
qg_convert_q_to_x_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 
12 use kinds
13 !$ use omp_lib
14 use qg_geom_mod
15 
16 implicit none
17 
18 private
20 ! ------------------------------------------------------------------------------
21 contains
22 ! ------------------------------------------------------------------------------
23 !> Convert potential vorticity to streamfunction
24 subroutine convert_q_to_x(geom,q,x_north,x_south,x)
25 
26 implicit none
27 
28 ! Passed variables
29 type(qg_geom),intent(in) :: geom !< Geometry
30 real(kind_real),intent(in) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
31 real(kind_real),intent(in) :: x_north(geom%nz) !< Streamfunction on northern wall
32 real(kind_real),intent(in) :: x_south(geom%nz) !< Streamfunction on southern wall
33 real(kind_real),intent(out) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
34 
35 ! Local variables
36 integer :: ix,iy,iz
37 real(kind_real) :: q_nobc(geom%nx,geom%ny,geom%nz),pinv_q(geom%nx,geom%ny,geom%nz),pinv_x(geom%nx,geom%ny,geom%nz)
38 
39 ! Subtract the beta term and the heating term
40 !$omp parallel do schedule(static) private(iy)
41 do iy=1,geom%ny
42  q_nobc(:,iy,:) = q(:,iy,:)-geom%bet(iy)
43  q_nobc(:,iy,1) = q_nobc(:,iy,1)-geom%heat(:,iy)
44 enddo
45 !$omp end parallel do
46 
47 ! Subtract the contribution from the boundaries
48 !$omp parallel do schedule(static) private(iz)
49 do iz=1,geom%nz
50  q_nobc(:,1,iz) = q_nobc(:,1,iz)-x_south(iz)/geom%deltay**2
51  q_nobc(:,geom%ny,iz) = q_nobc(:,geom%ny,iz)-x_north(iz)/geom%deltay**2
52 enddo
53 !$omp end parallel do
54 
55 ! Apply ff_pinv
56 !$omp parallel do schedule(static) private(iz,iy,ix)
57 do iz=1,geom%nz
58  do iy=1,geom%ny
59  do ix=1,geom%nx
60  pinv_q(ix,iy,iz) = sum(geom%f_pinv(iz,:)*q_nobc(ix,iy,:))
61  end do
62  end do
63 end do
64 !$omp end parallel do
65 
66 ! Solve Helmholz equation for each layer
67 do iz=1,geom%nz
68  call solve_helmholz(geom,geom%f_d(iz),pinv_q(:,:,iz),pinv_x(:,:,iz))
69 end do
70 
71 ! Apply ff_p
72 !$omp parallel do schedule(static) private(iz,iy,ix)
73 do iz=1,geom%nz
74  do iy=1,geom%ny
75  do ix=1,geom%nx
76  x(ix,iy,iz) = sum(geom%f_p(iz,:)*pinv_x(ix,iy,:))
77  end do
78  end do
79 end do
80 !$omp end parallel do
81 
82 end subroutine convert_q_to_x
83 ! ------------------------------------------------------------------------------
84 !> Convert potential vorticity to streamfunction - tangent linear
85 subroutine convert_q_to_x_tl(geom,q,x)
86 
87 implicit none
88 
89 ! Passed variables
90 type(qg_geom),intent(in) :: geom !< Geometry
91 real(kind_real),intent(in) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
92 real(kind_real),intent(out) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
93 
94 ! Local variables
95 integer :: ix,iy,iz
96 real(kind_real) :: pinv_q(geom%nx,geom%ny,geom%nz),pinv_x(geom%nx,geom%ny,geom%nz)
97 
98 ! Subtract the beta term and the heating term (=> TL of this is identity)
99 
100 ! Subtract the contribution from the boundaries (=> TL of this is identity)
101 
102 ! Apply ff_pinv
103 !$omp parallel do schedule(static) private(iz,iy,ix)
104 do iz=1,geom%nz
105  do iy=1,geom%ny
106  do ix=1,geom%nx
107  pinv_q(ix,iy,iz) = sum(geom%f_pinv(iz,:)*q(ix,iy,:))
108  end do
109  end do
110 end do
111 !$omp end parallel do
112 
113 ! Solve Helmholz equation for each layer
114 do iz=1,geom%nz
115  call solve_helmholz(geom,geom%f_d(iz),pinv_q(:,:,iz),pinv_x(:,:,iz))
116 end do
117 
118 ! Apply ff_p
119 !$omp parallel do schedule(static) private(iz,iy,ix)
120 do iz=1,geom%nz
121  do iy=1,geom%ny
122  do ix=1,geom%nx
123  x(ix,iy,iz) = sum(geom%f_p(iz,:)*pinv_x(ix,iy,:))
124  end do
125  end do
126 end do
127 !$omp end parallel do
128 
129 end subroutine convert_q_to_x_tl
130 ! ------------------------------------------------------------------------------
131 !> Convert potential vorticity to streamfunction - adjoint
132 subroutine convert_q_to_x_ad(geom,x,q)
133 
134 implicit none
135 
136 ! Passed variables
137 type(qg_geom),intent(in) :: geom !< Geometry
138 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
139 real(kind_real),intent(inout) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
140 
141 ! Local variables
142 integer :: ix,iy,iz
143 real(kind_real) :: pinv_q(geom%nx,geom%ny,geom%nz),pinv_x(geom%nx,geom%ny,geom%nz),qtmp(geom%nx,geom%ny,geom%nz)
144 
145 ! Apply ff_p
146 !$omp parallel do schedule(static) private(iz,iy,ix)
147 do iz=1,geom%nz
148  do iy=1,geom%ny
149  do ix=1,geom%nx
150  pinv_x(ix,iy,iz) = sum(geom%f_p(:,iz)*x(ix,iy,:))
151  end do
152  end do
153 end do
154 !$omp end parallel do
155 
156 ! Solve Helmholz equation for each layer
157 pinv_q = 0.0
158 do iz=1,geom%nz
159  call solve_helmholz_ad(geom,geom%f_d(iz),pinv_x(:,:,iz),pinv_q(:,:,iz))
160 end do
161 
162 ! Apply ff_pinv
163 !$omp parallel do schedule(static) private(iz,iy,ix)
164 do iz=1,geom%nz
165  do iy=1,geom%ny
166  do ix=1,geom%nx
167  qtmp(ix,iy,iz) = sum(geom%f_pinv(:,iz)*pinv_q(ix,iy,:))
168  end do
169  end do
170 end do
171 !$omp end parallel do
172 q = q+qtmp
173 
174 ! Subtract the beta term and the heating term (=> AD of this is identity)
175 
176 ! Subtract the contribution from the boundaries (=> AD of this is identity)
177 
178 end subroutine convert_q_to_x_ad
179 ! ------------------------------------------------------------------------------
180 end module qg_convert_q_to_x_mod
qg_geom_mod
Definition: qg_geom_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_geom_mod::qg_geom
Definition: qg_geom_mod.F90:26
qg_convert_q_to_x_mod
Definition: qg_convert_q_to_x_mod.F90:9
differential_solver_mod
Definition: qg_differential_solver_mod.F90:9
differential_solver_mod::solve_helmholz
subroutine, public solve_helmholz(geom, c, b, x)
Solve a Helmholz equation.
Definition: qg_differential_solver_mod.F90:28
qg_convert_q_to_x_mod::convert_q_to_x
subroutine, public convert_q_to_x(geom, q, x_north, x_south, x)
Convert potential vorticity to streamfunction.
Definition: qg_convert_q_to_x_mod.F90:25
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
differential_solver_mod::solve_helmholz_ad
subroutine, public solve_helmholz_ad(geom, c, x, b)
Solve a Helmholz equation - adjoint.
Definition: qg_differential_solver_mod.F90:87