OOPS
qg_convert_x_to_uv_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 ageom%ny jurisdiction.
8 
10 
11 use kinds
12 use qg_geom_mod
13 
14 implicit none
15 
16 private
18 ! ------------------------------------------------------------------------------
19 contains
20 ! ------------------------------------------------------------------------------
21 !> Convert streafunction to wind components
22 subroutine convert_x_to_uv(geom,x,x_north,x_south,u,v)
23 
24 implicit none
25 
26 ! Passed variables
27 type(qg_geom),intent(in) :: geom !< Geometry
28 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
29 real(kind_real),intent(in) :: x_north(geom%nz) !< Streamfunction on northern wall
30 real(kind_real),intent(in) :: x_south(geom%nz) !< Streamfunction on southern wall
31 real(kind_real),intent(out) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
32 real(kind_real),intent(out) :: v(geom%nx,geom%ny,geom%nz) !< Meridional wind
33 
34 ! Local variables
35 integer :: iz
36 
37 ! Zonal wind
38 !$omp parallel do schedule(static) private(iz)
39 do iz=1,geom%nz
40  u(:,2:geom%ny,iz) = 0.5*x(:,1:geom%ny-1,iz)/geom%deltay
41  u(:,1,iz) = 0.5*x_south(iz)/geom%deltay
42  u(:,1:geom%ny-1,iz) = u(:,1:geom%ny-1,iz)-0.5*x(:,2:geom%ny,iz)/geom%deltay
43  u(:,geom%ny,iz) = u(:,geom%ny,iz)-0.5*x_north(iz)/geom%deltay
44 enddo
45 !$omp end parallel do
46 
47 ! Meridional wind
48 !$omp parallel do schedule(static) private(iz)
49 do iz=1,geom%nz
50  v(1:geom%nx-1,:,iz) = 0.5*x(2:geom%nx,:,iz)/geom%deltax
51  v(geom%nx,:,iz) = 0.5*x(1,:,iz)/geom%deltax
52  v(2:geom%nx,:,iz) = v(2:geom%nx,:,iz)-0.5*x(1:geom%nx-1,:,iz)/geom%deltax
53  v(1,:,iz) = v(1,:,iz)-0.5*x(geom%nx,:,iz)/geom%deltax
54 enddo
55 !$omp end parallel do
56 
57 end subroutine convert_x_to_uv
58 ! ------------------------------------------------------------------------------
59 !> Convert streafunction to wind components - tangent Linear
60 subroutine convert_x_to_uv_tl(geom,x,u,v)
61 
62 implicit none
63 
64 ! Passed variables
65 type(qg_geom),intent(in) :: geom !< Geometry
66 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
67 real(kind_real),intent(out) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
68 real(kind_real),intent(out) :: v(geom%nx,geom%ny,geom%nz) !< Meridional wind
69 
70 ! Local variables
71 integer :: iz
72 
73 ! Zonal wind
74 !$omp parallel do schedule(static) private(iz)
75 do iz=1,geom%nz
76  u(:,2:geom%ny,iz) = 0.5*x(:,1:geom%ny-1,iz)/geom%deltay
77  u(:,1,iz) = 0.0
78  u(:,1:geom%ny-1,iz) = u(:,1:geom%ny-1,iz)-0.5*x(:,2:geom%ny,iz)/geom%deltay
79 enddo
80 !$omp end parallel do
81 
82 ! Meridional wind
83 !$omp parallel do schedule(static) private(iz)
84 do iz=1,geom%nz
85  v(1:geom%nx-1,:,iz) = 0.5*x(2:geom%nx,:,iz)/geom%deltax
86  v(geom%nx,:,iz) = 0.5*x(1,:,iz)/geom%deltax
87  v(2:geom%nx,:,iz) = v(2:geom%nx,:,iz)-0.5*x(1:geom%nx-1,:,iz)/geom%deltax
88  v(1,:,iz) = v(1,:,iz)-0.5*x(geom%nx,:,iz)/geom%deltax
89 enddo
90 !$omp end parallel do
91 
92 end subroutine convert_x_to_uv_tl
93 ! ------------------------------------------------------------------------------
94 !> Convert streafunction to wind components - adjoint
95 subroutine convert_x_to_uv_ad(geom,u,v,x)
96 
97 implicit none
98 
99 ! Passed variables
100 type(qg_geom),intent(in) :: geom !< Geometry
101 real(kind_real),intent(in) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
102 real(kind_real),intent(in) :: v(geom%nx,geom%ny,geom%nz) !< Meridional wind
103 real(kind_real),intent(inout) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
104 
105 ! Local variables
106 integer :: iz
107 
108 ! Zonal wind
109 do iz=1,geom%nz
110  x(:,2:geom%ny,iz) = x(:,2:geom%ny,iz)-0.5/geom%deltay*u(:,1:geom%ny-1,iz)
111  x(:,1:geom%ny-1,iz) = x(:,1:geom%ny-1,iz)+0.5/geom%deltay*u(:,2:geom%ny,iz)
112 enddo
113 
114 ! Meridional wind
115 do iz=1,geom%nz
116  x(geom%nx,:,iz) = x(geom%nx,:,iz)-0.5/geom%deltax*v(1,:,iz)
117  x(1:geom%nx-1,:,iz) = x(1:geom%nx-1,:,iz)-0.5/geom%deltax*v(2:geom%nx,:,iz)
118  x(1,:,iz) = x(1,:,iz)+0.5/geom%deltax*v(geom%nx,:,iz)
119  x(2:geom%nx,:,iz) = x(2:geom%nx,:,iz)+0.5/geom%deltax*v(1:geom%nx-1,:,iz)
120 enddo
121 
122 end subroutine convert_x_to_uv_ad
123 ! ------------------------------------------------------------------------------
124 end module qg_convert_x_to_uv_mod
qg_geom_mod
Definition: qg_geom_mod.F90:9
qg_geom_mod::qg_geom
Definition: qg_geom_mod.F90:26
qg_convert_x_to_uv_mod::convert_x_to_uv
subroutine, public convert_x_to_uv(geom, x, x_north, x_south, u, v)
Convert streafunction to wind components.
Definition: qg_convert_x_to_uv_mod.F90:23
qg_convert_x_to_uv_mod::convert_x_to_uv_tl
subroutine, public convert_x_to_uv_tl(geom, x, u, v)
Convert streafunction to wind components - tangent Linear.
Definition: qg_convert_x_to_uv_mod.F90:61
qg_convert_x_to_uv_mod
Definition: qg_convert_x_to_uv_mod.F90:9
qg_convert_x_to_uv_mod::convert_x_to_uv_ad
subroutine, public convert_x_to_uv_ad(geom, u, v, x)
Convert streafunction to wind components - adjoint.
Definition: qg_convert_x_to_uv_mod.F90:96