Cleanup of error codes in xmltools

NOTICE: if a tag is not found, ierr=1 is returned (and no longer ierr=-1)
This commit is contained in:
Paolo Giannozzi 2021-08-13 14:54:54 +02:00
parent e863308b10
commit d76276151c
3 changed files with 42 additions and 34 deletions

View File

@ -380,7 +380,7 @@ MODULE io_dyn_mat
!
IF (PRESENT(epsil)) THEN
CALL xmlr_opentag("DIELECTRIC_PROPERTIES", ierr)
IF (ierr == -1) THEN
IF (ierr == 1) THEN
IF (PRESENT(lrigid)) lrigid = .false.
IF (PRESENT(lraman)) lraman = .false.
epsil = 0.0_dp

View File

@ -45,7 +45,7 @@ CONTAINS
call xmlr_opentag ( 'qe_pp:pseudo', IERR = ierr )
if ( ierr == 0 ) then
v2 =.false.
else if ( ierr == -1 ) then
else if ( ierr == 1 ) then
rewind (iun)
call xmlr_opentag ( 'UPF', IERR = ierr )
if ( ierr == 0 ) then
@ -634,7 +634,7 @@ CONTAINS
! existing PP files may have pp_relbeta first, pp_relwfc later,
! but also the other way round - check that everything was right
!
if ( ierr == -1 .or. ierr > 0 ) then
if ( ierr > 0 ) then
ierr = -81
return
end if
@ -735,7 +735,8 @@ CONTAINS
CALL get_attr ('gipaw_data_format', upf%gipaw_data_format )
IF ( v2 ) THEN
CALL xmlr_opentag( 'PP_GIPAW_CORE_ORBITALS', IERR=ierr )
! ierr=-2 keeps track of case <PP_GIPAW_CORE_ORBITALS ... />
! ierr= 0 for case <PP_GIPAW_CORE_ORBITALS>...</PP_GIPAW_CORE_ORBITALS>
! ierr=-1 for case <PP_GIPAW_CORE_ORBITALS ... />
CALL get_attr ('number_of_core_orbitals', upf%gipaw_ncore_orbitals)
ELSE
CALL xmlr_readtag ('number_of_core_orbitals', upf%gipaw_ncore_orbitals)

View File

