TABLE OF CONTENTS
etsf_io_geometry_def
[ Top ] [ etsf_geometry ] [ Methods ]
NAME
etsf_io_geometry_def
FUNCTION
The given ETSF file must be opened and in define state (see etsf_io_low_set_define_mode() to change it). Then, all variable of the group are defined. All required dimensions must have already defined (see etsf_io_dims_def(). If some dimensions are missing, then the variable is not defined and no error are generated.
One can specify which variable may be splitted using the optional argument @split. For each associated array in this structure, variable with appropriated dimensions will use my_<something> instead of <something>.
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.
- k_dependent = (optional) use this argument to set the attribute flag 'k_dependent' to 'yes' or 'no' on variables that have it. If no variable in the group has the attribute 'k_dependent', this argument has no effect. The default value is .true. (which puts 'yes' in the file).
- flags = (optional) One can choose the variables of the group that will be defined (and disk allocated) using this flag. This is a sum of values taken from #FLAGS_VARIABLES.
- split <type(etsf_split)> = (optional) for each array associated in the type, the dimension used to declared the variables sizes will be 'my_/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_geometry_def(ncid, lstat, error_data, k_dependent, flags, & & split) !Arguments ------------------------------------ integer, intent(in) :: ncid logical, intent(out) :: lstat type(etsf_io_low_error), intent(out) :: error_data logical, optional, intent(in) :: k_dependent integer, optional, intent(in) :: flags type(etsf_split), optional, intent(in) :: split !Local variables------------------------------- character(len = *), parameter :: my_name = 'etsf_io_geometry_def' logical :: my_k_dependent integer :: my_flags type(etsf_split) :: my_split integer :: ivar type(split_dim_names) :: split_dims ! ************************************************************************* !DEBUG !write (*,*) 'etsf_io_geometry_def : enter' !ENDDEBUG ! Get values for optional arguments, set default. if (present(k_dependent)) then my_k_dependent = k_dependent else my_k_dependent = .true. end if if (present(flags)) then my_flags = flags else my_flags = etsf_geometry_all end if ! Consistency checks. if (my_flags < etsf_geometry_none .or. my_flags > etsf_geometry_all) then call etsf_io_low_error_set(error_data, ERROR_MODE_DEF, ERROR_TYPE_ARG, my_name, & & tgtname = "flags", errmess = "value out of bounds") lstat = .false. return end if if (iand(my_flags, etsf_geometry_space_group) /= 0) then call etsf_io_low_def_var(ncid, "space_group", & & etsf_io_low_integer, & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_primitive_vectors) /= 0) then call etsf_io_low_def_var(ncid, "primitive_vectors", & & etsf_io_low_double, & & (/ pad("number_of_cartesian_directions"), & & pad("number_of_vectors") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_red_sym_matrices) /= 0) then call etsf_io_low_def_var(ncid, "reduced_symmetry_matrices", & & etsf_io_low_integer, & & (/ pad("number_of_reduced_dimensions"), & & pad("number_of_reduced_dimensions"), & & pad("number_of_symmetry_operations") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if if (ivar >= 0) then ! Handle the symmorphic attribute. call etsf_io_low_write_att(ncid, ivar, & & "symmorphic", & & "yes", & & lstat, error_data = error_data) if (.not. lstat) then call etsf_io_low_error_update(error_data, my_name) return end if end if end if if (iand(my_flags, etsf_geometry_red_sym_trans) /= 0) then call etsf_io_low_def_var(ncid, "reduced_symmetry_translations", & & etsf_io_low_double, & & (/ pad("number_of_reduced_dimensions"), & & pad("number_of_symmetry_operations") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_atom_species) /= 0) then call etsf_io_low_def_var(ncid, "atom_species", & & etsf_io_low_integer, & & (/ pad("number_of_atoms") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_red_at_pos) /= 0) then call etsf_io_low_def_var(ncid, "reduced_atom_positions", & & etsf_io_low_double, & & (/ pad("number_of_reduced_dimensions"), & & pad("number_of_atoms") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_valence_charges) /= 0) then call etsf_io_low_def_var(ncid, "valence_charges", & & etsf_io_low_double, & & (/ pad("number_of_atom_species") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_atomic_numbers) /= 0) then call etsf_io_low_def_var(ncid, "atomic_numbers", & & etsf_io_low_double, & & (/ pad("number_of_atom_species") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_atom_species_names) /= 0) then call etsf_io_low_def_var(ncid, "atom_species_names", & & etsf_io_low_character, & & (/ pad("character_string_length"), & & pad("number_of_atom_species") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_chemical_symbols) /= 0) then call etsf_io_low_def_var(ncid, "chemical_symbols", & & etsf_io_low_character, & & (/ pad("symbol_length"), & & pad("number_of_atom_species") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if if (iand(my_flags, etsf_geometry_pseudo_types) /= 0) then call etsf_io_low_def_var(ncid, "pseudopotential_types", & & etsf_io_low_character, & & (/ pad("character_string_length"), & & pad("number_of_atom_species") /), & & lstat, ncvarid = ivar, error_data = error_data) ! We raise don't raise an error if a dimension is missing. if (.not. lstat .and. (error_data%access_mode_id /= ERROR_MODE_INQ .or. & & error_data%target_type_id /= ERROR_TYPE_DID)) then call etsf_io_low_error_update(error_data, my_name) return end if end if ! If we reach the end, then it should be OK. lstat = .true. !DEBUG !write (*,*) 'etsf_io_geometry_def : exit' !ENDDEBUG end subroutine etsf_io_geometry_def