Octopus
ps_cpi_file.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 debug_oct_m
23 use global_oct_m
27
28 implicit none
29
30 private
31
32 public :: &
36
38 type ps_cpi_file_t
39 ! Components are public by default
40 real(real64) :: zval
41 integer :: no_l_channels
42
43 integer :: nr
44 real(real64) :: a
45
46 real(real64), allocatable :: rofi(:)
47 real(real64), allocatable :: vps(:,:)
48 real(real64), allocatable :: rphi(:,:)
49
50 logical :: core_corrections
51 real(real64), allocatable :: chcore(:)
52 real(real64), allocatable, private :: d1chcore(:)
53 real(real64), allocatable, private :: d2chcore(:)
54 end type ps_cpi_file_t
55
56contains
57
58 ! ---------------------------------------------------------
59 subroutine ps_cpi_file_read(unit, psf)
60 integer, intent(in) :: unit
61 type(ps_cpi_file_t), intent(inout) :: psf
62
63 integer :: i, l, ios, idummy
64 real(real64) :: a, b, c, d
65
66 push_sub(ps_cpi_file_read)
67
68 read(unit, *) psf%zval, psf%no_l_channels
69 ! skip 10 lines
70 do i = 1, 10
71 read(unit, *)
72 end do
73
74 read(unit, *) psf%nr, psf%a
75
76 ! add extra point for the zero
77 psf%nr = psf%nr + 1
78
79 safe_allocate(psf%rofi (1:psf%nr))
80 safe_allocate(psf%vps (1:psf%nr, 1:psf%no_l_channels))
81 safe_allocate(psf%rphi (1:psf%nr, 1:psf%no_l_channels))
82
83 do l = 1, psf%no_l_channels
84 if (l /= 1) read(unit, *)
85
86 do i = 2, psf%nr
87 read(unit, *) idummy, psf%rofi(i), psf%rphi(i, l), psf%vps(i, l)
88 end do
89 end do
90
91 ! read core charge (if present)
92 read(unit, *, iostat=ios) a, b, c, d
93 if (ios == 0) then
94 psf%core_corrections = .true.
95
96 safe_allocate(psf%chcore (1:psf%nr))
97 safe_allocate(psf%d1chcore(1:psf%nr))
98 safe_allocate(psf%d2chcore(1:psf%nr))
99
100 psf% chcore(2) = b
101 psf%d1chcore(2) = c
102 psf%d2chcore(2) = d
103
104 do i = 3, psf%nr
105 read(unit, *) a, psf%chcore(i), psf%d1chcore(i), psf%d2chcore(i)
106 end do
107 else
108 psf%core_corrections = .false.
109 end if
110
111 ! add extra point at zero
112 psf%rofi(1) = m_zero
113 do l = 1, psf%no_l_channels
114 psf%vps(1, l) = first_point_extrapolate(psf%rofi, psf%vps(:, l))
115
116 psf%rphi(1, l) = m_zero
117 end do
118
119 if (psf%core_corrections) then
120 ! At this point, we use the normalization of the siesta format, where
121 ! psf%chcore(:) = 4*pi*\tilde{rho} r**2. As in the Fritz-Haber file we
122 ! have written 4*pi*\tilde{rho}, we multiply by r**2
123 psf%chcore(:) = psf%chcore(:) * psf%rofi(:)**2
124
125 psf%chcore(1) = first_point_extrapolate(psf%rofi, psf%chcore)
126 psf%d1chcore(1) = first_point_extrapolate(psf%rofi, psf%d1chcore)
127 psf%d2chcore(1) = first_point_extrapolate(psf%rofi, psf%d2chcore)
128 end if
129
130 ! WARNING: This should go away
131 psf%vps(:,:) = psf%vps(:,:)*m_two ! convert to Rydbergs
132
134 end subroutine ps_cpi_file_read
135
137 ! ---------------------------------------------------------
138 subroutine ps_cpi_file_end(psf)
139 type(ps_cpi_file_t), intent(inout) :: psf
142
143 safe_deallocate_a(psf%rofi)
144 safe_deallocate_a(psf%vps)
145 safe_deallocate_a(psf%rphi)
147 if (psf%core_corrections) then
148 safe_deallocate_a(psf%chcore)
149 safe_deallocate_a(psf%d1chcore)
150 safe_deallocate_a(psf%d2chcore)
151 end if
153 pop_sub(ps_cpi_file_end)
154 end subroutine ps_cpi_file_end
155
156end module ps_cpi_file_oct_m
157
158!! Local Variables:
159!! mode: f90
160!! coding: utf-8
161!! End:
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
subroutine, public ps_cpi_file_read(unit, psf)
subroutine, public ps_cpi_file_end(psf)
real(real64) function, public first_point_extrapolate(x, y, high_order)
Definition: ps_in_grid.F90:437
First, the contents of the file.
int true(void)