Octopus
kdotp_calc.F90
Go to the documentation of this file.
1!!! Copyright (C) 2008-2010 David Strubbe
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 comm_oct_m
23 use debug_oct_m
24 use global_oct_m
25 use grid_oct_m
29 use mesh_oct_m
33 use mpi_oct_m
37 use space_oct_m
41 use utils_oct_m
43
44 implicit none
45
46 private
47 public :: &
54
55contains
56
57! ---------------------------------------------------------
58 character(len=100) function kdotp_wfs_tag(dir, dir2) result(str)
59 integer, intent(in) :: dir
60 integer, optional, intent(in) :: dir2
61
62 push_sub(kdotp_wfs_tag)
63
64 str = "wfs_" // index2axis(dir)
65 if (present(dir2)) str = trim(str) // "_" // index2axis(dir2)
66
67 pop_sub(kdotp_wfs_tag)
68
69 end function kdotp_wfs_tag
70
71 ! ---------------------------------------------------------
80 subroutine calc_band_velocity(namespace, space, gr, st, hm, pert, velocity)
81 type(namespace_t), intent(in) :: namespace
82 class(space_t), intent(in) :: space
83 type(grid_t), intent(in) :: gr
84 type(states_elec_t), intent(in) :: st
85 type(hamiltonian_elec_t), intent(inout) :: hm
86 class(perturbation_t), intent(inout) :: pert
87 real(real64), intent(out) :: velocity(:,:,:)
88
89 integer :: ik, ist, idir, ib, ind
90 type(wfs_elec_t) :: pert_psib
91 complex(real64) :: dot(st%block_size)
92
93 push_sub(kdotp_calc_band_velocity)
94
95 velocity(:, :, :) = m_zero
96 if(states_are_real(st)) then
97 pop_sub(kdotp_calc_band_velocity)
98 return
99 end if
100
101 call profiling_in("CALC_BAND_VELOCITY")
102
103 do ik = st%d%kpt%start, st%d%kpt%end
104 do ib = st%group%block_start, st%group%block_end
105
106 call st%group%psib(ib, ik)%copy_to(pert_psib)
107
108 do idir = 1, space%periodic_dim
109 call pert%setup_dir(idir)
110 call pert%apply_batch(namespace, space, gr, hm, st%group%psib(ib, ik), pert_psib)
111
112 call zmesh_batch_dotp_vector(gr, st%group%psib(ib, ik), pert_psib, dot)
113
114 do ind = 1, pert_psib%nst
115 ist = st%group%psib(ib, ik)%ist(ind)
116 velocity(idir, ist, ik) = -aimag(dot(ind))
117 end do
118 end do
119 call pert_psib%end()
120 end do
121 end do
122
123 if (st%parallel_in_states .or. st%d%kpt%parallel) then
124 call comm_allreduce(st%st_kpt_mpi_grp, velocity)
125 end if
126
127 call profiling_out("CALC_BAND_VELOCITY")
128
129 pop_sub(kdotp_calc_band_velocity)
130 end subroutine calc_band_velocity
131
132#include "undef.F90"
133#include "real.F90"
134#include "kdotp_calc_inc.F90"
135
136#include "undef.F90"
137#include "complex.F90"
138#include "kdotp_calc_inc.F90"
139
140end module kdotp_calc_oct_m
141
142!! Local Variables:
143!! mode: f90
144!! coding: utf-8
145!! End:
real(real64), parameter, public m_zero
Definition: global.F90:187
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public calc_band_velocity(namespace, space, gr, st, hm, pert, velocity)
Computes the k-point and band-resolved velocity.
Definition: kdotp_calc.F90:174
subroutine, public dkdotp_add_occ(namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
add projection onto occupied states, by sum over states
Definition: kdotp_calc.F90:400
subroutine, public zkdotp_add_occ(namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
add projection onto occupied states, by sum over states
Definition: kdotp_calc.F90:639
subroutine, public dcalc_eff_mass_inv(namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
Computes the effective mass tensor.
Definition: kdotp_calc.F90:302
character(len=100) function, public kdotp_wfs_tag(dir, dir2)
Definition: kdotp_calc.F90:152
subroutine, public zcalc_eff_mass_inv(namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
Computes the effective mass tensor.
Definition: kdotp_calc.F90:541
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
subroutine, public zmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
pure logical function, public states_are_real(st)
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
character pure function, public index2axis(idir)
Definition: utils.F90:202