Octopus
states_elec_calc.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
22 use accel_oct_m
24 use batch_oct_m
26 use blas_oct_m
27 use blacs_oct_m
28 use iso_c_binding
29 use comm_oct_m
30 use debug_oct_m
34 use global_oct_m
35 use grid_oct_m
40 use math_oct_m
42 use mesh_oct_m
45 use mpi_oct_m
49 use pblas_oct_m
55 use space_oct_m
60 use types_oct_m
62
63 implicit none
64
65 private
66
67 public :: &
91
92 interface states_elec_rotate
94 end interface states_elec_rotate
95
96contains
97
98 ! ---------------------------------------------------------
100 !
101 subroutine states_elec_orthogonalize(st, namespace, mesh)
102 type(states_elec_t), intent(inout) :: st
103 type(namespace_t), intent(in) :: namespace
104 class(mesh_t), intent(in) :: mesh
105
106 integer :: ik
107
109
110 do ik = st%d%kpt%start, st%d%kpt%end
111 if (states_are_real(st)) then
112 call dstates_elec_orthogonalization_full(st, namespace, mesh, ik)
113 else
114 call zstates_elec_orthogonalization_full(st, namespace, mesh, ik)
115 end if
116 end do
117
119 end subroutine states_elec_orthogonalize
120
121 ! ---------------------------------------------------------
123 !
124 subroutine states_elec_calc_norms(grid, kpoints, st, norm_ks)
125 type(grid_t), intent(in) :: grid
126 type(kpoints_t), intent(in) :: kpoints
127 type(states_elec_t), intent(in) :: st
128 real(real64), contiguous, intent(out) :: norm_ks(:, :)
129
130 integer :: ik_ispin, iblock, minst, maxst
131
132 push_sub(states_elec_calc_norms)
133
134 ! If state or k is parallel, norm_ks is only partially touched, then gets summed over processes,
135 ! hence initialise to zero
136 norm_ks = m_zero
137
138 do ik_ispin = st%d%kpt%start, st%d%kpt%end
139 do iblock = st%group%block_start, st%group%block_end
140 minst = states_elec_block_min(st, iblock)
141 maxst = states_elec_block_max(st, iblock)
142 ! Compute norm, defering grid communication until outside the loops
143 call mesh_batch_nrm2(grid, st%group%psib(iblock, ik_ispin), norm_ks(minst:maxst, ik_ispin), reduce=.false.)
144 enddo
145 enddo
146
147 ! Reduction over domains
148 if (grid%parallel_in_domains) then
149 norm_ks = norm_ks**2
150 call grid%allreduce(norm_ks)
151 norm_ks = sqrt(norm_ks)
152 end if
153
154 ! Reduction over states and k-points
155 call comm_allreduce(st%st_kpt_mpi_grp, norm_ks)
156
158
159 end subroutine states_elec_calc_norms
160
161
162#include "undef.F90"
163#include "real.F90"
164#include "states_elec_calc_inc.F90"
165
166#include "undef.F90"
167#include "complex.F90"
168#include "states_elec_calc_inc.F90"
169#include "undef.F90"
170
171end module states_elec_calc_oct_m
172
173
174!! Local Variables:
175!! mode: f90
176!! coding: utf-8
177!! End:
double sqrt(double __x) __attribute__((__nothrow__
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 contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
Definition: global.F90:187
This module implements the underlying real-space grid.
Definition: grid.F90:117
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
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
Definition: mesh_batch.F90:177
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
Some general things and nomenclature:
Definition: par_vec.F90:171
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:131
pure logical function, public states_are_real(st)
subroutine, public zstates_elec_matrix(st1, st2, mesh, aa)
subroutine, public states_elec_orthogonalize(st, namespace, mesh)
Orthonormalizes nst orbitals in mesh (honours state parallelization).
subroutine, public zstates_elec_orthogonalize_single(st, mesh, nst, iqn, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, against_all)
ofthogonalize a single wave function against a set of states
subroutine, public dstates_elec_matrix(st1, st2, mesh, aa)
subroutine, public zstates_elec_calc_projections(st, gs_st, namespace, mesh, ik, proj, gs_nst)
This routine computes the projection between two set of states.
subroutine, public dstates_elec_orthogonalization_full(st, namespace, mesh, ik)
Orthonormalizes nst orbitals in mesh (honours state parallelization).
subroutine, public dstates_elec_rrqr_decomposition(st, namespace, mesh, nst, root, ik, jpvt)
Perform RRQR on the transpose states stored in the states object and return the pivot vector.
real(real64) function, public zstates_elec_residue(mesh, dim, hf, ee, ff)
subroutine, public zstates_elec_calc_orth_test(st, namespace, mesh, kpoints)
subroutine, public dstates_elec_orthogonalize_single_batch(st, mesh, nst, iqn, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, against_all)
orthogonalize a single wave function against a set of states
subroutine, public zstates_elec_rrqr_decomposition(st, namespace, mesh, nst, root, ik, jpvt)
Perform RRQR on the transpose states stored in the states object and return the pivot vector.
subroutine, public states_elec_calc_norms(grid, kpoints, st, norm_ks)
Compute the norms of the Kohn-Sham orbitals.
subroutine, public dstates_elec_calc_projections(st, gs_st, namespace, mesh, ik, proj, gs_nst)
This routine computes the projection between two set of states.
subroutine, public zstates_elec_orthogonalization(mesh, nst, dim, psi, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, gs_scheme)
Orthonormalizes phi to the nst orbitals psi.
subroutine, public dstates_elec_calc_orth_test(st, namespace, mesh, kpoints)
subroutine, public zstates_elec_calc_overlap(st, mesh, ik, overlap)
Computes the overlap matrix of the Kohn-Sham states with k-point ik.
subroutine, public dstates_elec_orthogonalization(mesh, nst, dim, psi, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, gs_scheme)
Orthonormalizes phi to the nst orbitals psi.
subroutine, public dstates_elec_orthogonalize_single(st, mesh, nst, iqn, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, against_all)
ofthogonalize a single wave function against a set of states
subroutine dstates_elec_rotate(st, namespace, mesh, uu, ik)
subroutine, public zstates_elec_orthogonalize_single_batch(st, mesh, nst, iqn, phi, normalize, mask, overlap, norm, Theta_fi, beta_ij, against_all)
orthogonalize a single wave function against a set of states
subroutine, public zstates_elec_orthogonalization_full(st, namespace, mesh, ik)
Orthonormalizes nst orbitals in mesh (honours state parallelization).
real(real64) function, public dstates_elec_residue(mesh, dim, hf, ee, ff)
subroutine, public dstates_elec_calc_overlap(st, mesh, ik, overlap)
Computes the overlap matrix of the Kohn-Sham states with k-point ik.
subroutine zstates_elec_rotate(st, namespace, mesh, uu, ik)
This module handles spin dimensions of the states and the k-point distribution.
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
This module provides routines for communicating states when using states parallelization.