OOPS
qg_differential_solver_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 
11 use iso_c_binding
12 use fckit_log_module, only: fckit_log
13 use fft_mod
14 use kinds
15 !$ use omp_lib
17 use qg_geom_mod
18 
19 implicit none
20 
21 private
23 ! ------------------------------------------------------------------------------
24 contains
25 ! ------------------------------------------------------------------------------
26 !> Solve a Helmholz equation
27 subroutine solve_helmholz(geom,c,b,x)
28 
29 implicit none
30 
31 ! Passed variables
32 type(qg_geom),intent(in) :: geom !< Geometry
33 real(kind_real),intent(in) :: c !< Coefficient in the linear operator
34 real(kind_real),intent(in) :: b(geom%nx,geom%ny) !< Right hand side
35 real(kind_real),intent(out) :: x(geom%nx,geom%ny) !< Solution
36 
37 ! Local variables
38 integer :: kx,iri,ix,iy
39 real(kind_real) :: v(geom%ny),w(geom%ny),z(geom%ny)
40 real(kind_real) :: bext(geom%nx+2,geom%ny),xext(geom%nx+2,geom%ny)
41 real(kind_real) :: am,bm
42 real(kind_real) :: tmp_in(geom%nx,geom%ny),tmp_out(geom%nx,geom%ny)
43 character(len=2014) :: record
44 
45 ! Transform
46 do iy=1,geom%ny
47  call fft_fwd(geom%nx,b(:,iy),bext(:,iy))
48 enddo
49 
50 ! Solve the tri-diagonal systems
51 do kx=0,geom%nx/2
52  ! am and bm parameters
53  am = c+2.0*(cos(2.0*real(kx,kind_real)*pi/real(geom%nx,kind_real))-1.0)/geom%deltax**2-2.0/geom%deltay**2
54  bm = 1.0/geom%deltay**2
55 
56  do iri=1,2
57  ! v and w parameters
58  v(1) = bm/am
59  do iy=2,geom%ny
60  w(iy) = am-bm*v(iy-1)
61  v(iy) = bm/w(iy)
62  enddo
63 
64  ! bext to z
65  z(1) = bext(2*kx+iri,1)/am
66  do iy=2,geom%ny
67  z(iy) = (bext(2*kx+iri,iy)-bm*z(iy-1))/w(iy)
68  enddo
69 
70  ! z to xext
71  xext(2*kx+iri,geom%ny) = z(geom%ny)
72  do iy=geom%ny-1,1,-1
73  xext(2*kx+iri,iy) = z(iy)-v(iy)*xext(2*kx+iri,iy+1)
74  enddo
75  enddo
76 enddo
77 
78 ! Transform back
79 do iy=1,geom%ny
80  call fft_inv(geom%nx,xext(:,iy),x(:,iy))
81 enddo
82 
83 end subroutine solve_helmholz
84 ! ------------------------------------------------------------------------------
85 !> Solve a Helmholz equation - adjoint
86 subroutine solve_helmholz_ad(geom,c,x,b)
87 
88 implicit none
89 
90 ! Passed variables
91 type(qg_geom),intent(in) :: geom !< Geometry
92 real(kind_real),intent(in) :: c !< Coefficient in the linear operator
93 real(kind_real),intent(in) :: x(geom%nx,geom%ny) !< Solution
94 real(kind_real),intent(inout) :: b(geom%nx,geom%ny) !< Right hand side
95 
96 ! Local variables
97 real(kind_real) :: v(geom%ny),w(geom%ny),z(geom%ny)
98 real(kind_real) :: bext(geom%nx+2,geom%ny),xext(geom%nx+2,geom%ny)
99 real(kind_real) :: btmp(geom%nx,geom%ny)
100 real(kind_real) :: am,bm
101 integer :: kx,iri,iy
102 
103 ! Initialization
104 bext = 0.0
105 
106 ! Transform back
107 do iy=geom%ny,1,-1
108  call fft_fwd(geom%nx,x(:,iy),xext(:,iy))
109 enddo
110 
111 ! Solve the tri-diagonal systems
112 do kx=geom%nx/2,0,-1
113  ! am and bm parameters
114  am = c+2.0*(cos(2.0*real(kx,kind_real)*pi/real(geom%nx,kind_real))-1.0)/geom%deltax**2-2.0/geom%deltay**2
115  bm = 1.0/geom%deltay**2
116 
117  do iri=2,1,-1
118  ! v and w parameters
119  v(1) = bm/am
120  do iy=2,geom%ny
121  w(iy) = am-bm*v(iy-1)
122  v(iy) = bm/w(iy)
123  enddo
124 
125  ! Initialization
126  z = 0.0
127 
128  ! z to xext
129  do iy=1,geom%ny-1
130  xext(2*kx+iri,iy+1) = xext(2*kx+iri,iy+1)-v(iy)*xext(2*kx+iri,iy)
131  z(iy) = z(iy)+xext(2*kx+iri,iy)
132  enddo
133  z(geom%ny) = z(geom%ny)+xext(2*kx+iri,geom%ny)
134 
135  ! bext to z
136  do iy=geom%ny,2,-1
137  bext(2*kx+iri,iy) = bext(2*kx+iri,iy)+z(iy)/w(iy)
138  z(iy-1) = z(iy-1)-bm*z(iy)/w(iy)
139  enddo
140  bext(2*kx+iri,1) = bext(2*kx+iri,1)+z(1)/am
141  enddo
142 enddo
143 
144 ! Transform
145 do iy=geom%ny,1,-1
146  call fft_inv(geom%nx,bext(:,iy),btmp(:,iy))
147 enddo
148 b = b+btmp
149 
150 end subroutine solve_helmholz_ad
151 ! ------------------------------------------------------------------------------
152 !> Horizontal Laplacian operator
153 subroutine laplacian_2d(geom,x,del2x)
154 
155 implicit none
156 
157 ! Passed variables
158 type(qg_geom),intent(in) :: geom !< Geometry
159 real(kind_real),intent(in) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
160 real(kind_real),intent(out) :: del2x(geom%nx,geom%ny,geom%nz) !< Result of applying Laplacian to x
161 
162 ! Local variables
163 integer :: iz
164 
165 !$omp parallel do schedule(static) private(iz)
166 do iz=1,geom%nz
167  ! 5-point laplacian
168  del2x(:,:,iz) = -2.0*x(:,:,iz)*(1.0/geom%deltax**2+1.0/geom%deltay**2)
169  del2x(1:geom%nx-1,:,iz) = del2x(1:geom%nx-1,:,iz)+x(2:geom%nx,:,iz)/geom%deltax**2
170  del2x(geom%nx,:,iz) = del2x(geom%nx,:,iz)+x(1,:,iz)/geom%deltax**2
171  del2x(2:geom%nx,:,iz) = del2x(2:geom%nx,:,iz)+x(1:geom%nx-1,:,iz)/geom%deltax**2
172  del2x(1,:,iz) = del2x(1,:,iz)+x(geom%nx,:,iz)/geom%deltax**2
173  del2x(:,1:geom%ny-1,iz) = del2x(:,1:geom%ny-1,iz)+x(:,2:geom%ny,iz)/geom%deltay**2
174  del2x(:,2:geom%ny,iz) = del2x(:,2:geom%ny,iz)+x(:,1:geom%ny-1,iz)/geom%deltay**2
175 enddo
176 !$omp end parallel do
177 
178 end subroutine laplacian_2d
179 ! ------------------------------------------------------------------------------
180 !> Horizontal Laplacian operator - adjoint
181 subroutine laplacian_2d_ad(geom,del2x,x)
182 
183 implicit none
184 
185 ! Passed variables
186 type(qg_geom),intent(in) :: geom !< Geometry
187 real(kind_real),intent(in) :: del2x(geom%nx,geom%ny,geom%nz) !< Result of applying Laplacian to x
188 real(kind_real),intent(inout) :: x(geom%nx,geom%ny,geom%nz) !< Streamfunction
189 
190 ! Local variables
191 integer :: iz
192 
193 do iz=1,geom%nz
194  ! 5-point laplacian
195  x(:,1:geom%ny-1,iz) = x(:,1:geom%ny-1,iz)+del2x(:,2:geom%ny,iz)/geom%deltay**2
196  x(:,2:geom%ny,iz) = x(:,2:geom%ny,iz)+del2x(:,1:geom%ny-1,iz)/geom%deltay**2
197  x(geom%nx,:,iz) = x(geom%nx,:,iz)+del2x(1,:,iz)/geom%deltax**2
198  x(1:geom%nx-1,:,iz) = x(1:geom%nx-1,:,iz)+del2x(2:geom%nx,:,iz)/geom%deltax**2
199  x(1,:,iz) = x(1,:,iz)+del2x(geom%nx,:,iz)/geom%deltax**2
200  x(2:geom%nx,:,iz) = x(2:geom%nx,:,iz)+del2x(1:geom%nx-1,:,iz)/geom%deltax**2
201  x(:,:,iz) = x(:,:,iz)-2.0*del2x(:,:,iz)*(1.0/geom%deltax**2+1.0/geom%deltay**2)
202 enddo
203 
204 end subroutine laplacian_2d_ad
205 ! ------------------------------------------------------------------------------
206 end module differential_solver_mod
subroutine, public solve_helmholz_ad(geom, c, x, b)
Solve a Helmholz equation - adjoint.
subroutine, public laplacian_2d_ad(geom, del2x, x)
Horizontal Laplacian operator - adjoint.
subroutine, public laplacian_2d(geom, x, del2x)
Horizontal Laplacian operator.
subroutine, public solve_helmholz(geom, c, b, x)
Solve a Helmholz equation.
Fortran module for Eigen FFTs.
Definition: fft_mod.F90:7
subroutine, public fft_inv(kk, pf, pg)
Definition: fft_mod.F90:66
subroutine, public fft_fwd(kk, pg, pf)
Definition: fft_mod.F90:49
real(kind_real), parameter, public pi
Pi.