Octopus
hgh_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 blas_oct_m
24 use comm_oct_m
25 use debug_oct_m
26 use global_oct_m
28 use mesh_oct_m
31 use ps_oct_m
34
35 implicit none
36
37 private
38 public :: &
47
48
49
51 private
52 integer :: n_s
53 real(real64), allocatable, public :: dp(:, :)
54 complex(real64), allocatable, public :: zp(:, :)
55 real(real64), public :: h(3, 3) = m_zero
56 real(real64), public :: k(3, 3) = m_zero
57 end type hgh_projector_t
58
59
60contains
61
62 ! ---------------------------------------------------------
63 subroutine hgh_projector_init(hgh_p, sm, reltyp, spec, l, lm, so_strength)
64 type(hgh_projector_t), intent(inout) :: hgh_p
65 type(submesh_t), intent(in) :: sm
66 integer, intent(in) :: reltyp
67 class(pseudopotential_t), intent(in) :: spec
68 integer, intent(in) :: l, lm
69 real(real64), intent(in) :: so_strength
70
71 integer :: i
72
73 push_sub(hgh_projector_init)
74
75 hgh_p%n_s = sm%np
76 if (reltyp == 0) then
77 safe_allocate(hgh_p%dp(1:hgh_p%n_s, 1:3))
78 do i = 1, 3
79 call pseudopotential_real_nl_projector(spec, hgh_p%n_s, sm%rel_x, sm%r, l, lm, i, hgh_p%dp(:, i))
80 end do
81 else
82 safe_allocate(hgh_p%zp(1:hgh_p%n_s, 1:3))
83 do i = 1, 3
84 call pseudopotential_nl_projector(spec, hgh_p%n_s, sm%rel_x, sm%r, l, lm, i, hgh_p%zp(:, i))
85 end do
86 end if
87
88 hgh_p%h(:, :) = spec%ps%h(l, :, :)
89 hgh_p%k(:, :) = spec%ps%k(l, :, :)*so_strength
90
91 pop_sub(hgh_projector_init)
92 end subroutine hgh_projector_init
93
94 ! ---------------------------------------------------------
95 subroutine hgh_projector_end(hgh_p)
96 type(hgh_projector_t), intent(inout) :: hgh_p
97
98 push_sub(hgh_projector_end)
99
100 safe_deallocate_a(hgh_p%dp)
101 safe_deallocate_a(hgh_p%zp)
102
103 pop_sub(hgh_projector_end)
104 end subroutine hgh_projector_end
105
106
107#include "undef.F90"
108#include "real.F90"
109#include "hgh_projector_inc.F90"
110
111#include "undef.F90"
112#include "complex.F90"
113#include "hgh_projector_inc.F90"
115end module hgh_projector_oct_m
116
117!! Local Variables:
118!! mode: f90
119!! coding: utf-8
120!! End:
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
real(real64), parameter, public m_zero
Definition: global.F90:187
subroutine, public zhgh_project_ket(hgh_p, ll, lmax, dim, reltype, uvpsi, ppsi)
THREADSAFE.
subroutine, public hgh_projector_end(hgh_p)
subroutine, public zhgh_project_bra(mesh, sm, hgh_p, dim, reltype, psi, uvpsi)
THREADSAFE.
subroutine, public hgh_projector_init(hgh_p, sm, reltyp, spec, l, lm, so_strength)
subroutine, public zhgh_project(mesh, sm, hgh_p, ll, lmax, dim, psi, ppsi, reltype)
zhgh_project calculates the action of the projector hgh_p on the psi wavefunction....
subroutine, public dhgh_project(mesh, sm, hgh_p, ll, lmax, dim, psi, ppsi, reltype)
dhgh_project calculates the action of the projector hgh_p on the psi wavefunction....
subroutine, public dhgh_project_bra(mesh, sm, hgh_p, dim, reltype, psi, uvpsi)
THREADSAFE.
subroutine, public dhgh_project_ket(hgh_p, ll, lmax, dim, reltype, uvpsi, ppsi)
THREADSAFE.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
Definition: ps.F90:114
subroutine, public pseudopotential_nl_projector(spec, np, x, r, l, lm, i, uV)
This routine returns the non-local projector, built using spherical harmonics.
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...