OOPS
qg_convert_x_to_q_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 streamfunction to potential vorticity
24 subroutine convert_x_to_q(geom,x,x_north,x_south,q)
25 
26 implicit none
27 
28 ! Passed variables
29 type(qg_geom),intent(in) :: geom !< Geometry
30 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
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) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
34 
35 ! Local variables
36 integer :: ix,iy,iz
37 real(kind_real) :: del2x(geom%nx,geom%ny,geom%nz)
38 real(kind_real) :: zz
39 
40 ! Laplacian of the streamfunction
41 call laplacian_2d(geom,x,del2x)
42 
43 ! Vertical differences
44 !$omp parallel do schedule(static) private(iz,iy,ix)
45 do iz=1,geom%nz
46  do iy=1,geom%ny
47  do ix=1,geom%nx
48  q(ix,iy,iz) = del2x(ix,iy,iz)+sum(geom%f(iz,:)*x(ix,iy,:))
49  end do
50  end do
51 end do
52 !$omp end parallel do
53 
54 ! Add the contribution from the boundaries
55 zz = 1.0 / (geom%deltay * geom%deltay)
56 !$omp parallel do schedule(static) private(iz)
57 do iz=1,geom%nz
58  q(:,1,iz) = q(:,1,iz)+x_south(iz)*zz
59  q(:,geom%ny,iz) = q(:,geom%ny,iz)+x_north(iz)*zz
60 enddo
61 !$omp end parallel do
62 
63 ! Add the beta term and the heating term
64 !$omp parallel do schedule(static) private(iy)
65 do iy=1,geom%ny
66  q(:,iy,:) = q(:,iy,:)+geom%bet(iy)
67  q(:,iy,1) = q(:,iy,1)+geom%heat(:,iy)
68 enddo
69 !$omp end parallel do
70 
71 end subroutine convert_x_to_q
72 ! ------------------------------------------------------------------------------
73 !> Convert streamfunction to potential vorticity - tangent linear
74 subroutine convert_x_to_q_tl(geom,x,q)
75 
76 implicit none
77 
78 ! Passed variables
79 type(qg_geom),intent(in) :: geom !< Geometry
80 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
81 real(kind_real),intent(inout) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
82 
83 ! Local variables
84 integer :: ix,iy,iz
85 real(kind_real) :: del2x(geom%nx,geom%ny,geom%nz)
86 
87 ! Laplacian of the streamfunction
88 call laplacian_2d(geom,x,del2x)
89 
90 ! Vertical differences
91 !$omp parallel do schedule(static) private(iz,iy,ix)
92 do iz=1,geom%nz
93  do iy=1,geom%ny
94  do ix=1,geom%nx
95  q(ix,iy,iz) = del2x(ix,iy,iz)+sum(geom%f(iz,:)*x(ix,iy,:))
96  end do
97  end do
98 end do
99 !$omp end parallel do
100 
101 ! Add the contribution from the boundaries (=> TL of this is identity)
102 
103 ! Add the beta term and the heating term (=> TL of this is identity)
104 
105 end subroutine convert_x_to_q_tl
106 ! ------------------------------------------------------------------------------
107 !> Convert streamfunction to potential vorticity - adjoint
108 subroutine convert_x_to_q_ad(geom,q,x)
109 
110 implicit none
111 
112 ! Passed variables
113 type(qg_geom),intent(in) :: geom !< Geometry
114 real(kind_real),intent(in) :: q(geom%nx,geom%ny,geom%nz) !< Potential vorticity
115 real(kind_real),intent(inout) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
116 
117 ! Local variables
118 integer :: ix,iy,iz
119 real(kind_real) :: del2x(geom%nx,geom%ny,geom%nz)
120 
121 ! Add the beta term and the heating term (=> AD of this is identity)
122 
123 ! Add the contribution from the boundaries (=> AD of this is identity)
124 
125 ! Initialization
126 del2x = 0.0
127 
128 ! Vertical differences
129 do iz=1,geom%nz
130  do iy=1,geom%ny
131  do ix=1,geom%nx
132  del2x(ix,iy,iz) = del2x(ix,iy,iz)+q(ix,iy,iz)
133  x(ix,iy,iz) = x(ix,iy,iz)+sum(geom%f(:,iz)*q(ix,iy,:))
134  end do
135  end do
136 end do
137 
138 ! Laplacian of the streamfunction
139 call laplacian_2d_ad(geom,del2x,x)
140 
141 end subroutine convert_x_to_q_ad
142 ! ------------------------------------------------------------------------------
143 end module qg_convert_x_to_q_mod
subroutine, public laplacian_2d_ad(geom, del2x, x)
Horizontal Laplacian operator - adjoint.
subroutine, public laplacian_2d(geom, x, del2x)
Horizontal Laplacian operator.
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.