TABLE OF CONTENTS


ABINIT/m_psxml2ab [ Modules ]

[ Top ] [ Modules ]

NAME

 m_psxml2ab

FUNCTION

  From a SIESTA XML format pseudopotential file
  convert to abinit internal datastructures for pspheader.

COPYRIGHT

 Copyright (C) 2005-2024 ABINIT group (MJV).
 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

OUTPUT

 psphead = psp information structure
 atmsymb = atomic symbol

SOURCE

 25 #if defined HAVE_CONFIG_H
 26 #include "config.h"
 27 #endif
 28 
 29 #include "abi_common.h"
 30 
 31 module m_psxml2ab
 32 
 33  use defs_basis
 34  use m_abicore
 35  use m_errors
 36 #ifdef HAVE_LIBPSML
 37  use m_psml
 38  use m_psml_api
 39 #endif
 40 
 41  use defs_datatypes, only : pspheader_type
 42  use m_fstrings,     only : yesno
 43 
 44 implicit none
 45 
 46 private
 47 
 48 #ifdef HAVE_LIBPSML
 49 public :: psxml2abheader
 50 !public :: psxml2abfull
 51 
 52 CONTAINS
 53 
 54 subroutine psxml2abheader(psxmlfile, psphead, atmsymb, creator, iwrite)
 55 
 56 !Arguments ------------------------------------
 57 !scalars
 58  integer, intent(in) :: iwrite
 59  character(len=fnlen), intent(in) :: psxmlfile
 60  type(pspheader_type),intent(inout) :: psphead
 61  character(len=3), intent(out) :: atmsymb
 62  character(len=30), intent(out) :: creator
 63 
 64 !Local variables-------------------------------
 65 !scalars
 66  character(len=500) :: flavor, frmt_str, label, message, relat, xc_name
 67  character(len=1) :: g1,g2
 68  integer :: dd,dm,dy
 69  integer :: il,lmm,ii
 70  integer :: nxc
 71  integer :: ishell, nvshells, shell_l, shell_n
 72  integer :: iproj, ll, ll_previous, nprojs, nprojsr, nprojso
 73  integer,parameter :: n1xccc_default=2501
 74  logical :: has_nlcc, has_spin
 75  real(dp) :: ekb
 76  type(ps_t) :: psxml
 77 !arrays
 78  integer, allocatable :: idx_sr(:), idx_so(:)
 79  real(dp),allocatable :: zeld(:),zelu(:)
 80 
 81 ! *********************************************************************
 82 
 83  call ps_destroy(psxml)
 84  call psml_reader(psxmlfile, psxml, debug=.true.)
 85 
 86  psphead%pspdat = 0
 87  call ps_Provenance_Get(psxml, 1, creator=creator, date=message)
 88  read (message(1:10), '(I4,A1,I2,A1,I2)', err=10) &
 89 &  dy, g1, dm, g2, dd
 90  psphead%pspdat = MODULO(dy,100) * 10000 + dm * 100 + dd
 91 
 92 10 continue
 93 
 94  call ps_PseudoAtomSpec_Get(psxml, &
 95 &  atomic_symbol=atmsymb, atomic_label=label, &
 96 &  atomic_number=psphead%znuclpsp, z_pseudo=psphead%zionpsp, &
 97 &  pseudo_flavor=flavor, relativity=relat, &
 98 &  spin_dft=has_spin, core_corrections=has_nlcc)
 99 
