TABLE OF CONTENTS


ABINIT/pawpsxml2ab [ Functions ]

[ Top ] [ Functions ]

NAME

 pawpsxml2ab

FUNCTION

  From a XML format pseudopotential file which has already been read in,
  convert to abinit internal datastructures.

COPYRIGHT

 Copyright (C) 2005-2018 ABINIT group (FJ).
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

 psxml  = pseudopotential data structure
 option = 1 for inpsphead call
          2 for pspatm call

OUTPUT

 pspheads data structure is filled

PARENTS

      inpspheads

CHILDREN

      pawpsp_read_header_xml,pawpsp_read_pawheader

SOURCE

 32 #if defined HAVE_CONFIG_H
 33 #include "config.h"
 34 #endif
 35 
 36 #include "abi_common.h"
 37 
 38 
 39 subroutine pawpsxml2ab( psxml, pspheads,option )
 40 
 41  use defs_basis
 42  use defs_datatypes
 43  use m_profiling_abi
 44  use m_errors
 45  use m_pawpsp,only: pawpsp_read_header_xml,pawpsp_read_pawheader
 46  use m_pawxmlps, only : paw_setup_t
 47 
 48 !This section has been created automatically by the script Abilint (TD).
 49 !Do not modify the following lines by hand.
 50 #undef ABI_FUNC
 51 #define ABI_FUNC 'pawpsxml2ab'
 52 !End of the abilint section
 53 
 54  implicit none
 55 
 56 !Arguments ------------------------------------
 57 !scalars
 58  type(paw_setup_t),intent(in) :: psxml
 59  type(pspheader_type),intent(inout) :: pspheads !vz_i
 60  integer ,intent(in) ::option
 61 
 62 !Local variables-------------------------------
 63 !scalars
 64  integer :: ii,il,lloc,mesh_size_ae,mesh_size_proj
 65  real(dp) :: r2well
 66 ! character(len=100) :: xclibxc
 67 ! character(len=500) :: message
 68 !arrays
 69 
 70 ! *********************************************************************
 71 
 72  if(option==1.or.option==2) then
 73    call pawpsp_read_header_xml(lloc,pspheads%lmax,pspheads%pspcod,&
 74 &   pspheads%pspxc,psxml,r2well,pspheads%zionpsp,pspheads%znuclpsp)
 75 
 76    call pawpsp_read_pawheader(pspheads%pawheader%basis_size,&
 77 &   pspheads%lmax,pspheads%pawheader%lmn_size,&
 78 &   pspheads%pawheader%l_size,pspheads%pawheader%mesh_size,&
 79 &   pspheads%pawheader%pawver,psxml,pspheads%pawheader%rpaw,&
 80 &   pspheads%pawheader%rshp,pspheads%pawheader%shape_type)
 81 
 82    pspheads%nproj=0
 83    do il=0,pspheads%lmax
 84      do ii=1,pspheads%pawheader%basis_size
 85        if(psxml%valence_states%state(ii)%ll==il) pspheads%nproj(il)=pspheads%nproj(il)+1
 86      end do
 87    end do
 88    pspheads%nprojso=0
 89    pspheads%pspdat=27061961
 90    pspheads%pspso=0
 91    pspheads%xccc=1
 92    pspheads%title=psxml%atom%symbol
 93 !
 94  end if
 95 
 96  if (option==2) then
 97    write(std_out,*) "paw_setup version= ",psxml%version
 98    write(std_out,*)"atom symbol= ",psxml%atom%symbol,"Z= ",psxml%atom%znucl,&
 99 &   "core= ",psxml%atom%zion,"valence= ",psxml%atom%zval
