TABLE OF CONTENTS
etsf_io_wavedata_put
NAME
etsf_io_wavedata_put
FUNCTION
Write data related to the given group in an opened ETSF file (it must be in write mode, use etsf_io_low_set_write_mode() to change it). Only associated pointers of argument @folder will be accessed. If any errors occurs it returns with @lstat = .false..
INPUTS
- ncid = integer returned by an 'open' NetCDF call. The file can be either in define or write mode. This status can be changed by the call.
- folder <type(etsf_wavedata)> = an allocated structure with pointers on allocated areas in memory. These areas will be read or written if the pointer is associated, if not, the variable will be ignored. It is possible to access to specific dimensions of a variable using the <short_var_name>__kpoint_access or <short_var_name>__spin_access of this @folder structure. The <short_var_name>__number_of_<something> can also been set if only a subpart in one dimension should be accessed (this is possible when the specifications have been declared with a max_number_of_<something>.
OUTPUT
- lstat = return .true. if all the actions succeed, if not the status of the file is undefined.
- error_data <type(etsf_io_low_error)> = contains the details of the error is @lstat is false.
NOTES
This file has been automatically generated by the autogen_subroutines.py script. Any change you would bring to it will systematically be overwritten.
SOURCE
subroutine etsf_io_wavedata_put(ncid, folder, lstat, error_data) !Arguments ------------------------------------ integer, intent(intent) :: ncid type(etsf_wavedata), intent(intent) :: folder logical, intent(out) :: lstat type(etsf_io_low_error), intent(out) :: error_data !Local variables------------------------------- character(len = *), parameter :: my_name = 'etsf_io_wavedata_put' integer,allocatable :: varid(:) integer,allocatable :: start(:) integer,allocatable :: count(:) integer :: len character(etsf_charlen) :: flag ! ************************************************************************* !DEBUG !write (*,*) 'etsf_io_wavedata_put : enter' !ENDDEBUG allocate(varid(4)) ! Begin by putting the file in write mode. call etsf_io_low_set_write_mode(ncid, lstat, error_data = error_data) if (.not. lstat) return if (associated(folder%basis_set)) then call etsf_io_low_write_var(ncid, "basis_set", & & folder%basis_set, etsf_charlen, & & lstat, ncvarid = varid(1), & & error_data = error_data) if (.not. lstat) return end if if (associated(folder%kinetic_energy_cutoff)) then call etsf_io_low_write_var(ncid, "kinetic_energy_cutoff", & & folder%kinetic_energy_cutoff, & & lstat, ncvarid = varid(2), & & error_data = error_data) if (.not. lstat) return end if if (associated(folder%number_of_coefficients)) then call etsf_io_low_write_var(ncid, "number_of_coefficients", & & folder%number_of_coefficients, & & lstat, ncvarid = varid(3), & & error_data = error_data) if (.not. lstat) return end if if (etsf_io_low_var_associated(folder%reduced_coordinates_of_plane_waves)) then ! Handle the k_dependent attribute. call etsf_io_low_read_att(ncid, "reduced_coordinates_of_plane_waves", & & "k_dependent", & & etsf_charlen, flag, & & lstat, error_data = error_data) if (.not. lstat) return if (flag(1:2) == "no") then allocate(start(2), count(2)) else allocate(start(3), count(3)) end if start(:) = 1 count(:) = 0 if (flag(1:3) == "yes" .and. & & folder%red_coord_pw__kpoint_access /= etsf_no_sub_access) then start(3) = folder%red_coord_pw__kpoint_access count(3) = 1 end if count(2) = folder%red_coord_pw__number_of_coefficients call etsf_io_low_write_var(ncid, "reduced_coordinates_of_plane_waves", & & folder%reduced_coordinates_of_plane_waves, & & lstat, ncvarid = varid(4), & & error_data = error_data, start = start, count = count) deallocate(start, count) if (.not. lstat) return end if deallocate(varid) !DEBUG !write (*,*) 'etsf_io_wavedata_put : exit' !ENDDEBUG end subroutine etsf_io_wavedata_put