OOPS
qg_convert_x_to_u_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 zonal wind
22 subroutine convert_x_to_u(geom,x,x_north,x_south,u)
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(inout) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
32 
33 ! Local variables
34 integer :: iz
35 
36 !$omp parallel do schedule(static) private(iz)
37 do iz=1,geom%nz
38  u(:,2:geom%ny,iz) = 0.5*x(:,1:geom%ny-1,iz)/geom%deltay
39  u(:,1,iz) = 0.5*x_south(iz)/geom%deltay
40  u(:,1:geom%ny-1,iz) = u(:,1:geom%ny-1,iz)-0.5*x(:,2:geom%ny,iz)/geom%deltay
41  u(:,geom%ny,iz) = u(:,geom%ny,iz)-0.5*x_north(iz)/geom%deltay
42 enddo
43 !$omp end parallel do
44 
45 end subroutine convert_x_to_u
46 ! ------------------------------------------------------------------------------
47 !> Convert streafunction to zonal wind - tangent Linear
48 subroutine convert_x_to_u_tl(geom,x,u)
49 
50 implicit none
51 
52 ! Passed variables
53 type(qg_geom),intent(in) :: geom !< Geometry
54 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
55 real(kind_real),intent(out) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
56 
57 ! Local variables
58 integer :: iz
59 
60 !$omp parallel do schedule(static) private(iz)
61 do iz=1,geom%nz
62  u(:,2:geom%ny,iz) = 0.5*x(:,1:geom%ny-1,iz)/geom%deltay
63  u(:,1,iz) = 0.0
64  u(:,1:geom%ny-1,iz) = u(:,1:geom%ny-1,iz)-0.5*x(:,2:geom%ny,iz)/geom%deltay
65 enddo
66 !$omp end parallel do
67 
68 end subroutine convert_x_to_u_tl
69 ! ------------------------------------------------------------------------------
70 !> Convert streafunction to zonal wind - adjoint
71 subroutine convert_x_to_u_ad(geom,u,x)
72 
73 implicit none
74 
75 ! Passed variables
76 type(qg_geom),intent(in) :: geom !< Geometry
77 real(kind_real),intent(in) :: u(geom%nx,geom%ny,geom%nz) !< Zonal wind
78 real(kind_real),intent(inout) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
79 
80 ! Local variables
81 integer :: iz
82 
83 do iz=1,geom%nz
84  x(:,2:geom%ny,iz) = x(:,2:geom%ny,iz)-0.5/geom%deltay*u(:,1:geom%ny-1,iz)
85  x(:,1:geom%ny-1,iz) = x(:,1:geom%ny-1,iz)+0.5/geom%deltay*u(:,2:geom%ny,iz)
86 enddo
87 
88 end subroutine convert_x_to_u_ad
89 ! ------------------------------------------------------------------------------
90 end module qg_convert_x_to_u_mod
subroutine, public convert_x_to_u(geom, x, x_north, x_south, u)
Convert streafunction to zonal wind.
subroutine, public convert_x_to_u_tl(geom, x, u)
Convert streafunction to zonal wind - tangent Linear.
subroutine, public convert_x_to_u_ad(geom, u, x)
Convert streafunction to zonal wind - adjoint.