Octopus
ps_fhi.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 ps_fhi_oct_m
22 use atomic_oct_m
23 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
30 use ps_cpi_oct_m
33
34 implicit none
35
36 private
37 public :: &
38 ps_fhi_t, &
40 ps_fhi_end, &
42
44 type ps_fhi_t
45 ! Components are public by default
46 type(ps_fhi_file_t), allocatable, private :: fhi_file
47 type(ps_cpi_file_t), allocatable, private :: cpi_file
48 type(ps_in_grid_t), allocatable :: ps_grid
49 type(valconf_t), allocatable, private :: conf
50 end type ps_fhi_t
51
52contains
53
54 ! ---------------------------------------------------------
55 subroutine ps_fhi_init(ps_fhi, filename, namespace)
56 type(ps_fhi_t), intent(inout) :: ps_fhi
57 character(len=*), intent(in) :: filename
58 type(namespace_t), intent(in) :: namespace
59
60 integer :: iunit
61 logical :: found
62
63 push_sub(ps_fhi_init)
64
65 safe_allocate(ps_fhi%fhi_file)
66 safe_allocate(ps_fhi%cpi_file)
67 safe_allocate(ps_fhi%ps_grid)
68 safe_allocate(ps_fhi%conf)
69
70 inquire(file = filename, exist = found)
71
72 if (.not. found) then
73 call messages_write("Pseudopotential file '" // trim(filename) // "' not found")
74 call messages_fatal(namespace=namespace)
75 end if
76
77 iunit = io_open(filename, action='read', form='formatted', status='old')
78 call ps_fhi_file_read(iunit, ps_fhi%fhi_file, namespace)
79 call ps_cpi_file_read(iunit, ps_fhi%cpi_file)
80 call io_close(iunit)
81
82 call ps_cpi_file_to_grid(ps_fhi%cpi_file, ps_fhi%ps_grid)
83
84 pop_sub(ps_fhi_init)
85 end subroutine ps_fhi_init
86
87
88 ! ---------------------------------------------------------
89 subroutine ps_fhi_end(ps_fhi)
90 type(ps_fhi_t), intent(inout) :: ps_fhi
91
92 push_sub(ps_fhi_end)
93
94 safe_deallocate_a(ps_fhi%fhi_file)
95
96 call ps_cpi_file_end(ps_fhi%cpi_file)
97 safe_deallocate_a(ps_fhi%cpi_file)
98 safe_deallocate_a(ps_fhi%conf)
99
100 call ps_in_grid_end(ps_fhi%ps_grid)
101 safe_deallocate_a(ps_fhi%ps_grid)
102
103 pop_sub(ps_fhi_end)
104 end subroutine ps_fhi_end
105
106
107 ! ---------------------------------------------------------
108 subroutine ps_fhi_process(ps_fhi, lmax, lloc, namespace)
109 type(ps_fhi_t), intent(inout) :: ps_fhi
110 integer, intent(in) :: lmax, lloc
111 type(namespace_t), intent(in) :: namespace
112
113 push_sub(ps_fhi_process)
115 if (lmax /= ps_fhi%fhi_file%lmax) then
116 message(1) = "Inconsistency in pseudopotential :"
117 write(message(2),'(a,i2,a,i2)') " Input file says lmax = ", lmax, &
118 " but ps file says lmax = ", ps_fhi%fhi_file%lmax
119 call messages_warning(2, namespace=namespace)
120 end if
121 if (lloc /= ps_fhi%fhi_file%lloc) then
122 message(1) = "Inconsistency in pseudopotential :"
123 write(message(2),'(a,i2,a,i2)') " Input file says lloc = ", lloc, &
124 " but ps file says lloc = ", ps_fhi%fhi_file%lloc
125 call messages_warning(2, namespace=namespace)
126 end if
127
128 ! check norm of rphi
129 call ps_in_grid_check_rphi(ps_fhi%ps_grid, namespace)
130
131 ! Fix the local potential. Final argument is the core radius
132 call ps_in_grid_vlocal(ps_fhi%ps_grid, lloc, m_three, namespace)
133
134 ! Calculate kb cosines and norms
135 call ps_in_grid_kb_cosines(ps_fhi%ps_grid, lloc)
136
137 ! Define the KB-projector cut-off radii
138 call ps_in_grid_cutoff_radii(ps_fhi%ps_grid, lloc)
140 ! Calculate KB-projectors
141 call ps_in_grid_kb_projectors(ps_fhi%ps_grid)
143 pop_sub(ps_fhi_process)
144 end subroutine ps_fhi_process
145
146end module ps_fhi_oct_m
147
148!! Local Variables:
149!! mode: f90
150!! coding: utf-8
151!! End:
real(real64), parameter, public m_three
Definition: global.F90:190
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
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 ps_cpi_file_read(unit, psf)
subroutine, public ps_cpi_file_end(psf)
subroutine, public ps_cpi_file_to_grid(cpi_file, ps_grid)
Definition: ps_cpi.F90:197
subroutine, public ps_fhi_file_read(unit, psf, namespace)
subroutine, public ps_fhi_end(ps_fhi)
Definition: ps_fhi.F90:183
subroutine, public ps_fhi_init(ps_fhi, filename, namespace)
Definition: ps_fhi.F90:149
subroutine, public ps_fhi_process(ps_fhi, lmax, lloc, namespace)
Definition: ps_fhi.F90:202
subroutine, public ps_in_grid_kb_projectors(ps)
KB-projectors kb = (vps - vlocal) |phi> * dknorm.
Definition: ps_in_grid.F90:284
subroutine, public ps_in_grid_end(ps)
Definition: ps_in_grid.F90:208
subroutine, public ps_in_grid_cutoff_radii(ps, lloc)
Definition: ps_in_grid.F90:358
subroutine, public ps_in_grid_kb_cosines(ps, lloc)
KB-cosines and KB-norms: dkbcos stores the KB "cosines:" || (v_l - v_local) phi_l ||^2 / < (v_l - v_l...
Definition: ps_in_grid.F90:313
subroutine, public ps_in_grid_vlocal(ps, l_loc, rcore, namespace)
Definition: ps_in_grid.F90:238
subroutine, public ps_in_grid_check_rphi(ps, namespace)
checks normalization of the pseudo wavefunctions
Definition: ps_in_grid.F90:412
remember that the FHI format is basically the CPI format with a header
Definition: ps_fhi.F90:137