Neko  0.8.1
A portable framework for high-order spectral element flow simulations
cg_sx.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021, 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 cg_sx
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, add2s1, abscmp
44  implicit none
45  private
46 
48  type, public, extends(ksp_t) :: sx_cg_t
49  real(kind=rp), allocatable :: w(:)
50  real(kind=rp), allocatable :: r(:)
51  real(kind=rp), allocatable :: p(:)
52  real(kind=rp), allocatable :: z(:)
53  contains
54  procedure, pass(this) :: init => sx_cg_init
55  procedure, pass(this) :: free => sx_cg_free
56  procedure, pass(this) :: solve => sx_cg_solve
57  end type sx_cg_t
58 
59 contains
60 
62  subroutine sx_cg_init(this, n, max_iter, M, rel_tol, abs_tol)
63  class(sx_cg_t), intent(inout) :: this
64  class(pc_t), optional, intent(inout), target :: M
65  integer, intent(in) :: n
66  integer, intent(in) :: max_iter
67  real(kind=rp), optional, intent(inout) :: rel_tol
68  real(kind=rp), optional, intent(inout) :: abs_tol
69 
70  call this%free()
71 
72  allocate(this%w(n))
73  allocate(this%r(n))
74  allocate(this%p(n))
75  allocate(this%z(n))
76 
77  if (present(m)) then
78  this%M => m
79  end if
80 
81  if (present(rel_tol) .and. present(abs_tol)) then
82  call this%ksp_init(max_iter, rel_tol, abs_tol)
83  else if (present(rel_tol)) then
84  call this%ksp_init(max_iter, rel_tol=rel_tol)
85  else if (present(abs_tol)) then
86  call this%ksp_init(max_iter, abs_tol=abs_tol)
87  else
88  call this%ksp_init(max_iter)
89  end if
90 
91  end subroutine sx_cg_init
92 
94  subroutine sx_cg_free(this)
95  class(sx_cg_t), intent(inout) :: this
96 
97  call this%ksp_free()
98 
99  if (allocated(this%w)) then
100  deallocate(this%w)
101  end if
102 
103  if (allocated(this%r)) then
104  deallocate(this%r)
105  end if
106 
107  if (allocated(this%p)) then
108  deallocate(this%p)
109  end if
110 
111  if (allocated(this%z)) then
112  deallocate(this%z)
113  end if
114 
115  nullify(this%M)
116 
117  end subroutine sx_cg_free
118 
120  function sx_cg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
121  class(sx_cg_t), intent(inout) :: this
122  class(ax_t), intent(inout) :: ax
123  type(field_t), intent(inout) :: x
124  integer, intent(in) :: n
125  real(kind=rp), dimension(n), intent(inout) :: f
126  type(coef_t), intent(inout) :: coef
127  type(bc_list_t), intent(inout) :: blst
128  type(gs_t), intent(inout) :: gs_h
129  type(ksp_monitor_t) :: ksp_results
130  integer, optional, intent(in) :: niter
131  real(kind=rp), parameter :: one = 1.0
132  real(kind=rp), parameter :: zero = 0.0
133  integer :: i, iter, max_iter
134  real(kind=rp) :: rnorm, rtr, rtr0, rtz2, rtz1
135  real(kind=rp) :: beta, pap, alpha, alphm, norm_fac
136 
137  if (present(niter)) then
138  max_iter = niter
139  else
140  max_iter = this%max_iter
141  end if
142  norm_fac = one / sqrt(coef%volume)
143 
144  rtz1 = one
145  do i = 1, n
146  x%x(i,1,1,1) = 0.0_rp
147  this%p(i) = 0.0_rp
148  this%r(i) = f(i)
149  end do
150 
151  rtr = glsc3(this%r, coef%mult, this%r, n)
152  rnorm = sqrt(rtr)*norm_fac
153  ksp_results%res_start = rnorm
154  ksp_results%res_final = rnorm
155  ksp_results%iter = 0
156  if(abscmp(rnorm, zero)) return
157 
158  do iter = 1, max_iter
159  call this%M%solve(this%z, this%r, n)
160  rtz2 = rtz1
161  rtz1 = glsc3(this%r, coef%mult, this%z, n)
162 
163  beta = rtz1 / rtz2
164  if (iter .eq. 1) beta = zero
165  call add2s1(this%p, this%z, beta, n)
166 
167  call ax%compute(this%w, this%p, coef, x%msh, x%Xh)
168  call gs_h%op(this%w, n, gs_op_add)
169  call bc_list_apply(blst, this%w, n)
170 
171  pap = glsc3(this%w, coef%mult, this%p, n)
172 
173  alpha = rtz1 / pap
174  alphm = -alpha
175  do i = 1, n
176  x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * this%p(i)
177  this%r(i) = this%r(i) + alphm * this%w(i)
178  end do
179 
180  rtr = glsc3(this%r, coef%mult, this%r, n)
181  if (iter .eq. 1) rtr0 = rtr
182  rnorm = sqrt(rtr) * norm_fac
183  if (rnorm .lt. this%abs_tol) then
184  exit
185  end if
186  end do
187  ksp_results%res_final = rnorm
188  ksp_results%iter = iter
189  end function sx_cg_solve
190 
191 end module cg_sx
192 
193 
Defines a Matrix-vector product.
Definition: ax.f90:34
Defines a boundary condition.
Definition: bc.f90:34
Defines various Conjugate Gradient methods.
Definition: cg_sx.f90:34
subroutine sx_cg_init(this, n, max_iter, M, rel_tol, abs_tol)
Initialise a standard PCG solver.
Definition: cg_sx.f90:63
subroutine sx_cg_free(this)
Deallocate a standard PCG solver.
Definition: cg_sx.f90:95
type(ksp_monitor_t) function sx_cg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter)
Standard PCG solve.
Definition: cg_sx.f90:121
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 add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition: math.f90:574
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Krylov preconditioner.
Definition: precon.f90:34
Base type for a matrix-vector product providing .
Definition: ax.f90:43
A list of boundary conditions.
Definition: bc.f90:102
Standard preconditioned conjugate gradient method (SX version)
Definition: cg_sx.f90:48
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