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