TABLE OF CONTENTS
etsf_io_file_contents
[ Top ] [ etsf_io_file ] [ Methods ]
NAME
etsf_io_file_contents
FUNCTION
This is a high level routine to inquire a file and get the specifications it matches. It is usefull to know if the file is a valid density file...
INPUTS
- file_name = a list of path where input files can be found.
OUTPUT
- read_flags = a serie of flags the file matches. These flags are defined in the module etsf_io (see ETSF_IO_VALIDITY_FLAGS). It is an addition of all matching flags.
- errors = an array of size etsf_nspecs_data. For each flag that missed, it gives the error why it missed.
- lstat = return .false. if something make the file unreadable.
- error_data <type(etsf_io_low_error)> = contains the details of the error is @lstat is false.
SOURCE
subroutine etsf_io_file_contents(read_flags, errors, file_name, lstat, error_data) integer, intent(out) :: read_flags type(etsf_io_low_error), intent(out) :: errors(etsf_nspecs_data) character(len = *), intent(in) :: file_name logical, intent(out) :: lstat type(etsf_io_low_error), intent(out) :: error_data character(len = *), parameter :: me = "etsf_io_file_contents" integer :: ncid read_flags = etsf_specs_none call etsf_io_low_open_read(ncid, trim(file_name), lstat, & & error_data = error_data) if (.not. lstat) then call etsf_io_low_error_update(error_data, me) return end if call etsf_io_file_check_dielectric_function_data(ncid, lstat, errors(1)) if (lstat) then read_flags = read_flags + etsf_dielectric_function_data else call etsf_io_low_error_update(errors(1), me) end if call etsf_io_file_check_wavefunctions_data(ncid, lstat, errors(2)) if (lstat) then read_flags = read_flags + etsf_wavefunctions_data else call etsf_io_low_error_update(errors(2), me) end if call etsf_io_file_check_scalar_field_data(ncid, lstat, errors(3)) if (lstat) then read_flags = read_flags + etsf_scalar_field_data else call etsf_io_low_error_update(errors(3), me) end if call etsf_io_file_check_crystallographic_data(ncid, lstat, errors(4)) if (lstat) then read_flags = read_flags + etsf_crystallographic_data else call etsf_io_low_error_update(errors(4), me) end if ! We close the file after the definitions. call etsf_io_low_close(ncid, lstat, error_data = error_data) if (.not. lstat) then call etsf_io_low_error_update(error_data, me) return end if end subroutine etsf_io_file_contents