45 integer,
public,
parameter :: &
46 READ_COORDS_ERR = 0, &
54 integer,
public,
parameter :: &
55 XYZ_FLAGS_RESIDUE = 1, &
60 character(len=15) :: label
61 real(real64),
allocatable :: x(:)
63 real(real64) :: charge
64 character(len=3) :: residue
73 type(read_coords_atom),
allocatable :: atom(:)
75 real(real64),
allocatable :: latvec(:,:)
83 type(read_coords_info),
intent(out) :: gf
87 gf%source = read_coords_err
97 type(read_coords_info),
intent(inout) :: gf
104 safe_deallocate_a(gf%atom(ia)%x)
106 safe_deallocate_a(gf%atom)
108 safe_deallocate_a(gf%latvec)
116 character(len=*),
intent(in) :: what
117 type(read_coords_info),
intent(inout) :: gf
118 class(space_t),
intent(in) :: space
119 type(namespace_t),
intent(in) :: namespace
121 integer :: ia, ncol, iunit, jdir, int_one, nsteps, istep, step_to_use, periodic_dim
123 character(len=256) :: str
159 gf%flags = ior(gf%flags, xyz_flags_residue)
165 message(1) =
"Reading " // trim(what) //
" from " // trim(str)
168 iunit =
io_open(str, namespace, action=
'read')
194 message(1) =
"Reading " // trim(what) //
" from " // trim(str)
201 write(
message(1),
'(a,i6)')
"Invalid number of atoms ", gf%n
207 safe_allocate(gf%atom(1:gf%n))
210 safe_allocate(gf%atom(ia)%x(1:space%dim))
211 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x
237 message(1) =
"Reading " // trim(what) //
" from " // trim(str)
240 iunit =
io_open(str, namespace, action=
'read', status=
'old')
242 read(iunit,
'(a256)') str
243 if (str(1:9) ==
'ANIMSTEPS')
then
244 read(str(10:), *) nsteps
255 call parse_variable(namespace,
'XSFCoordinatesAnimStep', 1, step_to_use)
256 if (step_to_use < 1)
then
257 message(1) =
"XSFCoordinatesAnimStep must be > 0."
259 else if (step_to_use > nsteps)
then
260 write(
message(1),
'(a,i9)')
"XSFCoordinatesAnimStep must be <= available number of steps ", nsteps
263 write(
message(1),
'(a,i9,a,i9)')
'Using animation step ', step_to_use,
' out of ', nsteps
271 select case (trim(str))
283 write(
message(1),
'(3a)')
'Line in file was "', trim(str),
'" instead of CRYSTAL/SLAB/POLYMER/MOLECULE.'
286 if (periodic_dim /= space%periodic_dim)
then
287 message(1) =
"Periodicity in XSF input is incompatible with the value of PeriodicDimensions."
291 do istep = 1, step_to_use - 1
293 do jdir = 1, space%dim
297 read(iunit, *) gf%n, int_one
303 read(iunit,
'(a256)') str
304 if (str(1:7) /=
'PRIMVEC')
then
305 write(
message(1),
'(3a)')
'Line in file was "', trim(str),
'" instead of "PRIMVEC".'
309 read(str(8:), *) istep
310 if (istep /= step_to_use)
then
311 write(
message(1),
'(a,i9,a,i9)')
'Found PRIMVEC index ', istep,
' instead of ', step_to_use
316 safe_allocate(gf%latvec(1:space%dim, 1:space%dim))
318 do jdir = 1, space%dim
319 read(iunit, *) gf%latvec(1:space%dim, jdir)
323 read(iunit,
'(a256)') str
324 if (str(1:9) /=
'PRIMCOORD')
then
325 write(
message(1),
'(3a)')
'Line in file was "', trim(str),
'" instead of "PRIMCOORD".'
329 read(str(10:), *) istep
330 if (istep /= step_to_use)
then
331 write(
message(1),
'(a,i9,a,i9)')
'Found PRIMCOORD index ', istep,
' instead of ', step_to_use
336 read(iunit, *) gf%n, int_one
338 write(
message(1),
'(a,i6)')
"Invalid number of atoms ", gf%n
341 if (int_one /= 1)
then
342 write(
message(1),
'(a,i6,a)')
'Number in file was ', int_one,
' instead of 1.'
345 safe_allocate(gf%atom(1:gf%n))
349 safe_allocate(gf%atom(ia)%x(1:space%dim))
350 read(iunit,*) gf%atom(ia)%label, gf%atom(ia)%x(1:space%dim)
380 if (
parse_block(namespace, trim(what), blk) == 0)
then
388 message(1) =
"Reading " // trim(what) //
" from " // trim(what) //
" block"
391 safe_allocate(gf%atom(1:gf%n))
395 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2))
then
396 write(
message(1),
'(3a,i2)')
'Error in block ', what,
' line #', ia
401 safe_allocate(gf%atom(ia)%x(1:space%dim))
402 do jdir = 1, space%dim
405 if (ncol == space%dim + 2)
then
408 gf%atom(ia)%move = .
true.
431 if (trim(what) ==
'Coordinates')
then
432 if (
parse_block(namespace,
'Reduced'//trim(what), blk) == 0)
then
440 message(1) =
"Reading " // trim(what) //
" from Reduced" // trim(what) //
" block"
443 safe_allocate(gf%atom(1:gf%n))
447 if ((ncol < space%dim + 1) .or. (ncol > space%dim + 2))
then
448 write(
message(1),
'(3a,i2)')
'Error in block ', what,
' line #', ia
453 safe_allocate(gf%atom(ia)%x(1:space%dim))
454 do jdir = 1, space%dim
457 if (ncol == space%dim + 2)
then
460 gf%atom(ia)%move = .
true.
485 logical,
intent(inout) :: done
492 message(1) =
'Multiple definitions of '//trim(what)//
' in the input file.'
504 character(len=*),
intent(in) :: what
505 integer,
intent(in) :: iunit
508 character(len=80) :: record
509 character(len=6) :: record_name
519 read(iunit,
'(a80)', err=990,
end=990) record
520 read(record,
'(a6)') record_name
521 if (trim(record_name) ==
'ATOM' .or. trim(record_name) ==
'HETATM')
then
527 safe_allocate(gf%atom(1:gf%n))
533 read(iunit,
'(a80)', err=991,
end=991) record
534 read(record,
'(a6)') record_name
535 if (trim(record_name) ==
'ATOM' .or. trim(record_name) ==
'HETATM')
then
536 read(record,
'(12x,a4,1x,a3)') gf%atom(na)%label, gf%atom(na)%residue
537 call str_trim(gf%atom(na)%label)
538 gf%atom(na)%label = gf%atom(na)%label(1:1)
539 call str_trim(gf%atom(na)%residue)
541 read(record,
'(30x,3f8.3)') gf%atom(na)%x(1:3)
545 gf%atom(na)%x = units_to_atomic(unit_angstrom, gf%atom(na)%x)
real(real64), parameter, public m_zero
subroutine, public io_close(iunit, grp)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
type(namespace_t), public global_namespace
logical function, public parse_is_defined(namespace, name)
subroutine, public parse_block_logical(blk, l, c, res)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
integer, parameter, public read_coords_xsf
integer, parameter, public xyz_flags_charge
integer, parameter, public read_coords_reduced
subroutine, public read_coords_init(gf)
integer, parameter, public xyz_flags_move
subroutine, public read_coords_end(gf)
integer, parameter, public read_coords_xyz
subroutine read_coords_read_pdb(what, iunit, gf)
subroutine, public read_coords_read(what, gf, space, namespace)
integer, parameter, public read_coords_pdb
integer, parameter, public read_coords_inp
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine check_duplicated(done)