100  psphead%pspcod = 9
101 
102 ! impose libxc coding for pspxc
103  psphead%pspxc = 0
104  call ps_ExchangeCorrelation_Get(psxml, n_libxc_functionals=nxc)
105  do ii=1, min(2,nxc)
106    call ps_LibxcFunctional_Get(psxml, ii, code=dd)
107    psphead%pspxc = dd + psphead%pspxc*1000
108  end do
109  psphead%pspxc = -psphead%pspxc
110 
111  call ps_ValenceConfiguration_Get(psxml, nshells=nvshells)
112 
113  if (iwrite == 1) then
114    write (message,'(a,a)') '- psxml2ab: ps_PseudoFlavor ', trim(message)
115 !  call wrtout(ab_out,  message,'COLL')
116    call wrtout(std_out,  message,'COLL')
117    write (message,'(a,I5)') '- psxml2ab: ps_NLibxcFunctionals ', nxc
118 !  call wrtout(ab_out,  message,'COLL')
119    call wrtout(std_out,  message,'COLL')
120    write (message,'(a,I5)') '- psxml2ab: ps_NValenceShells ', nvshells
121 !  call wrtout(ab_out,  message,'COLL')
122    call wrtout(std_out,  message,'COLL')
123  end if
124 
125  ABI_MALLOC(zeld, (nvshells))
126  ABI_MALLOC(zelu, (nvshells))
127  zeld = zero
128  zelu = zero
129  frmt_str="(a"
130  do ishell = 1, nvshells
131    if (has_spin) then
132      call ps_ValenceShell_Get(psxml, ishell, n=shell_n, l=shell_l, &
133 &    occ_up=zelu(ishell), occ_down=zeld(ishell))
134    else
135      call ps_ValenceShell_Get(psxml, ishell, n=shell_n, l=shell_l, &
136 &    occupation=zeld(ishell))
137    end if
138    if (iwrite == 1) then
139      write (message,'(a,I5)') '- psxml2ab: ps_ValenceShellN ', shell_n
140 !    call wrtout(ab_out,  message,'COLL')
141      call wrtout(std_out,  message,'COLL')
142      write (message,'(a,I5)') '- psxml2ab: ps_ValenceShellL ', shell_l
143 !    call wrtout(ab_out,  message,'COLL')
144      call wrtout(std_out,  message,'COLL')
145    end if
146    frmt_str = trim(frmt_str) // ", F10.3"
147  end do
148  frmt_str = trim(frmt_str) // ")"
149  if (iwrite == 1) then
150    write (message,frmt_str) '- psxml2ab: zeld ', zeld
151 !  call wrtout(ab_out,  message,'COLL')
152    call wrtout(std_out,  message,'COLL')
153    write (message,frmt_str) '- psxml2ab: zelu ', zelu
154 !  call wrtout(ab_out,  message,'COLL')
155    call wrtout(std_out,  message,'COLL')
156  end if
157 
158 !    Find the number of projectors per angular momentum shell
159  nprojs = 0
160  nprojsr = 0
161  nprojso = 0
162  psphead%nproj(:)=0
163  call ps_NonlocalProjectors_Filter(psxml, set=SET_NONREL, number=nprojs)
164  if (nprojs > 0) then
165    call ps_NonlocalProjectors_Filter(psxml, set=SET_NONREL, indexes=idx_sr)
166    do iproj = 1, nprojs
167      if (iwrite == 1) then
168        write (message,'(a,2I5)') '- psxml2ab: iproj, idx for nonrel ', iproj, idx_sr(iproj)
169 !      call wrtout(ab_out,  message,'COLL')
170        call wrtout(std_out,  message,'COLL')
171      end if
172      call ps_Projector_Get(psxml, idx_sr(iproj), l=il)
173      psphead%nproj(il) = psphead%nproj(il) + 1
174    end do
175  else
176    call ps_NonlocalProjectors_Filter(psxml, set=SET_SREL, number=nprojsr)
177    if (nprojsr > 0) then
178      call ps_NonlocalProjectors_Filter(psxml, set=SET_SREL, indexes=idx_sr)
179      do iproj = 1, nprojsr
180        if (iwrite == 1) then
181          write (message,'(a,2I5)') '- psxml2ab: iproj, idx for srel ', iproj, idx_sr(iproj)
182 !        call wrtout(ab_out,  message,'COLL')
183          call wrtout(std_out,  message,'COLL')
184        end if
185        call ps_Projector_Get(psxml, idx_sr(iproj), l=il)
186        psphead%nproj(il) = psphead%nproj(il) + 1
187      end do
188    else
189      ABI_BUG('Your psml potential should have either scalar- or non- relativistic projectors')
190    end if
191  end if
192 
193  psphead%nprojso(:)=0
194  call ps_NonlocalProjectors_Filter(psxml, set=SET_SO, number=nprojso, &
195 &  indexes=idx_so)
196  do iproj = 1, nprojso
197    if (iwrite == 1) then
198      write (message,'(a,2I5)') '- psxml2ab: iproj, idx for soc ', iproj, idx_so(iproj)
199 !    call wrtout(ab_out,  message,'COLL')
200      call wrtout(std_out,  message,'COLL')
201    end if
202    call ps_Projector_Get(psxml, idx_so(iproj), l=il)
203    psphead%nprojso(il) = psphead%nprojso(il) + 1
204  end do
205  if (iwrite == 1) then
206    write (message,'(a,5I5)') '- psxml2ab: nproj ', psphead%nproj
207 !  call wrtout(ab_out,  message,'COLL')
208    call wrtout(std_out,  message,'COLL')
209    write (message,'(a,5I5)') '- psxml2ab: nprojso ', psphead%nprojso
210 !  call wrtout(ab_out,  message,'COLL')
211    call wrtout(std_out,  message,'COLL')
212  end if
213 
214  if( has_nlcc) then
215    psphead%xccc  = n1xccc_default
216  else
217    psphead%xccc  = 0
218  end if
219 
220  psphead%pspso = 0
221  if (sum(abs(psphead%nprojso(:))) > 0) psphead%pspso = 2
222 
223 
224  psphead%lmax = 0
225  if (iwrite == 1) then
226    write (message,'(a,I5)') '- psxml2ab: ps_Number_of_Projectors not relativistic ',&
227 &        nprojs
228    call wrtout(ab_out,  message,'COLL')
229    call wrtout(std_out,  message,'COLL')
230    write (message,'(a,I5)') '- psxml2ab: ps_Number_of_Projectors scalar relativistic ', &
231 &        nprojsr
232    call wrtout(ab_out,  message,'COLL')
233    call wrtout(std_out,  message,'COLL')
234  end if
235  ll_previous=-1
236  do iproj = 1, max(nprojs, nprojsr)
237    call ps_Projector_Get(psxml, idx_sr(iproj), l=ll)
238    if (iwrite == 1) then
239      if(ll/=ll_previous)then
240        write (message,'(a,I5)') '- psxml2ab: ps_Projector_L ', ll
241        call wrtout(ab_out,  message,'COLL')
242        call wrtout(std_out,  message,'COLL')
243        ll_previous=ll
244      endif
245      call ps_Projector_Get(psxml, idx_sr(iproj), ekb=ekb)
246      write (message,'(a,E20.10)') '- psxml2ab: ps_Projector_Ekb ', ekb
247      call wrtout(ab_out,  message,'COLL')
248      call wrtout(std_out,  message,'COLL')
249    end if
250    psphead%lmax = max( psphead%lmax, ll)
251  end do
252 
253  if (iwrite == 1) then
254    write (message,'(a,I5)') '- psxml2ab: ps_Number_of_Projectors SOC ', nprojso
255    call wrtout(ab_out,  message,'COLL')
256    call wrtout(std_out,  message,'COLL')
257  end if
258  ll_previous=-1
259  do iproj = 1, nprojso
260    call ps_Projector_Get(psxml, idx_so(iproj), l=ll)
261    if (iwrite == 1) then
262      if(ll/=ll_previous)then
263        write (message,'(a,I5)') '- psxml2ab: ps_Projector_L ', ll
264        call wrtout(ab_out,  message,'COLL')
265        call wrtout(std_out,  message,'COLL')
266        ll_previous=ll
267      endif
268      call ps_Projector_Get(psxml, idx_so(iproj), ekb=ekb)
269      write (message,'(a,E20.10)') '- psxml2ab: ps_Projector_Ekb ', ekb
270      call wrtout(ab_out,  message,'COLL')
271      call wrtout(std_out,  message,'COLL')
272    end if
273    psphead%lmax = max( psphead%lmax, ll)
274  end do
275 
276  lmm     = max( nprojs, nprojsr)
277  lmm     = max( lmm, nprojso)
278 
279  write(std_out,'(a,f5.1,a,i4,a,i4)' ) '- psxml2ab:  read the values zionpsp=',&
280 & psphead%zionpsp,' , pspcod=',psphead%pspcod,' , lmax=',psphead%lmax
281 
282  if ( iwrite .eq. 1 ) then
283 
284    write(message,'(a,a)') &
285 &   '- psxml2ab: Atomic Label:                      ', &
286 &  trim(label)
287    call wrtout(ab_out,  message,'COLL')
288    call wrtout(std_out,  message,'COLL')
289 
290    write(message,'(a,f12.5)') &
291 &   '- psxml2ab: Atomic Number:                     ', &
292 &   psphead%znuclpsp
293    call wrtout(ab_out,  message,'COLL')
294    call wrtout(std_out,  message,'COLL')
295 
296    write(message,'(a,f12.5)') &
297 &   '- psxml2ab: Valence charge:                    ', &
298 &   psphead%zionpsp
299    call wrtout(ab_out,  message,'COLL')
300    call wrtout(std_out,  message,'COLL')
301 
302    write(message,'(a,a)') &
303 &   '- psxml2ab: Pseudopotential generator code :    ', &
304 &   creator
305    call wrtout(ab_out,  message,'COLL')
306    call wrtout(std_out,  message,'COLL')
307 
308    write(message,'(a,i8)') &
309 &   '- psxml2ab: Date of pseudopotential generation: ', &
310 &   psphead%pspdat
311    call wrtout(ab_out,  message,'COLL')
312    call wrtout(std_out,  message,'COLL')
313 
314    write(message,'(a,a)') &
315 &   '- psxml2ab: Pseudopotential flavor:             ', &
316 &   trim(flavor)
317    call wrtout(ab_out,  message,'COLL')
318    call wrtout(std_out,  message,'COLL')
319 
320    do ii = 1, nxc
321      call ps_LibxcFunctional_Get(psxml, ii, name=xc_name, code=dd)
322      write(message,'(a,I4,2a," (",I4,")")') &
323 &     '- psxml2ab: Exchange-correlation functional ',ii,' :    ', &
324 &     trim(xc_name), dd
325      call wrtout(ab_out,  message,'COLL')
326      call wrtout(std_out,  message,'COLL')
327    end do
328 
329    write(message,'(2a)') &
330 &   '- psxml2ab: Relativistically generated pseudopotential (not necessarily SOC!):   ', &
331 &   trim(relat)
332    call wrtout(ab_out,  message,'COLL')
333    call wrtout(std_out,  message,'COLL')
334 
335    write(message,'(2a)') &
336 &   '- psxml2ab: Spin-polarized pseudopotential:     ', &
337 &   yesno(has_spin)
338    call wrtout(ab_out,  message,'COLL')
339    call wrtout(std_out,  message,'COLL')
340 
341    select case(has_nlcc)
342      case(.true.)
343        write(message, '(a)' ) &
344 &       '- psxml2ab: XC core correction read in from XML file.'
345      case(.false.)
346        write(message, '(a)' ) &
347 &       '- psxml2ab: No core corrections.'
348        call wrtout(ab_out,message,'COLL')
349        call wrtout(std_out,  message,'COLL')
350    end select
351  end if
352 
353  ABI_FREE(zeld)
354  ABI_FREE(zelu)
355 
356  if (allocated(idx_sr)) then
357    ABI_FREE_NOCOUNT(idx_sr)
358  end if
359  if (allocated(idx_so)) then
360    ABI_FREE_NOCOUNT(idx_so)
361  end if
362 
363  call ps_destroy(psxml)
364 
365 end subroutine psxml2abheader

ABINIT/psml_die [ Functions ]

[ Top ] [ Functions ]

NAME

 psml_die

FUNCTION

  auxiliary die function for calling libPSML. Needed at link time
  allows calling software to decide how fatal the PSML die call actually is.

COPYRIGHT

 Copyright (C) 2005-2024 ABINIT group (MJV).
 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

   str = string with error message

OUTPUT

SOURCE

395 subroutine psml_die(str)
396 
397   use m_errors
398   character(len=*), intent(in) :: str
399 
400   ABI_BUG(str)
401 
402 end subroutine psml_die