Octopus
ps_fhi_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
26
27 implicit none
28
29 private
30
31 public :: &
34
35 ! First, the contents of the file.
36 type ps_fhi_file_t
37 ! Components are public by default
38
39 ! This is the general ABINIT header
40 character(len=256), private :: title
41 real(real64), private :: znucl ! charge of the nucleus
42 real(real64), private :: zion ! valence charge
43 integer, private :: pspdat ! date of creation of PP (DDMMYY)
44 integer, private :: pspcod ! code for the pseudopotential (6 for .fhi)
45 integer, private :: pspxc ! exchange-correlation used to generate the psp
46 integer :: lmax ! Maximum l to use
47 integer :: lloc ! local part of the pseudo to use
48 integer, private :: mmax ! Maximum number of points in real space grid
49 real(real64), private :: r2well ! ??
50
51 ! this is specific for FHI
52 real(real64), private :: rchrg ! The core charge becomes zero beyond rchrg
53 real(real64), private :: fchrg !
54 real(real64), private :: qchrg !
55 end type ps_fhi_file_t
56
57contains
58
59 ! ---------------------------------------------------------
60 subroutine ps_fhi_file_read(unit, psf, namespace)
61 integer, intent(in) :: unit
62 type(ps_fhi_file_t), intent(inout) :: psf
63 type(namespace_t), intent(in) :: namespace
64
65 character(len=3) :: line
66
67 push_sub(read_file_data)
68
69 read(unit, *) psf%title
70 read(unit, *) psf%znucl, psf%zion, psf%pspdat
71 read(unit, *) psf%pspcod, psf%pspxc, psf%lmax, psf%lloc, psf%mmax, psf%r2well
72
73 if (psf%pspcod /= 6) then
74 message(1) = "Inconsistency in pseudopotential file:"
75 write(message(2),'(a,i2)') " expecting pspcod = 6, but found ", psf%pspcod
76 call messages_fatal(2, namespace=namespace)
77 end if
78
79 read(unit, '(a3)') line
80 if (line /= '4--') then ! have non-local core corrections
81 backspace(unit)
82 read(unit, *) psf%rchrg, psf%fchrg, psf%qchrg
83 else
84 psf%rchrg = m_zero
85 psf%fchrg = m_zero
86 psf%qchrg = m_zero
87 end if
88
89 ! skip next 3 lines (#5, #6, #7)
90 read(unit, '(a3)') line
91 read(unit, '(a3)') line
92 read(unit, '(a3)') line
93
94 pop_sub(read_file_data)
95 end subroutine ps_fhi_file_read
96
97end module ps_fhi_file_oct_m
98
99!! Local Variables:
100!! mode: f90
101!! coding: utf-8
102!! End:
real(real64), parameter, public m_zero
Definition: global.F90:187
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_fhi_file_read(unit, psf, namespace)