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