Octopus
subspace.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
21module subspace_oct_m
22 use accel_oct_m
24 use batch_oct_m
26 use blas_oct_m
27 use blacs_oct_m
29 use comm_oct_m
30 use debug_oct_m
31#ifdef HAVE_ELPA
32 use elpa
33#endif
34 use global_oct_m
37 use, intrinsic :: iso_fortran_env
40 use math_oct_m
41 use mesh_oct_m
45 use parser_oct_m
46 use pblas_oct_m
53 use types_oct_m
56
57 implicit none
58
59 private
60
61 public :: &
62 subspace_t, &
65
68 type subspace_t
69 private
70 integer :: method
71 contains
72 procedure :: init => subspace_init
73 end type subspace_t
74
75
76contains
77
79 subroutine subspace_init(this, namespace, st)
80 class(subspace_t), intent(out) :: this
81 type(namespace_t), intent(in) :: namespace
82 type(states_elec_t), intent(in) :: st
83
84 integer :: default
85
86 push_sub(subspace_init)
87
88 !%Variable SubspaceDiagonalization
89 !%Type integer
90 !%Default standard
91 !%Section SCF::Eigensolver
92 !%Description
93 !% Selects the method to perform subspace diagonalization. The
94 !% default is <tt>standard</tt>, unless states parallelization is used,
95 !% when the default is <tt>scalapack</tt>.
96 !% Note that this variable is not parsed in the case of the evolution eigensolver.
97 !%Option none 0
98 !% No subspace diagonalization. WARNING: this will generally give incorrect results.
99 !%Option standard 1
100 !% The standard routine. Can be used with domain parallelization but not
101 !% state parallelization.
102 !%Option scalapack 3
103 !% State-parallelized version using ScaLAPACK. (Requires that
104 !% Octopus was compiled with ScaLAPACK support.)
105 !%End
106
107 default = option__subspacediagonalization__standard
108
109#ifdef HAVE_SCALAPACK
110 if (st%parallel_in_states) default = option__subspacediagonalization__scalapack
111#endif
112
113 call parse_variable(namespace, 'SubspaceDiagonalization', default, this%method)
115 if (.not. varinfo_valid_option('SubspaceDiagonalization', this%method)) then
116 call messages_input_error(namespace, 'SubspaceDiagonalization')
117 end if
118
119 call messages_print_var_option('SubspaceDiagonalization', this%method, namespace=namespace)
120
121 ! some checks for ingenious users
122 if (this%method == option__subspacediagonalization__scalapack) then
123#ifndef HAVE_MPI
124 message(1) = 'The scalapack subspace diagonalization can only be used with MPI parallelization.'
125 call messages_fatal(1, only_root_writes = .true., namespace=namespace)
126#else
127#ifndef HAVE_SCALAPACK
128 message(1) = 'The scalapack subspace diagonalization requires scalapack.'
129 call messages_fatal(1, only_root_writes = .true., namespace=namespace)
130#endif
131 if (st%dom_st_mpi_grp%size == 1) then
132 message(1) = 'The scalapack subspace diagonalization is designed to be used with domain or state parallelization.'
133 call messages_warning(1, namespace=namespace)
134 end if
135
136 if (st%d%kpt%parallel) then
137 message(1) = 'Currently the scalapack subspace diagonalization does not use k-point parallelization.'
138 call messages_warning(1, namespace=namespace)
139 end if
140#endif
141 end if
142
143 pop_sub(subspace_init)
144 end subroutine subspace_init
145
146#include "undef.F90"
147#include "real.F90"
148#include "subspace_inc.F90"
149
150#include "undef.F90"
151#include "complex.F90"
152#include "subspace_inc.F90"
153
154end module subspace_oct_m
155
156!! Local Variables:
157!! mode: f90
158!! coding: utf-8
159!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module contains interfaces for BLACS routines Interfaces are from http:
Definition: blacs.F90:27
This module provides the BLACS processor grid.
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:131
This module provides routines for communicating states when using states parallelization.
subroutine subspace_init(this, namespace, st)
Initialize Subspace type.
Definition: subspace.F90:173
subroutine, public dsubspace_diag(this, namespace, mesh, st, hm, ik, eigenval, diff, nonortho)
Diagonalises the Hamiltonian in the subspace defined by the states.
Definition: subspace.F90:304
subroutine, public zsubspace_diag(this, namespace, mesh, st, hm, ik, eigenval, diff, nonortho)
Diagonalises the Hamiltonian in the subspace defined by the states.
Definition: subspace.F90:915
int true(void)