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