Neko  0.8.1
A portable framework for high-order spectral element flow simulations
bicgstab.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module bicgstab
35  use num_types, only: rp
37  use precon, only : pc_t
38  use ax_product, only : ax_t
39  use field, only : field_t
40  use coefs, only : coef_t
41  use gather_scatter, only : gs_t, gs_op_add
42  use bc, only : bc_list_t, bc_list_apply
43  use math, only : glsc3, rzero, copy, neko_eps, add2s2, x_update, &
45  use utils, only : neko_error
46  implicit none
47  private
48 
50  type, public, extends(ksp_t) :: bicgstab_t
51  real(kind=rp), allocatable :: p(:)
52  real(kind=rp), allocatable :: p_hat(:)
53  real(kind=rp), allocatable :: r(:)
54  real(kind=rp), allocatable :: s(:)
55  real(kind=rp), allocatable :: s_hat(:)
56  real(kind=rp), allocatable :: t(:)
57  real(kind=rp), allocatable :: v(:)
58  contains
60  procedure, pass(this) :: init => bicgstab_init
62  procedure, pass(this) :: free => bicgstab_free
64  procedure, pass(this) :: solve => bicgstab_solve
65  end type bicgstab_t
66 
67 contains
68 
70  subroutine bicgstab_init(this, n, max_iter, M, rel_tol, abs_tol)
71  class(bicgstab_t), intent(inout) :: this
72  class(pc_t), optional, intent(inout), target :: M
73  integer, intent(in) :: n
74  integer, intent(in) :: max_iter
75  real(kind=rp), optional, intent(inout) :: rel_tol
76  real(kind=rp), optional, intent(inout) :: abs_tol
77 
78 
79  call this%free()
80 
81  allocate(this%p(n))
82  allocate(this%p_hat(n))
83  allocate(this%r(n))
84  allocate(this%s(n))
85  allocate(this%s_hat(n))
86  allocate(this%t(n))
87  allocate(this%v(n))
88  if (present(m)) then
89  this%M => m
90  end if
91 
92  if (present(rel_tol) .and. present(abs_tol)) then
93  call this%ksp_init(max_iter, rel_tol, abs_tol)
94  else if (present(rel_tol)) then
95  call this%ksp_init(max_iter, rel_tol=rel_tol)
96  else if (present(abs_tol)) then
97  call this%ksp_init(max_iter, abs_tol=abs_tol)
98  else
99  call this%ksp_init(max_iter)
100  end if
101 
102  end subroutine bicgstab_init
103 
105  subroutine bicgstab_free(this)
106  class(bicgstab_t), intent(inout) :: this
107 
108  call this%ksp_free()
109 
110  if (allocated(this%v)) then
111  deallocate(this%v)
112  end if
113 
114  if (allocated(this%r)) then
115  deallocate(this%r)
116  end if
117 
118  if (allocated(this%t)) then
119  deallocate(this%t)
120  end if
121 
122  if (allocated(this%p)) then
123  deallocate(this%p)
124  end if
125 
126  if (allocated(this%p_hat)) then
127  deallocate(this%p_hat)
128  end if
129 
130  if (allocated(this%s)) then
131  deallocate(this%s)
132  end if
133 
134  if (allocated(this%s_hat)) then
135  deallocate(this%s_hat)
136  end if
137 
138  nullify(this%M)
139 
140 
141  end subroutine bicgstab_free
142 
144  function bicgstab_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
145  class(bicgstab_t), intent(inout) :: this
146  class(ax_t), intent(inout) :: ax
147  type(field_t), intent(inout) :: x
148  integer, intent(in) :: n
149  real(kind=rp), dimension(n), intent(inout) :: f
150  type(coef_t), intent(inout) :: coef
151  type(bc_list_t), intent(inout) :: blst
152  type(gs_t), intent(inout) :: gs_h
153  type(ksp_monitor_t) :: ksp_results
154  integer, optional, intent(in) :: niter
155  integer :: iter, max_iter
156  real(kind=rp) :: rnorm, rtr, norm_fac, gamma
157  real(kind=rp) :: beta, alpha, omega, rho_1, rho_2
158 
159  if (present(niter)) then
160  max_iter = niter
161  else
162  max_iter = this%max_iter
163  end if
164  norm_fac = 1.0_rp / sqrt(coef%volume)
165 
166  associate(r => this%r, t => this%t, s => this%s, v => this%v, p => this%p, &
167  s_hat => this%s_hat, p_hat => this%p_hat)
168 
169  call rzero(x%x, n)
170  call copy(r, f, n)
171 
172  rtr = sqrt(glsc3(r, coef%mult, r, n))
173  rnorm = rtr * norm_fac
174  gamma = rnorm * this%rel_tol
175  ksp_results%res_start = rnorm
176  ksp_results%res_final = rnorm
177  ksp_results%iter = 0
178  if(abscmp(rnorm, 0.0_rp)) return
179  do iter = 1, max_iter
180 
181  rho_1 = glsc3(r, coef%mult, f ,n)
182 
183  if (abs(rho_1) .lt. neko_eps) then
184  call neko_error('Bi-CGStab rho failure')
185  end if
186 
187  if (iter .eq. 1) then
188  call copy(p, r, n)
189  else
190  beta = (rho_1 / rho_2) * (alpha / omega)
191  call p_update(p, r, v, beta, omega, n)
192  end if
193 
194  call this%M%solve(p_hat, p, n)
195  call ax%compute(v, p_hat, coef, x%msh, x%Xh)
196  call gs_h%op(v, n, gs_op_add)
197  call bc_list_apply(blst, v, n)
198  alpha = rho_1 / glsc3(f, coef%mult, v, n)
199  call copy(s, r, n)
200  call add2s2(s, v, -alpha, n)
201  rtr = glsc3(s, coef%mult, s, n)
202  rnorm = sqrt(rtr) * norm_fac
203  if (rnorm .lt. this%abs_tol .or. rnorm .lt. gamma) then
204  call add2s2(x%x, p_hat, alpha,n)
205  exit
206  end if
207 
208  call this%M%solve(s_hat, s, n)
209  call ax%compute(t, s_hat, coef, x%msh, x%Xh)
210  call gs_h%op(t, n, gs_op_add)
211  call bc_list_apply(blst, t, n)
212  omega = glsc3(t, coef%mult, s, n) &
213  / glsc3(t, coef%mult, t, n)
214  call x_update(x%x, p_hat, s_hat, alpha, omega, n)
215  call copy(r, s, n)
216  call add2s2(r, t, -omega, n)
217 
218  rtr = glsc3(r, coef%mult, r, n)
219  rnorm = sqrt(rtr) * norm_fac
220  if (rnorm .lt. this%abs_tol .or. rnorm .lt. gamma) then
221  exit
222  end if
223 
224  if (omega .lt. neko_eps) then
225  call neko_error('Bi-CGstab omega failure')
226  end if
227  rho_2 = rho_1
228 
229  end do
230  end associate
231  ksp_results%res_final = rnorm
232  ksp_results%iter = iter
233  end function bicgstab_solve
234 
235 end module bicgstab
236 
237 
Defines a Matrix-vector product.
Definition: ax.f90:34
Defines a boundary condition.
Definition: bc.f90:34
Defines various Bi-Conjugate Gradient Stabilized methods.
Definition: bicgstab.f90:34
subroutine bicgstab_init(this, n, max_iter, M, rel_tol, abs_tol)
Constructor.
Definition: bicgstab.f90:71
subroutine bicgstab_free(this)
Deallocate a standard BiCGSTAB solver.
Definition: bicgstab.f90:106
type(ksp_monitor_t) function bicgstab_solve(this, Ax, x, f, n, coef, blst, gs_h, niter)
Bi-Conjugate Gradient Stabilized method solve.
Definition: bicgstab.f90:145
Coefficients.
Definition: coef.f90:34
Defines a field.
Definition: field.f90:34
Gather-scatter.
Implements the base abstract type for Krylov solvers plus helper types.
Definition: krylov.f90:34
integer, parameter, public ksp_max_iter
Maximum number of iters.
Definition: krylov.f90:50
Definition: math.f90:60
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition: math.f90:810
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:777
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:211
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition: math.f90:67
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:167
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: math.f90:589
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:762
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Krylov preconditioner.
Definition: precon.f90:34
Utilities.
Definition: utils.f90:35
Base type for a matrix-vector product providing .
Definition: ax.f90:43
A list of boundary conditions.
Definition: bc.f90:102
Standard preconditioned Bi-Conjugate Gradient Stabilized method.
Definition: bicgstab.f90:50
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
Type for storing initial and final residuals in a Krylov solver.
Definition: krylov.f90:55
Base abstract type for a canonical Krylov method, solving .
Definition: krylov.f90:65
Defines a canonical Krylov preconditioner.
Definition: precon.f90:40