Octopus
kb_projector.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 atom_oct_m
23 use comm_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use grid_oct_m
27 use mesh_oct_m
30 use ps_oct_m
31 use pseudo_oct_m
35
36 implicit none
37
38 private
39 public :: &
48
50 private
51 integer :: n_s
52 integer, public :: n_c
53 real(real64), allocatable, public :: p(:,:)
54 real(real64), allocatable, public :: e(:)
55 end type kb_projector_t
56
57
58contains
59
60 ! ---------------------------------------------------------
61 subroutine kb_projector_init(kb_p, sm, pseudo, l, lm)
62 type(kb_projector_t), intent(inout) :: kb_p
63 type(submesh_t), intent(in) :: sm
64 class(pseudopotential_t), intent(in) :: pseudo
65 integer, intent(in) :: l, lm
66
67 integer :: ic
68
69 push_sub(kb_projector_init)
70
71 kb_p%n_s = sm%np
72 kb_p%n_c = pseudo%ps%projectors_per_l(l+1)
73
74 safe_allocate(kb_p%p (1:kb_p%n_s, 1:max(kb_p%n_c,2)))
75 kb_p%p = m_zero
76 safe_allocate(kb_p%e (1:max(kb_p%n_c,2)))
77 kb_p%e = m_zero
78
79 do ic = 1, kb_p%n_c
80 call pseudopotential_real_nl_projector(pseudo, sm%np, sm%rel_x, sm%r, l, lm, ic, kb_p%p(:, ic))
81
82 kb_p%e(ic) = pseudo%ps%h(l, ic, ic)
83 end do
84
85 ! Average potential for fully relativistic pseudopotentials with SOC
86 ! This is needed to get SOStrength=0 and RelativisticCorrections = 0 to be identical
87 if (pseudo%ps%relativistic_treatment == proj_j_dependent .and. l /= 0) then
88 do ic = 0, kb_p%n_c/2-1
89 kb_p%e(2*ic+1) = kb_p%e(2*ic+1)*real(l+1, real64) /real(2*l+1, real64)
90 kb_p%e(2*ic+2) = kb_p%e(2*ic+2)*real(l, real64) /real(2*l+1, real64)
91 end do
92 end if
93
94 pop_sub(kb_projector_init)
95 end subroutine kb_projector_init
96
97 ! ---------------------------------------------------------
98 subroutine kb_projector_end(kb_p)
99 type(kb_projector_t), intent(inout) :: kb_p
100
101 push_sub(kb_projector_end)
102
103 safe_deallocate_a(kb_p%p)
104 safe_deallocate_a(kb_p%e)
105
106 pop_sub(kb_projector_end)
107 end subroutine kb_projector_end
108
109#include "undef.F90"
110#include "real.F90"
111#include "kb_projector_inc.F90"
112
113#include "undef.F90"
114#include "complex.F90"
115#include "kb_projector_inc.F90"
116
117end module kb_projector_oct_m
118
119!! Local Variables:
120!! mode: f90
121!! coding: utf-8
122!! 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 zkb_project_ket(kb_p, dim, uvpsi, psi)
THREADSAFE.
subroutine, public zkb_project_bra(mesh, sm, kb_p, dim, psi, uvpsi)
THREADSAFE.
subroutine, public dkb_project(mesh, sm, kb_p, dim, psi, ppsi)
dkb_project calculates the action of the projector kb_p on the psi wavefunction. The action of the pr...
subroutine, public kb_projector_end(kb_p)
subroutine, public dkb_project_bra(mesh, sm, kb_p, dim, psi, uvpsi)
THREADSAFE.
subroutine, public zkb_project(mesh, sm, kb_p, dim, psi, ppsi)
zkb_project calculates the action of the projector kb_p on the psi wavefunction. The action of the pr...
subroutine, public dkb_project_ket(kb_p, dim, uvpsi, psi)
THREADSAFE.
subroutine, public kb_projector_init(kb_p, sm, pseudo, l, lm)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
Definition: ps.F90:114
integer, parameter, public proj_j_dependent
Fully-relativistic j-dependent pseudopotentials.
Definition: ps.F90:171
subroutine, public pseudopotential_real_nl_projector(spec, np, x, r, l, lm, i, uV)
This routine returns the non-local projector and its derivative, built using real spherical harmonics...