100    write(std_out,*) "xc_functional  xc_type=",psxml%xc_functional%functionaltype,&
101 &   "xc_name= ",psxml%xc_functional%name
102    write(std_out,*)"generator  type=",psxml%generator%gen,"name= ",psxml%generator%name
103    write(std_out,*)"PAW_radius rpaw=",psxml%rpaw
104    write(std_out,*)"valence_states"
105    do ii=1,psxml%valence_states%nval
106      write(std_out,*)"state n=",psxml%valence_states%state(ii)%nn,&
107 &     "l= ",psxml%valence_states%state(ii)%ll,&
108 &     "f= ",psxml%valence_states%state(ii)%ff,&
109 &     "rc= ",psxml%valence_states%state(ii)%rc,&
110 &     "e= ",psxml%valence_states%state(ii)%ee,&
111 &     "id= ",psxml%valence_states%state(ii)%id
112    end do
113    do ii=1,psxml%ngrid
114      write(std_out,*)"radial_grid  eq= ",psxml%radial_grid(ii)%eq,&
115 &     "a= ",psxml%radial_grid(ii)%aa,&
116 &     "n= ",psxml%radial_grid(ii)%nn,&
117 &     "d= ",psxml%radial_grid(ii)%dd,&
118 &     "b= ",psxml%radial_grid(ii)%bb,&
119 &     "istart= ",psxml%radial_grid(ii)%istart,&
120 &     "iend= ",psxml%radial_grid(ii)%iend,&
121 &     "id= ",psxml%radial_grid(ii)%id
122    end do
123    write(std_out,*)"shape_function  type= ",psxml%shape_function%gtype,&
124 &   "rc= ",psxml%shape_function%rc,&
125 &   "lamb",psxml%shape_function%lamb
126    if(psxml%ae_core_density%tread) then
127      write(std_out,*)"ae_core_density grid=  ",psxml%ae_core_density%grid
128      write(std_out,*)psxml%ae_core_density%data
129    end if
130    if(psxml%pseudo_core_density%tread) then
131      write(std_out,*)"pseudo_core_density grid= ",psxml%pseudo_core_density%grid
132      write(std_out,*)psxml%pseudo_core_density%data
133    end if
134    if(psxml%pseudo_valence_density%tread) then
135      write(std_out,*)"pseudo_valence_densit grid= ",psxml%pseudo_valence_density%grid
136      write(std_out,*)psxml%pseudo_valence_density%data
137    end if
138    if(psxml%zero_potential%tread) then
139      write(std_out,*)"zero_potential grid= ",psxml%zero_potential%grid
140      write(std_out,*)psxml%zero_potential%data
141    end if
142    if(psxml%LDA_minus_half_potential%tread) then
143      write(std_out,*)"LDA_minus_half_potential grid= ",psxml%LDA_minus_half_potential%grid
144      write(std_out,*)psxml%LDA_minus_half_potential%data
145    end if
146    if(psxml%ae_core_kinetic_energy_density%tread) then
147      write(std_out,*)"ae_core_kinetic_energy_density grid= ",psxml%ae_core_kinetic_energy_density%grid
148      write(std_out,*)psxml%ae_core_kinetic_energy_density%data
149    end if
150    if(psxml%pseudo_core_kinetic_energy_density%tread) then
151      write(std_out,*)"pseudo_core_kinetic_energy_density grid= ",psxml%pseudo_core_kinetic_energy_density%grid
152      write(std_out,*)psxml%pseudo_core_kinetic_energy_density%data
153    end if
154    if(psxml%kresse_joubert_local_ionic_potential%tread) then
155      write(std_out,*)"kresse_joubert_local_ionic_potential grid =",psxml%kresse_joubert_local_ionic_potential%grid
156      write(std_out,*)psxml%kresse_joubert_local_ionic_potential%data
157    end if
158    if(psxml%blochl_local_ionic_potential%tread) then
159      write(std_out,*)"blochl_local_ionic_potential grid= ",psxml%blochl_local_ionic_potential%grid
160      write(std_out,*)psxml%blochl_local_ionic_potential%data
161    end if
162    do ii=1,psxml%ngrid
163      if(trim(psxml%ae_partial_wave(1)%grid)==trim(psxml%radial_grid(ii)%id)) then
164        mesh_size_ae=psxml%radial_grid(ii)%iend-psxml%radial_grid(ii)%istart+1
165      end if
166    end do
167    do ii=1,psxml%ngrid
168      if(trim(psxml%projector_function(1)%grid)==trim(psxml%radial_grid(ii)%id)) then
169        mesh_size_proj=psxml%radial_grid(ii)%iend-psxml%radial_grid(ii)%istart+1
170      end if
171    end do
172    do ii=1,psxml%valence_states%nval
173      write(std_out,*)"ae_partial_wave state= ",psxml%ae_partial_wave(ii)%state,&
174 &     "grid= ",psxml%ae_partial_wave(ii)%grid
175      write(std_out,*)psxml%ae_partial_wave(ii)%data(1:mesh_size_ae)
176      write(std_out,*)"pseudo_partial_wave state= ",psxml%pseudo_partial_wave(ii)%state,&
177 &     "grid= ",psxml%pseudo_partial_wave(ii)%grid
178      write(std_out,*)psxml%pseudo_partial_wave(ii)%data(1:mesh_size_ae)
179      write(std_out,*)"projector_function state= ",psxml%projector_function(ii)%state,&
180 &     "grid= ",psxml%projector_function(ii)%grid
181      write(std_out,*)psxml%projector_function(ii)%data(1:mesh_size_proj)
182    end do
183    write(std_out,*)"kinetic_energy_differences"
184    write(std_out,*)psxml%kinetic_energy_differences%data
185  end if
186 
187 
188 end subroutine pawpsxml2ab