Octopus
linear_solver.F90
Go to the documentation of this file.
1!! Copyright (C) 2004 E.S. Kadantsev, M. Marques
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use batch_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
29 use, intrinsic :: iso_fortran_env
32 use mesh_oct_m
39 use parser_oct_m
42 use smear_oct_m
44 use space_oct_m
47
48 implicit none
49
50 private
51
52 public :: &
62
64 private
65 integer, public :: solver
66 type(preconditioner_t), public :: pre
67 integer :: max_iter
68 type(multigrid_t), allocatable :: mgrid
69 end type linear_solver_t
70
71
73 private
74 type(namespace_t), pointer :: namespace
75 type(linear_solver_t), pointer :: ls
76 type(hamiltonian_elec_t), pointer :: hm
77 type(mesh_t), pointer :: mesh
78 type(states_elec_t), pointer :: st
79 integer :: ist
80 integer :: ik
81 real(real64) :: dshift
82 complex(real64) :: zshift
84
85 type(linear_solver_args_t) :: args
86
87contains
88
89 ! ---------------------------------------------------------
90 subroutine linear_solver_init(this, namespace, gr, states_are_real, mc, space)
91 type(linear_solver_t), intent(out) :: this
92 type(namespace_t), intent(in) :: namespace
93 type(grid_t), intent(inout) :: gr
94 logical, intent(in) :: states_are_real
95 type(multicomm_t), intent(in) :: mc
96 class(space_t), intent(in) :: space
97
98 integer :: fsolver
99 integer :: defsolver
100
101 push_sub(linear_solver_init)
102
103 !%Variable LinearSolver
104 !%Type integer
105 !%Default qmr_symmetric
106 !%Section Linear Response::Solver
107 !%Description
108 !% Method for solving linear equations, which occur for Sternheimer linear
109 !% response and OEP. The solvers vary in speed, reliability (ability to
110 !% converge), and domain of applicability. QMR solvers are most reliable.
111 !%Option bicgstab 4
112 !% Biconjugate gradients stabilized. Slower than <tt>cg</tt>, but more reliable.
113 !% General matrices.
114 !%Option cg 5
115 !% Conjugate gradients. Fast but unreliable. Hermitian matrices only
116 !% (no eta in Sternheimer).
117 !%Option multigrid 7
118 !% Multigrid solver, currently only Gauss-Jacobi (experimental).
119 !% Slow, but fairly reliable. General matrices.
120 !%Option qmr_symmetric 81
121 !% Quasi-minimal residual solver, for (complex) symmetric matrices. [Real symmetric
122 !% is equivalent to Hermitian.] Slightly slower than <tt>bicgstab</tt> but more reliable.
123 !% For Sternheimer, must be real wavefunctions, but can have eta.
124 !%Option qmr_symmetrized 82
125 !% Quasi-minimal residual solver, using the symmetrized form <math>A^\dagger A x = A^\dagger y</math> instead of
126 !% <math>A x = y</math>. Reliable but very slow. General matrices.
127 !%Option qmr_dotp 83
128 !% Quasi-minimal residual solver, for Hermitian matrices, using the
129 !% symmetric algorithm with conjugated dot product (experimental). Slightly slower than <tt>bicgstab</tt>
130 !% but more reliable. Can always be used in Sternheimer.
131 !%Option qmr_general 84
132 !% Quasi-minimal residual solver, for general matrices, using the
133 !% most general form of the algorithm. Slow and unreliable.
134 !%Option sos 9
135 !% Sum over states: the Sternheimer equation is solved by using
136 !% the explicit solution in terms of the ground-state
137 !% wavefunctions. You need unoccupied states to use this method.
138 !% Unlike the other methods, may not give the correct answer.
139 !%Option idrs 11
140 !% This is the "Induced Dimension Reduction", IDR(s) (for s=4). IDR(s) is a robust and efficient short recurrence
141 !% Krylov subspace method for solving large nonsymmetric systems of linear equations. It is described in
142 !% [Peter Sonneveld and Martin B. van Gijzen, SIAM J. Sci. Comput. 31, 1035 (2008)]. We have adapted the code
143 !% released by M. B. van Gizjen [http://ta.twi.tudelft.nl/nw/users/gijzen/IDR.html].
144 !%End
145
146 if (conf%devel_version) then
147 defsolver = option__linearsolver__qmr_dotp
148 else
149 if (states_are_real) then
150 defsolver = option__linearsolver__qmr_symmetric
151 ! in this case, it is equivalent to LS_QMR_DOTP
152 else
153 defsolver = option__linearsolver__qmr_symmetrized
154 end if
155 end if
156 call parse_variable(namespace, "LinearSolver", defsolver, fsolver)
157
158 ! set up pointer for dot product and norm in QMR solvers
159 call mesh_init_mesh_aux(gr)
161 !the last 2 digits select the linear solver
162 this%solver = mod(fsolver, 100)
164 call preconditioner_init(this%pre, namespace, gr, mc, space)
165
166 !%Variable LinearSolverMaxIter
167 !%Type integer
168 !%Default 1000
169 !%Section Linear Response::Solver
170 !%Description
171 !% Maximum number of iterations the linear solver does, even if
172 !% convergence is not achieved.
173 !%End
174 call parse_variable(namespace, "LinearSolverMaxIter", 1000, this%max_iter)
176 call messages_print_with_emphasis(msg='Linear Solver', namespace=namespace)
178 ! solver
179 select case (this%solver)
180 case (option__linearsolver__cg)
181 message(1)='Linear Solver: Conjugate Gradients'
182
183 case (option__linearsolver__bicgstab)
184 message(1)='Linear Solver: Biconjugate Gradients Stabilized'
186 case (option__linearsolver__idrs)
187 message(1)='Linear Solver: IDRS'
188
189 case (option__linearsolver__multigrid)
190 message(1)='Multigrid (currently only Gauss-Jacobi - EXPERIMENTAL)'
191
192 case (option__linearsolver__qmr_symmetric)
193 message(1)='Linear Solver: Quasi-Minimal Residual, for symmetric matrix'
194
195 case (option__linearsolver__qmr_symmetrized)
196 message(1)='Linear Solver: Quasi-Minimal Residual, for symmetrized matrix'
197
198 case (option__linearsolver__qmr_dotp)
199 message(1)='Linear Solver: Quasi-Minimal Residual, symmetric with conjugated dot product'
200
201 case (option__linearsolver__qmr_general)
202 message(1)='Linear Solver: Quasi-Minimal Residual, general algorithm'
203
204 case (option__linearsolver__sos)
205 message(1)='Linear Solver: Sum-over-States'
206 end select
207
208 call messages_info(1, namespace=namespace)
209
210 call messages_print_with_emphasis(namespace=namespace)
211
212 if (this%solver == option__linearsolver__multigrid) then
213 call messages_experimental("Multigrid linear solver")
214
215 safe_allocate(this%mgrid)
216 call multigrid_init(this%mgrid, namespace, space, gr, gr%der, gr%stencil, mc)
217 end if
218
219 if (this%solver == option__linearsolver__qmr_dotp) then
220 call messages_experimental("QMR solver (symmetric with conjugated dot product)")
221 end if
222
223 pop_sub(linear_solver_init)
224 end subroutine linear_solver_init
225
226 ! ---------------------------------------------------------
227 subroutine linear_solver_end(this)
228 type(linear_solver_t), intent(inout) :: this
229 this%solver = -1
230
231 call preconditioner_end(this%pre)
232 if (allocated(this%mgrid)) then
233 call multigrid_end(this%mgrid)
234 safe_deallocate_a(this%mgrid)
235 end if
236
237 end subroutine linear_solver_end
238
239
240 ! ---------------------------------------------------------
241 integer function linear_solver_ops_per_iter(this) result(n)
242 type(linear_solver_t), intent(inout) :: this
243
244 select case (this%solver)
245 case (option__linearsolver__bicgstab)
246 n = 2
247 case default ! LS_CG, LS_MULTIGRID, LS_QMR, LS_SOS
248 n = 1
249 end select
250
251 end function linear_solver_ops_per_iter
252
253 ! ----------------------------------------------------------
254
255 subroutine linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
256 type(namespace_t), intent(in) :: namespace
257 character(len=*), intent(in) :: old_prefix
258 character(len=*), intent(in) :: new_prefix
259
260 call messages_obsolete_variable(namespace, trim(old_prefix)//"LinearSolver", trim(new_prefix)//"LinearSolver")
261 call messages_obsolete_variable(namespace, trim(old_prefix)//"LinearSolverMaxIter", trim(new_prefix)//"LinearSolverMaxIter")
262
263 call preconditioner_obsolete_variables(namespace, old_prefix, new_prefix)
264
266
267#include "undef.F90"
268
269#include "real.F90"
270
271#include "linear_solver_inc.F90"
272
273#include "undef.F90"
274
275#include "complex.F90"
276#include "linear_solver_inc.F90"
277
278
279
280end module linear_solver_oct_m
281
282!! Local Variables:
283!! mode: f90
284!! coding: utf-8
285!! End:
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
type(conf_t), public conf
Global instance of Octopus configuration.
Definition: global.F90:180
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dlinear_solver_solve_hxey_batch(this, namespace, hm, mesh, st, xb, yb, shift, tol, residue, iter_used, occ_response, use_initial_guess)
subroutine, public zlinear_solver_solve_hxey(this, namespace, hm, mesh, st, ist, ik, x, y, shift, tol, residue, iter_used, occ_response)
This subroutine calculates the solution of (H + shift) x = y Typically shift = - eigenvalue + omega.
subroutine, public dlinear_solver_solve_hxey(this, namespace, hm, mesh, st, ist, ik, x, y, shift, tol, residue, iter_used, occ_response)
This subroutine calculates the solution of (H + shift) x = y Typically shift = - eigenvalue + omega.
integer function, public linear_solver_ops_per_iter(this)
subroutine, public linear_solver_end(this)
subroutine, public linear_solver_init(this, namespace, gr, states_are_real, mc, space)
subroutine, public zlinear_solver_solve_hxey_batch(this, namespace, hm, mesh, st, xb, yb, shift, tol, residue, iter_used, occ_response, use_initial_guess)
subroutine, public linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:118
This module defines various routines, operating on mesh functions.
subroutine, public mesh_init_mesh_aux(mesh)
Initialise a pointer to the grid/mesh, that is globally exposed, such that low level mesh operations ...
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:904
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1029
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1069
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public multigrid_end(mgrid)
Definition: multigrid.F90:516
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
Definition: multigrid.F90:187
subroutine, public preconditioner_end(this)
subroutine, public preconditioner_obsolete_variables(namespace, old_prefix, new_prefix)
subroutine, public preconditioner_init(this, namespace, gr, mc, space)
This module is intended to contain "only mathematical" functions and procedures.
Definition: solvers.F90:117