Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
34module 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_list, only : bc_list_t
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 procedure, pass(this) :: solve_coupled => sx_cg_solve_coupled
58 end type sx_cg_t
59
60contains
61
63 subroutine sx_cg_init(this, n, max_iter, M, rel_tol, abs_tol, monitor)
64 class(sx_cg_t), intent(inout) :: this
65 class(pc_t), optional, intent(in), target :: M
66 integer, intent(in) :: n
67 integer, intent(in) :: max_iter
68 real(kind=rp), optional, intent(in) :: rel_tol
69 real(kind=rp), optional, intent(in) :: abs_tol
70 logical, optional, intent(in) :: monitor
71
72 call this%free()
73
74 allocate(this%w(n))
75 allocate(this%r(n))
76 allocate(this%p(n))
77 allocate(this%z(n))
78
79 if (present(m)) then
80 this%M => m
81 end if
82
83 if (present(rel_tol) .and. present(abs_tol) .and. present(monitor)) then
84 call this%ksp_init(max_iter, rel_tol, abs_tol, monitor = monitor)
85 else if (present(rel_tol) .and. present(abs_tol)) then
86 call this%ksp_init(max_iter, rel_tol, abs_tol)
87 else if (present(monitor) .and. present(abs_tol)) then
88 call this%ksp_init(max_iter, abs_tol = abs_tol, monitor = monitor)
89 else if (present(rel_tol) .and. present(monitor)) then
90 call this%ksp_init(max_iter, rel_tol, monitor = monitor)
91 else if (present(rel_tol)) then
92 call this%ksp_init(max_iter, rel_tol = rel_tol)
93 else if (present(abs_tol)) then
94 call this%ksp_init(max_iter, abs_tol = abs_tol)
95 else if (present(monitor)) then
96 call this%ksp_init(max_iter, monitor = monitor)
97 else
98 call this%ksp_init(max_iter)
99 end if
100
101 end subroutine sx_cg_init
102
104 subroutine sx_cg_free(this)
105 class(sx_cg_t), intent(inout) :: this
106
107 call this%ksp_free()
108
109 if (allocated(this%w)) then
110 deallocate(this%w)
111 end if
112
113 if (allocated(this%r)) then
114 deallocate(this%r)
115 end if
116
117 if (allocated(this%p)) then
118 deallocate(this%p)
119 end if
120
121 if (allocated(this%z)) then
122 deallocate(this%z)
123 end if
124
125 nullify(this%M)
126
127 end subroutine sx_cg_free
128
130 function sx_cg_solve(this, Ax, x, f, n, coef, blst, gs_h, niter) result(ksp_results)
131 class(sx_cg_t), intent(inout) :: this
132 class(ax_t), intent(in) :: ax
133 type(field_t), intent(inout) :: x
134 integer, intent(in) :: n
135 real(kind=rp), dimension(n), intent(in) :: f
136 type(coef_t), intent(inout) :: coef
137 type(bc_list_t), intent(inout) :: blst
138 type(gs_t), intent(inout) :: gs_h
139 type(ksp_monitor_t) :: ksp_results
140 integer, optional, intent(in) :: niter
141 real(kind=rp), parameter :: one = 1.0
142 real(kind=rp), parameter :: zero = 0.0
143 integer :: i, iter, max_iter
144 real(kind=rp) :: rnorm, rtr, rtr0, rtz2, rtz1
145 real(kind=rp) :: beta, pap, alpha, alphm, norm_fac
146
147 if (present(niter)) then
148 max_iter = niter
149 else
150 max_iter = this%max_iter
151 end if
152 norm_fac = one / sqrt(coef%volume)
153
154 rtz1 = one
155 do i = 1, n
156 x%x(i,1,1,1) = 0.0_rp
157 this%p(i) = 0.0_rp
158 this%r(i) = f(i)
159 end do
160
161 rtr = glsc3(this%r, coef%mult, this%r, n)
162 rnorm = sqrt(rtr)*norm_fac
163 ksp_results%res_start = rnorm
164 ksp_results%res_final = rnorm
165 ksp_results%iter = 0
166 if(abscmp(rnorm, zero)) return
167 call this%monitor_start('CG')
168 do iter = 1, max_iter
169 call this%M%solve(this%z, this%r, n)
170 rtz2 = rtz1
171 rtz1 = glsc3(this%r, coef%mult, this%z, n)
172
173 beta = rtz1 / rtz2
174 if (iter .eq. 1) beta = zero
175 call add2s1(this%p, this%z, beta, n)
176
177 call ax%compute(this%w, this%p, coef, x%msh, x%Xh)
178 call gs_h%op(this%w, n, gs_op_add)
179 call blst%apply_scalar(this%w, n)
180
181 pap = glsc3(this%w, coef%mult, this%p, n)
182
183 alpha = rtz1 / pap
184 alphm = -alpha
185 do i = 1, n
186 x%x(i,1,1,1) = x%x(i,1,1,1) + alpha * this%p(i)
187 this%r(i) = this%r(i) + alphm * this%w(i)
188 end do
189
190 rtr = glsc3(this%r, coef%mult, this%r, n)
191 if (iter .eq. 1) rtr0 = rtr
192 rnorm = sqrt(rtr) * norm_fac
193 call this%monitor_iter(iter, rnorm)
194 if (rnorm .lt. this%abs_tol) then
195 exit
196 end if
197 end do
198 call this%monitor_stop()
199 ksp_results%res_final = rnorm
200 ksp_results%iter = iter
201 ksp_results%converged = this%is_converged(iter, rnorm)
202 end function sx_cg_solve
203
205 function sx_cg_solve_coupled(this, Ax, x, y, z, fx, fy, fz, &
206 n, coef, blstx, blsty, blstz, gs_h, niter) result(ksp_results)
207 class(sx_cg_t), intent(inout) :: this
208 class(ax_t), intent(in) :: ax
209 type(field_t), intent(inout) :: x
210 type(field_t), intent(inout) :: y
211 type(field_t), intent(inout) :: z
212 integer, intent(in) :: n
213 real(kind=rp), dimension(n), intent(in) :: fx
214 real(kind=rp), dimension(n), intent(in) :: fy
215 real(kind=rp), dimension(n), intent(in) :: fz
216 type(coef_t), intent(inout) :: coef
217 type(bc_list_t), intent(inout) :: blstx
218 type(bc_list_t), intent(inout) :: blsty
219 type(bc_list_t), intent(inout) :: blstz
220 type(gs_t), intent(inout) :: gs_h
221 type(ksp_monitor_t), dimension(3) :: ksp_results
222 integer, optional, intent(in) :: niter
223
224 ksp_results(1) = this%solve(ax, x, fx, n, coef, blstx, gs_h, niter)
225 ksp_results(2) = this%solve(ax, y, fy, n, coef, blsty, gs_h, niter)
226 ksp_results(3) = this%solve(ax, z, fz, n, coef, blstz, gs_h, niter)
227
228 end function sx_cg_solve_coupled
229
230end module cg_sx
231
232
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines various Conjugate Gradient methods.
Definition cg_sx.f90:34
type(ksp_monitor_t) function, dimension(3) sx_cg_solve_coupled(this, ax, x, y, z, fx, fy, fz, n, coef, blstx, blsty, blstz, gs_h, niter)
Standard PCG coupled solve.
Definition cg_sx.f90:207
subroutine sx_cg_free(this)
Deallocate a standard PCG solver.
Definition cg_sx.f90:105
subroutine sx_cg_init(this, n, max_iter, m, rel_tol, abs_tol, monitor)
Initialise a standard PCG solver.
Definition cg_sx.f90:64
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:131
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:51
Definition math.f90:60
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:894
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:657
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 allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:46
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:55
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:68
Defines a canonical Krylov preconditioner.
Definition precon.f90:40