@ -70,6 +70,21 @@ MODULE xmltools
! utility functions
PUBLIC :: xml_protect, i2c, l2c, r2c
!
! Error codes returned by xmlr_opentag / xml_readtag:
! -1 tag with no value (e.g. <tag attr="val"/>) found (no error)
! 0 tag found and read (no error)
! 1 tag not found
! 2 error parsing file
! 3 line too long
! 4 too many levels of tags
!
! Error codes returned by xmlw_opentag / xml_writetag:
! 0 tag open and/or written (no error)
! 1 cannot write to unit "xmlunit"
! 2 tag name too long
! 3 wrong number of values for attributes
! 4 too many levels of tag
!
INTERFACE xmlr_readtag
MODULE PROCEDURE readtag_c, readtag_r, readtag_l, readtag_i, &
readtag_iv, readtag_rv, readtag_rm, readtag_rt, &
@ -296,13 +311,8 @@ CONTAINS
! name required, character: tag name
! On output: the tag is left open, ready for addition of data -
! the tag must be subsequently closed with close_xml_tag
! If ierr is present, the following value is returned:
! ierr = 0 normal execution
! ierr = 1 cannot write to unit "xmlunit"
! ierr = 2 tag name too long
! ierr = 3 too many tag levels
! ierr =10 wrong number of values for attributes
! If absent, the above error messages are printed.
! If ierr is present, the error code set in write_tag_and_attr is returned
! If ierr is absent, the above error code is reprinted on output
!
CHARACTER(LEN=*), INTENT(IN) :: name
INTEGER, INTENT(OUT),OPTIONAL :: ierr
@ -545,6 +555,7 @@ CONTAINS
!
CHARACTER(LEN=*), INTENT(IN) :: name
INTEGER :: ierr
! See list of error codes in the header of this file
!
LOGICAL :: have_list, have_vals
INTEGER :: i, la, lv, n1a,n2a, n1v, n2v
@ -555,7 +566,7 @@ CONTAINS
END IF
!
IF ( nlevel+1 > maxlevel ) THEN
ierr = 3
ierr = 4
RETURN
END IF
nlevel = nlevel+1
@ -574,7 +585,7 @@ CONTAINS
!
! attributes (if present)
!
ierr = 10
ierr = 3
if ( allocated (attrlist) ) then
WRITE (xmlunit, "(A)", ADVANCE='no', ERR=10) attrlist
deallocate (attrlist)
@ -917,10 +928,6 @@ CONTAINS
character(len=*), intent(in) :: tag
character(len=*), intent(out):: cval
integer, intent(out), optional :: ierr
! 0: tag found and read
!-1: tag not found
! 1: error parsing file
! 2: error in arguments
!
integer :: i, j, lt, ll
character(len=1) :: endtag
@ -929,8 +936,11 @@ CONTAINS
!
cval = ''
if ( eot < 0 ) then
! print *, 'end of file reached, tag not found'
if ( present(ierr) ) ierr =-1
if ( .not. present(ierr) ) then
print *, 'end of file reached, tag not found'
else
ierr = 1
end if
return
else if ( eot == 0 ) then
! print *, 'tag found, no value to read on line'
@ -955,8 +965,11 @@ CONTAINS
lt = len_trim(tag)
endtag = adjustl( line(j+i+1+lt:) )
if ( endtag /= '>' ) then
! print *, 'tag ',trim(tag),' not correctly closed'
if (present(ierr)) ierr = 1
if ( .not.present(ierr)) then
print *, 'tag ',trim(tag),' not correctly closed'
else
ierr = 2
endif
else
! end of tag found, read value (if any) and exit
if ( i > 1 ) cval = trim(cval) // adjustl(trim(line(j:j+i-2)))
@ -990,13 +1003,7 @@ CONTAINS
!
character(len=*), intent(in) :: tag
integer, intent(out), optional :: ierr
! 0: tag found and read
!-1: tag not found
!-2: tag with no value (e.g. <tag attr="val"/>) found
! 1: error parsing file
! 2: line too long
! 3: too many levels of tags
!
! See list of error codes in the header of this file
integer :: stat, ntry, ll, lt, i, j, j0
! stat=-1: in comment (not actually used)
! stat= 0: begin
@ -1017,7 +1024,7 @@ CONTAINS
ll = len_trim(line)
if ( ll == maxline ) then
print *, 'xmlr_opentag: severe error, line too long'
if (present(ierr)) ierr = 2
if (present(ierr)) ierr = 3
return
end if
! j is the current scan position
@ -1070,7 +1077,7 @@ CONTAINS
j0= j
else if ( line(j:j+1) == '/>' ) then
! <tag ... /> found : return
if (present(ierr)) ierr =-2
if (present(ierr)) ierr =-1
! eot = 0: tag with no value found
eot = 0
!
@ -1084,7 +1091,7 @@ CONTAINS
nlevel = nlevel+1
IF ( nlevel > maxlevel ) THEN
print *, 'xmlr_opentag: severe error, too many levels'
if (present(ierr)) ierr = 3
if (present(ierr)) ierr = 4
else
open_tags(nlevel) = trim(tag)
#if defined ( __debug )
@ -1131,7 +1138,7 @@ CONTAINS
!
10 if ( stat == 0 ) then
if ( present(ierr) ) then
ierr =-1
ierr = 1
! quick-and-dirty pseudo-fix to deal with tags not found:
! rewind and try again - will work if the desired tag is
! found above the current position (and nowhere else)
@ -1142,7 +1149,7 @@ CONTAINS
end if
else
print *, 'xmlr_opentag: severe parsing error'
if ( present(ierr) ) ierr = 1
if ( present(ierr) ) ierr = 2
end if
!
end subroutine xmlr_opentag
@ -1178,7 +1185,7 @@ CONTAINS
ll = len_trim(line)
if ( ll == maxline ) then
print *, 'Fatal error: line too long'
if (present(ierr)) ierr = 1
if (present(ierr)) ierr = 2
return
end if
! j is the current scan position