FE-Project
Loading...
Searching...
No Matches
scale_gmres.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
14#include "scaleFElib.h"
16 !-----------------------------------------------------------------------------
17 !
18 !++ used modules
19 !
20 use scale_precision
21 use scale_prc
22
23 !-----------------------------------------------------------------------------
24 implicit none
25 private
26
27 !-----------------------------------------------------------------------------
28 !
29 !++ Public type
30 !
31
32 type, public :: gmres
33 real(rp) :: eps
34 real(rp) :: eps0
35 integer :: n
36 integer :: m
37 integer :: m_out
38 real(rp), allocatable :: v(:,:)
39 real(rp), allocatable :: hj(:)
40 real(rp), allocatable :: g(:)
41 real(rp), allocatable :: r(:,:)
42 real(rp), allocatable :: co(:)
43 real(rp), allocatable :: si(:)
44 real(rp), allocatable :: y(:)
45 contains
46 procedure, public :: init => gmres_init
47 procedure, public :: iterate_pre => gmres_iterate_pre
48 procedure, public :: iterate_step_j => gmres_iterate_step_j
49 procedure, public :: iterate_post => gmres_iterate_post
50 procedure, public :: final => gmres_final
51 end type
52
53 !-----------------------------------------------------------------------------
54 !
55 !++ Public parameters & variables
56 !
57
58 !-----------------------------------------------------------------------------
59 !
60 !++ Private procedure
61 !
62
63 !-----------------------------------------------------------------------------
64 !
65 !++ Private parameters & variables
66 !
67
68 !-----------------------------------------------------------------------------
69
70 ! Private procedure
71 !
72
73 ! Private variable
74 !
75
76contains
77
78subroutine gmres_init( this, N, m, eps, eps0 )
79
80 implicit none
81 class(gmres), intent(inout) :: this
82 integer, intent(in) :: N
83 integer, intent(in) :: m
84 real(RP), intent(in) :: eps
85 real(RP), intent(in) :: eps0
86
87 !--------------------------------------
88
89 this%N = n
90 this%m = m
91 this%EPS = eps
92 this%EPS0 = eps0
93
94 allocate( this%v(n,m+1) )
95 allocate( this%hj(m+1) )
96 allocate( this%g(m+1) )
97 allocate( this%r(m+1,m) )
98 allocate( this%co(m), this%si(m) )
99 allocate( this%y(m) )
100
101 return
102end subroutine gmres_init
103
104!OCL SERIAL
105subroutine gmres_iterate_pre( this, b, w0, is_converged )
106 implicit none
107 class(gmres), intent(inout) :: this
108 real(RP), intent(in) :: b(this%N)
109 real(RP), intent(in) :: w0(this%N)
110 logical, intent(out) :: is_converged
111
112 integer :: i
113 !--------------------------------------
114
115 this%g(1) = 0.0_rp
116 do i=1, this%N
117 this%v(i,1) = b(i) - w0(i)
118 this%g(1) = this%g(1) + this%v(i,1)**2
119 end do
120 this%g(1) = sqrt(this%g(1))
121
122 if ( this%g(1) < this%EPS0 * this%N ) then
123 is_converged = .true.
124 return
125 else
126 is_converged = .false.
127 end if
128
129 do i=1, this%N
130 this%v(i,1) = this%v(i,1) / this%g(1)
131 end do
132
133 this%m_out = min(this%m, this%N)
134 return
135end subroutine gmres_iterate_pre
136
137!OCL SERIAL
138subroutine gmres_iterate_step_j( this, j, wj, is_converged )
139 implicit none
140 class(gmres), intent(inout) :: this
141 integer, intent(in) :: j
142 real(RP), intent(inout) :: wj(this%N)
143 logical, intent(out) :: is_converged
144
145 integer :: i
146 real(RP) :: tmp1, tmp2
147 !--------------------------------------
148
149 do i=1, j
150 this%hj(i) = sum( wj(:) * this%v(:,i) )
151 wj(:) = wj(:) - this%hj(i) * this%v(:,i)
152 end do
153 this%hj(j+1) = sqrt( sum(wj(:)**2) )
154
155 ! if ( abs(this%hj(j+1)) < this%EPS0 ) then
156 ! this%m_out = j
157 ! is_converged = .true.
158 ! return
159 ! else
160 this%v(:,j+1) = wj(:) / this%hj(j+1)
161 ! end if
162
163 this%r(1,j) = this%hj(1)
164 do i=1, j-1
165 tmp1 = this%co(i) * this%r(i,j) + this%si(i) * this%hj(i+1)
166 tmp2 = - this%si(i) * this%r(i,j) + this%co(i) * this%hj(i+1)
167 this%r(i ,j) = tmp1
168 this%r(i+1,j) = tmp2
169 end do
170
171 tmp1 = 1.0_rp / sqrt(this%r(j,j)**2 + this%hj(j+1)**2)
172 this%co(j) = tmp1 * this%r(j,j)
173 this%si(j) = tmp1 * this%hj(j+1)
174
175 this%g(j+1) = - this%si(j) * this%g(j)
176 this%g(j) = this%co(j) * this%g(j)
177
178 this%r(j,j) = this%co(j) * this%r(j,j) + this%si(j) * this%hj(j+1)
179 this%r(j+1,j) = 0.0_rp
180
181 if ( abs(this%g(j+1)) < this%EPS ) then
182 this%m_out = j
183 is_converged = .true.
184 else
185 is_converged = .false.
186 end if
187
188 return
189end subroutine gmres_iterate_step_j
190
191!OCL SERIAL
192subroutine gmres_iterate_post( this, x )
193 implicit none
194 class(gmres), intent(inout) :: this
195 real(RP), intent(inout) :: x(this%N)
196
197 integer :: i, j
198 !--------------------------------------
199
200 do j=this%m_out, 1, -1
201 this%y(j) = this%g(j)
202 do i=j+1, this%m_out
203 this%y(j) = this%y(j) - this%r(j,i) * this%y(i)
204 end do
205 this%y(j) = this%y(j) / this%r(j,j)
206 end do
207
208 x(:) = x(:) + matmul(this%v(:,1:this%m_out), this%y(1:this%m_out))
209
210 return
211end subroutine gmres_iterate_post
212
213subroutine gmres_final( this )
214 implicit none
215 class(gmres), intent(inout) :: this
216 !--------------------------------------
217
218 if ( allocated(this%v) ) then
219 deallocate( this%v )
220 deallocate( this%hj )
221 deallocate( this%g )
222 deallocate( this%r )
223 deallocate( this%co, this%si )
224 deallocate( this%y )
225 end if
226
227 return
228end subroutine gmres_final
229
230end module scale_gmres
module common / GMRES
subroutine gmres_init(this, n, m, eps, eps0)