TABLE OF CONTENTS
- ABINIT/m_crystal_io
- m_crystal_io/crystal_from_hdr
- m_crystal_io/crystal_ncwrite
- m_crystal_io/crystal_ncwrite_path
ABINIT/m_crystal_io [ Modules ]
NAME
m_crystal_io
FUNCTION
Module containing the methods used to do IO on crystal objects.
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group (MG, YP, DC) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
PARENTS
CHILDREN
SOURCE
21 #if defined HAVE_CONFIG_H 22 #include "config.h" 23 #endif 24 25 #include "abi_common.h" 26 27 MODULE m_crystal_io 28 29 use defs_basis 30 use m_abicore 31 use m_xmpi 32 use m_errors 33 use m_crystal 34 use m_atomdata 35 use m_nctk 36 use iso_c_binding 37 #ifdef HAVE_NETCDF 38 use netcdf 39 #endif 40 41 use m_fstrings, only : sjoin, yesno 42 use m_io_tools, only : file_exists 43 use defs_abitypes, only : hdr_type 44 45 implicit none 46 47 private
m_crystal_io/crystal_from_hdr [ Functions ]
[ Top ] [ m_crystal_io ] [ Functions ]
NAME
crystal_from_hdr
FUNCTION
Initializes a crystal_t data type starting from the abinit header.
INPUTS
hdr<hdr_type>=the abinit header timrev ==2 => take advantage of time-reversal symmetry ==1 ==> do not use time-reversal symmetry remove_inv [optional]= if .TRUE. the inversion symmetry is removed from the set of operations even if it is present in the header
OUTPUT
cryst<crystal_t>= the data type filled with data reported in the abinit header
TODO
Add information on the use of time-reversal in the Abinit header.
PARENTS
cut3d,eph,fold2Bloch,gstate,m_ddk,m_dvdb,m_ioarr,m_iowf,m_wfd,m_wfk mlwfovlp_qp,mrgscr,setup_bse,setup_screening,setup_sigma,wfk_analyze
CHILDREN
atomdata_from_znucl,crystal_init
SOURCE
87 subroutine crystal_from_hdr(cryst,hdr,timrev,remove_inv) 88 89 90 !This section has been created automatically by the script Abilint (TD). 91 !Do not modify the following lines by hand. 92 #undef ABI_FUNC 93 #define ABI_FUNC 'crystal_from_hdr' 94 !End of the abilint section 95 96 implicit none 97 98 !Arguments ------------------------------------ 99 type(hdr_type),intent(in) :: hdr 100 type(crystal_t),intent(out) :: cryst 101 integer,intent(in) :: timrev 102 logical,optional,intent(in) :: remove_inv 103 104 !Local variables------------------------------- 105 integer :: space_group 106 logical :: rinv,use_antiferro 107 ! ********************************************************************* 108 109 rinv=.FALSE.; if (PRESENT(remove_inv)) rinv=remove_inv 110 use_antiferro = (hdr%nspden==2.and.hdr%nsppol==1) 111 112 ! Consistency check 113 ABI_CHECK(any(timrev == [1, 2]),"timrev should be in (1|2)") 114 if (use_antiferro) then 115 ABI_CHECK(ANY(hdr%symafm==-1),"Wrong nspden, nsppol, symafm.") 116 end if 117 118 space_group=0 !FIXME not known 119 120 call crystal_init(hdr%amu,cryst,space_group,hdr%natom,hdr%npsp,hdr%ntypat,hdr%nsym,hdr%rprimd,hdr%typat,hdr%xred,& 121 & hdr%zionpsp,hdr%znuclpsp,timrev,use_antiferro,rinv,hdr%title,& 122 & symrel=hdr%symrel,tnons=hdr%tnons,symafm=hdr%symafm) ! Optional 123 124 end subroutine crystal_from_hdr
m_crystal_io/crystal_ncwrite [ Functions ]
[ Top ] [ m_crystal_io ] [ Functions ]
NAME
crystal_ncwrite
FUNCTION
Output system geometry to a file, using the NETCDF file format and ETSF I/O. Data are taken from the crystal_t object.
INPUTS
cryst<crystal_t>=Object defining the unit cell and its symmetries. ncid=NC file handle.
OUTPUT
Only writing
NOTES
Alchemy not treated, since crystal should be initialized at the beginning of the run.
PARENTS
anaddb,eig2tot,exc_spectra,dfpt_looppert,m_haydock,m_phonons,m_shirley outscfcv,sigma
CHILDREN
atomdata_from_znucl
SOURCE
156 integer function crystal_ncwrite(cryst, ncid) result(ncerr) 157 158 159 !This section has been created automatically by the script Abilint (TD). 160 !Do not modify the following lines by hand. 161 #undef ABI_FUNC 162 #define ABI_FUNC 'crystal_ncwrite' 163 !End of the abilint section 164 165 implicit none 166 167 !Arguments ------------------------------------ 168 !scalars 169 integer,intent(in) :: ncid 170 type(crystal_t),intent(in) :: cryst 171 172 #ifdef HAVE_NETCDF 173 !Local variables------------------------------- 174 !scalars 175 integer :: itypat 176 character(len=500) :: msg 177 character(len=etsfio_charlen) :: symmorphic 178 type(atomdata_t) :: atom 179 !arrays 180 character(len=2) :: symbols(cryst%ntypat) 181 character(len=80) :: psp_desc(cryst%ntypat),symbols_long(cryst%ntypat) 182 183 ! ************************************************************************* 184 185 ! @crystal_t 186 187 ! TODO alchemy not treated correctly 188 if (isalchemical(cryst)) then 189 write(msg,"(3a)")& 190 "Alchemical crystals are not fully supported by the netcdf format",ch10,& 191 "Important parameters (e.g. znucl, symbols) are not written with the correct value" 192 MSG_WARNING(msg) 193 end if 194 195 symmorphic = yesno(isymmorphic(cryst)) 196 197 ! Define dimensions. 198 ncerr = nctk_def_dims(ncid, [ & 199 nctkdim_t("complex", 2), nctkdim_t("symbol_length", 2),& 200 nctkdim_t("character_string_length", 80), nctkdim_t("number_of_cartesian_directions", 3),& 201 nctkdim_t("number_of_reduced_dimensions", 3), nctkdim_t("number_of_vectors", 3),& 202 nctkdim_t("number_of_atoms", cryst%natom), nctkdim_t("number_of_atom_species", cryst%ntypat),& 203 nctkdim_t("number_of_symmetry_operations", cryst%nsym)], defmode=.True.) 204 NCF_CHECK(ncerr) 205 206 ! Define variables 207 NCF_CHECK(nctk_def_iscalars(ncid, [character(len=nctk_slen) :: "space_group"])) 208 209 ncerr = nctk_def_arrays(ncid, [ & 210 ! Atomic structure and symmetry operations 211 nctkarr_t("primitive_vectors", "dp", "number_of_cartesian_directions, number_of_vectors"), & 212 nctkarr_t("reduced_symmetry_matrices", "int", & 213 "number_of_reduced_dimensions, number_of_reduced_dimensions, number_of_symmetry_operations"), & 214 nctkarr_t("reduced_symmetry_translations", "dp", "number_of_reduced_dimensions, number_of_symmetry_operations"), & 215 nctkarr_t("atom_species", "int", "number_of_atoms"), & 216 nctkarr_t("reduced_atom_positions", "dp", "number_of_reduced_dimensions, number_of_atoms"), & 217 nctkarr_t("atomic_numbers", "dp", "number_of_atom_species"), & 218 nctkarr_t("atom_species_names", "char", "character_string_length, number_of_atom_species"), & 219 nctkarr_t("chemical_symbols", "char", "symbol_length, number_of_atom_species"), & 220 nctkarr_t('atomic_mass_units', "dp", "number_of_atom_species"), & 221 ! Atomic information. 222 nctkarr_t("valence_charges", "dp", "number_of_atom_species"), & ! NB: This variable is not written if alchemical 223 nctkarr_t("pseudopotential_types", "char", "character_string_length, number_of_atom_species") & 224 ]) 225 NCF_CHECK(ncerr) 226 227 ! Some variables require the "symmorphic" attribute. 228 NCF_CHECK(nf90_put_att(ncid, vid("reduced_symmetry_matrices"), "symmorphic", symmorphic)) 229 NCF_CHECK(nf90_put_att(ncid, vid("reduced_symmetry_translations"), "symmorphic", symmorphic)) 230 231 ! At this point we have an ETSF-compliant file. Add additional data for internal use in abinit. 232 ! TODO add spinat. 233 ncerr = nctk_def_arrays(ncid, nctkarr_t('symafm', "int", "number_of_symmetry_operations")) 234 NCF_CHECK(ncerr) 235 236 ! Set-up atomic symbols. 237 do itypat=1,cryst%ntypat 238 call atomdata_from_znucl(atom,cryst%znucl(itypat)) 239 symbols(itypat) = atom%symbol 240 write(symbols_long(itypat),'(a2,a78)') symbols(itypat),REPEAT(CHAR(0),78) 241 write(psp_desc(itypat),'(2a)') & 242 & cryst%title(itypat)(1:MIN(80,LEN_TRIM(cryst%title(itypat)))),REPEAT(CHAR(0),MAX(0,80-LEN_TRIM(cryst%title(itypat)))) 243 end do 244 245 ! Write data. 246 NCF_CHECK(nctk_set_datamode(ncid)) 247 NCF_CHECK(nf90_put_var(ncid, vid("space_group"), cryst%space_group)) 248 NCF_CHECK(nf90_put_var(ncid, vid("primitive_vectors"), cryst%rprimd)) 249 NCF_CHECK(nf90_put_var(ncid, vid("reduced_symmetry_matrices"), cryst%symrel)) 250 NCF_CHECK(nf90_put_var(ncid, vid("reduced_symmetry_translations"), cryst%tnons)) 251 NCF_CHECK(nf90_put_var(ncid, vid("atom_species"), cryst%typat)) 252 NCF_CHECK(nf90_put_var(ncid, vid("reduced_atom_positions"), cryst%xred)) 253 NCF_CHECK(nf90_put_var(ncid, vid("atomic_numbers"), cryst%znucl(1:cryst%ntypat))) 254 NCF_CHECK(nf90_put_var(ncid, vid("atom_species_names"), symbols_long)) 255 NCF_CHECK(nf90_put_var(ncid, vid("chemical_symbols"), symbols)) 256 NCF_CHECK(nf90_put_var(ncid, vid('atomic_mass_units'), cryst%amu)) 257 NCF_CHECK(nf90_put_var(ncid, vid("pseudopotential_types"), psp_desc)) 258 if (cryst%npsp == cryst%ntypat) then 259 NCF_CHECK(nf90_put_var(ncid, vid("valence_charges"), cryst%zion)) 260 end if 261 262 NCF_CHECK(nf90_put_var(ncid, vid("symafm"), cryst%symafm)) 263 264 #else 265 MSG_ERROR("netcdf library not available") 266 #endif 267 268 contains 269 integer function vid(vname) 270 271 272 !This section has been created automatically by the script Abilint (TD). 273 !Do not modify the following lines by hand. 274 #undef ABI_FUNC 275 #define ABI_FUNC 'vid' 276 !End of the abilint section 277 278 character(len=*),intent(in) :: vname 279 vid = nctk_idname(ncid, vname) 280 end function vid 281 282 end function crystal_ncwrite
m_crystal_io/crystal_ncwrite_path [ Functions ]
[ Top ] [ m_crystal_io ] [ Functions ]
NAME
crystal_ncwrite_path
FUNCTION
Output system geometry to a file, using the NETCDF file format and ETSF I/O.
INPUTS
crystal<crystal_t>=Object defining the unit cell and its symmetries. path=filename
OUTPUT
Only writing
PARENTS
CHILDREN
SOURCE
307 integer function crystal_ncwrite_path(crystal, path) result(ncerr) 308 309 310 !This section has been created automatically by the script Abilint (TD). 311 !Do not modify the following lines by hand. 312 #undef ABI_FUNC 313 #define ABI_FUNC 'crystal_ncwrite_path' 314 !End of the abilint section 315 316 implicit none 317 318 !Arguments ------------------------------------ 319 !scalars 320 character(len=*),intent(in) :: path 321 type(crystal_t),intent(in) :: crystal 322 323 #ifdef HAVE_NETCDF 324 !Local variables------------------------------- 325 !scalars 326 integer :: ncid 327 328 ! ************************************************************************* 329 330 ncerr = nf90_noerr 331 if (file_exists(path)) then 332 NCF_CHECK(nctk_open_modify(ncid, path, xmpi_comm_self)) 333 else 334 ncerr = nctk_open_create(ncid, path, xmpi_comm_self) 335 NCF_CHECK_MSG(ncerr, sjoin("creating:", path)) 336 end if 337 338 NCF_CHECK(crystal_ncwrite(crystal, ncid)) 339 NCF_CHECK(nf90_close(ncid)) 340 #endif 341 342 end function crystal_ncwrite_path