46#if defined(HAVE_NETCDF)
82#if defined(HAVE_NETCDF)
89 integer,
parameter,
private :: &
94 character(len=3),
parameter :: &
95 index2label(3) = (/
're ',
'im ',
'abs' /)
119 what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
120 type(namespace_t),
intent(in) :: namespace
121 class(space_t),
intent(in) :: space
122 logical,
intent(inout) :: what(MAX_OUTPUT_TYPES)
123 integer(int64),
intent(out) :: how(0:MAX_OUTPUT_TYPES)
124 integer,
intent(out) :: output_interval(0:MAX_OUTPUT_TYPES)
125 character(len=*),
optional,
intent(in) :: what_tag_in
126 character(len=*),
optional,
intent(in) :: how_tag_in
127 character(len=*),
optional,
intent(in) :: output_interval_tag_in
128 logical,
optional,
intent(in) :: ignore_error
131 integer :: ncols, nrows, iout, column_index, what_i
132 integer(int64) :: what_no_how(12)
133 character(len=80) :: what_tag
134 character(len=80) :: how_tag
135 character(len=80) :: output_interval_tag
136 character(len=80) :: output_column_marker
140 output_interval = 0_8
144 output_interval_tag =
optional_default(output_interval_tag_in,
'OutputInterval')
479 option__output__matrix_elements, &
480 option__output__berkeleygw, &
481 option__output__dos, &
482 option__output__tpa, &
483 option__output__mmb_den, &
484 option__output__j_flow, &
485 option__output__occ_matrices, &
486 option__output__effectiveu, &
487 option__output__magnetization, &
488 option__output__kanamoriu, &
489 option__output__stress, &
490 option__output__jdos &
493 if (
parse_block(namespace, what_tag, blk) == 0)
then
495 do iout = 0, nrows - 1
511 what(what_i) = .
true.
512 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
513 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
514 .or. (what_tag /=
'Output'))
then
522 else if (ncols == 2)
then
535 what(what_i) = .
true.
536 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
537 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
538 .or. (what_tag /=
'Output'))
then
540 if (how(what_i) == 0)
call parse_variable(namespace, how_tag, 0, how(what_i))
567 what(what_i) = .
true.
568 do column_index = 0, 1
570 if (output_column_marker ==
'output_interval')
then
572 else if (output_column_marker ==
'output_format')
then
573 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
574 .or. (what_tag /=
'Output'))
then
580 else if (len_trim(output_column_marker) /= 0)
then
585 if (output_interval(what_i) == 0)
then
586 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
588 if (how(what_i) == 0)
then
589 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
590 .or. (what_tag /=
'Output'))
then
609 call parse_variable(namespace, output_interval_tag, 50, output_interval(0))
613 do what_i = lbound(what, 1), ubound(what, 1)
614 if (what_tag ==
'Output')
then
615 if (what(what_i) .and. (.not. any(what_no_how == what_i)))
then
620 if (how(what_i) == 0 .and. .not.
optional_default(ignore_error, .false.))
then
621 write(
message(1),
'(a)')
'Must specify output method with variable OutputFormat.'
626 if (space%dim == 1)
then
627 if (
bitand(how(what_i), option__outputformat__axis_y) /= 0)
then
628 message(1) =
"OutputFormat = axis_y not available with Dimensions = 1."
631 if (
bitand(how(what_i), option__outputformat__plane_z) /= 0)
then
632 message(1) =
"OutputFormat = plane_z not available with Dimensions = 1."
635 if (
bitand(how(what_i), option__outputformat__xcrysden) /= 0)
then
636 message(1) =
"OutputFormat = xcrysden not available with Dimensions = 1."
641 if (space%dim <= 2)
then
642 if (
bitand(how(what_i), option__outputformat__axis_z) /= 0)
then
643 message(1) =
"OutputFormat = axis_z not available with Dimensions <= 2."
646 if (
bitand(how(what_i), option__outputformat__plane_x) /= 0)
then
647 message(1) =
"OutputFormat = plane_x not available with Dimensions <= 2."
650 if (
bitand(how(what_i), option__outputformat__plane_y) /= 0)
then
651 message(1) =
"OutputFormat = plane_y not available with Dimensions <= 2."
654 if (
bitand(how(what_i), option__outputformat__integrate_xy) /= 0)
then
655 message(1) =
"OutputFormat = integrate_xy not available with Dimensions <= 2."
658 if (
bitand(how(what_i), option__outputformat__integrate_xz) /= 0)
then
659 message(1) =
"OutputFormat = integrate_xz not available with Dimensions <= 2."
662 if (
bitand(how(what_i), option__outputformat__integrate_yz) /= 0)
then
663 message(1) =
"OutputFormat = integrate_yz not available with Dimensions <= 2."
666 if (
bitand(how(what_i), option__outputformat__dx) /= 0)
then
667 message(1) =
"OutputFormat = dx not available with Dimensions <= 2."
670 if (
bitand(how(what_i), option__outputformat__cube) /= 0)
then
671 message(1) =
"OutputFormat = cube not available with Dimensions <= 2."
674 if (
bitand(how(what_i), option__outputformat__poscar) /= 0)
then
675 message(1) =
"OutputFormat = POSCAR not available with Dimensions <= 2."
680 if (space%periodic_dim<=2)
then
681 if (
bitand(how(what_i), option__outputformat__poscar) /= 0)
then
682 message(1) =
"OutputFormat = POSCAR not available with periodic_dim <= 2."
687#if !defined(HAVE_NETCDF)
688 if (
bitand(how(what_i), option__outputformat__netcdf) /= 0)
then
689 message(1) =
'Octopus was compiled without NetCDF support.'
690 message(2) =
'It is not possible to write output in NetCDF format.'
694#if !defined(HAVE_ETSF_IO)
695 if (
bitand(how(what_i), option__outputformat__etsf) /= 0)
then
696 message(1) =
'Octopus was compiled without ETSF_IO support.'
697 message(2) =
'It is not possible to write output in ETSF format.'
706 if (output_interval(what_i) < 0)
then
707 message(1) =
"OutputInterval must be >= 0."
721 character(len=*),
intent(in) :: where
726 if (index(
where,
"AxisX") /= 0) how = ior(how, option__outputformat__axis_x)
727 if (index(
where,
"AxisY") /= 0) how = ior(how, option__outputformat__axis_y)
728 if (index(
where,
"AxisZ") /= 0) how = ior(how, option__outputformat__axis_z)
729 if (index(
where,
"PlaneX") /= 0) how = ior(how, option__outputformat__plane_x)
730 if (index(
where,
"PlaneY") /= 0) how = ior(how, option__outputformat__plane_y)
731 if (index(
where,
"PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
732 if (index(
where,
"IntegrateXY") /= 0) how = ior(how, option__outputformat__integrate_xy)
733 if (index(
where,
"IntegrateXZ") /= 0) how = ior(how, option__outputformat__integrate_xz)
734 if (index(
where,
"IntegrateYZ") /= 0) how = ior(how, option__outputformat__integrate_yz)
735 if (index(
where,
"PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
736 if (index(
where,
"DX") /= 0) how = ior(how, option__outputformat__dx)
737 if (index(
where,
"XCrySDen") /= 0) how = ior(how, option__outputformat__xcrysden)
738 if (index(
where,
"Binary") /= 0) how = ior(how, option__outputformat__binary)
739 if (index(
where,
"MeshIndex") /= 0) how = ior(how, option__outputformat__mesh_index)
740 if (index(
where,
"XYZ") /= 0) how = ior(how, option__outputformat__xyz)
741#if defined(HAVE_NETCDF)
742 if (index(
where,
"NETCDF") /= 0) how = ior(how, option__outputformat__netcdf)
744 if (index(
where,
"Cube") /= 0) how = ior(how, option__outputformat__cube)
745 if (index(
where,
"VTK") /= 0) how = ior(how, option__outputformat__vtk)
756 character(len=*),
intent(in) :: fname
757 real(real64),
intent(in) :: pos(:,:)
758 character(len=*),
intent(in) :: species(:)
759 character(len=*),
optional,
intent(in) :: header
761 integer :: iunit, n_atoms, iatom, min_dim
762 real(real64) :: position(3)
766 n_atoms =
size(pos, 2)
767 assert(
size(species) == n_atoms)
769 iunit =
io_open(trim(fname)//
'.xyz', namespace, action=
'write', position=
'asis')
771 write(iunit,
'(i6)') n_atoms
772 if (
present(header))
then
773 write(iunit,
'(a)') trim(adjustl(header))
779 min_dim = min(
size(pos, 1), 3)
780 position(1:3) = 0.0_real64
783 do iatom = 1, n_atoms
784 position(1:min_dim) = pos(1:min_dim, iatom)
785 write(iunit,
'(a10, 3(1x, f11.6))') trim(adjustl(species(iatom))), &
801 character(len=*),
intent(in) :: dir
802 character(len=*),
intent(in) :: fname
803 class(
space_t),
intent(in) :: space
805 real(real64),
intent(in) :: pos(:,:)
806 type(
atom_t),
intent(in) :: atoms(:)
807 class(
box_t),
intent(in) :: box
811 character(len=:),
allocatable :: info_str, dir_fname
812 character(len=LABEL_LEN),
allocatable :: species(:)
818 info_str = trim(space%short_info()) //
'; ' // trim(box%short_info(
unit_angstrom))
819 if (space%is_periodic())
then
820 info_str = info_str //
'; ' // trim(latt%short_info(
unit_angstrom))
823 safe_allocate(species(1:
size(atoms)))
824 do ia = 1,
size(atoms)
825 species(ia) = atoms(ia)%get_label()
828 dir_fname = trim(dir) //
'/' // trim(fname)
831 safe_deallocate_a(species)
838 character(len=*),
intent(in) :: dir, fname
839 class(
space_t),
intent(in) :: space
841 real(real64),
intent(in) :: pos(:,:)
842 type(
atom_t),
intent(in) :: atoms(:)
843 class(
mesh_t),
intent(in) :: mesh
845 real(real64),
optional,
intent(in) :: total_forces(:,:)
848 real(real64),
allocatable :: forces(:,:)
855 iunit =
io_open(trim(dir)//
'/'//trim(fname)//
'.xsf', namespace, action=
'write', position=
'asis')
857 if (
present(total_forces))
then
858 safe_allocate(forces(1:space%dim, 1:
size(atoms)))
861 safe_deallocate_a(forces)
875 integer,
intent(in) :: iunit
876 class(
space_t),
intent(in) :: space
878 real(real64),
intent(in) :: pos(:,:)
879 type(
atom_t),
intent(in) :: atoms(:)
880 class(
mesh_t),
intent(in) :: mesh
881 real(real64),
optional,
intent(in) :: forces(:, :)
882 integer,
optional,
intent(in) :: index
884 integer :: idir, idir2, iatom, index_, natoms
885 character(len=7) :: index_str
886 real(real64) :: offset(space%dim)
887 real(real64) :: rlattice(space%dim, space%dim)
891 if (
present(index))
then
892 write(index_str,
'(a,i6)')
' ', index
895 write(index_str,
'(a)')
''
898 natoms =
size(pos, dim=2)
904 offset(1:space%dim) = latt%red_to_cart(spread(-
m_half, 1, space%dim))
906 do idir = space%periodic_dim + 1, space%dim
907 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
910 if (space%is_periodic())
then
911 if (index_ == 1)
then
912 select case (space%periodic_dim)
914 write(iunit,
'(a)')
'CRYSTAL'
916 write(iunit,
'(a)')
'SLAB'
918 write(iunit,
'(a)')
'POLYMER'
922 write(iunit,
'(a)')
'PRIMVEC'//trim(index_str)
925 rlattice = latt%rlattice
926 do idir = space%periodic_dim+1, space%dim
927 rlattice(:,idir) = rlattice(:,idir)*
m_two*mesh%box%bounding_box_l(idir)
930 do idir = 1, space%dim
934 write(iunit,
'(a)')
'PRIMCOORD'//trim(index_str)
935 write(iunit,
'(i10, a)') natoms,
' 1'
937 write(iunit,
'(a)')
'ATOMS'//trim(index_str)
942 write(iunit,
'(a10, 3f12.6)', advance=
'no') trim(atoms(iatom)%label), &
944 if (
present(forces))
then
945 write(iunit,
'(5x, 3f12.6)', advance=
'no') forces(:, iatom)
957 integer,
intent(in) :: iunit
958 class(
space_t),
intent(in) :: space
960 real(real64),
intent(in) :: pos(:,:)
961 type(
atom_t),
intent(in) :: atoms(:)
962 class(
mesh_t),
intent(in) :: mesh
963 real(real64),
intent(in) :: centers(:, :)
964 integer,
intent(in) :: supercell(:)
965 real(real64),
optional,
intent(in) :: extra_atom(:)
967 integer :: idir, idir2, iatom, index_
968 character(len=7) :: index_str
969 real(real64) :: offset(3)
970 integer :: irep, nreplica
974 write(index_str,
'(a)')
''
977 nreplica = product(supercell(1:space%dim))
982 offset(1:space%dim) = latt%red_to_cart(spread(-
m_half, 1, space%dim))
983 offset(1:3) = offset(1:3) + centers(1:3,1)
985 do idir = space%periodic_dim + 1, 3
986 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
989 if(space%is_periodic())
then
991 select case(space%periodic_dim)
993 write(iunit,
'(a)')
'CRYSTAL'
995 write(iunit,
'(a)')
'SLAB'
997 write(iunit,
'(a)')
'POLYMER'
1001 write(iunit,
'(a)')
'PRIMVEC'//trim(index_str)
1003 do idir = 1, space%dim
1005 latt%rlattice(idir2, idir)*supercell(idir)), idir2 = 1, space%dim)
1008 write(iunit,
'(a)')
'PRIMCOORD'//trim(index_str)
1009 if(.not.
present(extra_atom))
then
1010 write(iunit,
'(i10, a)')
size(atoms)*nreplica,
' 1'
1012 write(iunit,
'(i10, a)')
size(atoms)*nreplica+1,
' 1'
1015 write(iunit,
'(a)')
'ATOMS'//trim(index_str)
1019 do irep = 1, nreplica
1021 do iatom = 1,
size(atoms)
1022 write(iunit,
'(a10, 3f12.6)', advance=
'no') trim(atoms(iatom)%label), &
1024 - offset(idir)), idir = 1, space%dim)
1028 write(iunit,
'(a10, 3f12.6)', advance=
'no')
'X', &
1036#if defined(HAVE_NETCDF)
1038 subroutine ncdf_error(func, status, filename, namespace, ierr)
1039 character(len=*),
intent(in) ::
func
1040 integer,
intent(in) :: status
1041 character(len=*),
intent(in) :: filename
1043 integer,
intent(inout) :: ierr
1045 push_sub(ncdf_error)
1047 if (status == nf90_noerr)
then
1052 write(
message(1),
'(3a)')
"NETCDF error in function '" , trim(
func) ,
"'"
1053 write(
message(2),
'(3a)')
"(reading/writing ", trim(filename) ,
")"
1054 write(
message(3),
'(6x,a,a)')
'Error code = ', trim(nf90_strerror(status))
1059 end subroutine ncdf_error
1064 real(real64),
intent(in) :: in(:, :, :)
1065 real(real64),
intent(out) :: out(:, :, :)
1066 integer :: ix, iy, iz
1070 do ix = lbound(in, 1), ubound(in, 1)
1071 do iy = lbound(in, 2), ubound(in, 2)
1072 do iz = lbound(in, 3), ubound(in, 3)
1073 out(iz, iy, ix) = in(ix, iy, iz)
1083#include "io_function_inc.F90"
1086#include "complex.F90"
1087#include "io_function_inc.F90"
real(real64) function func(r1, rn, n, a)
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
real(real64), parameter, public m_two
real(real64), parameter, public m_half
This module implements the index, used for the mesh points.
subroutine zio_function_output_vector_bz(how, dir, fname, namespace, space, kpt, kpoints, ff, unit, ierr, grp, root)
subroutine, public write_canonicalized_xyz_file(dir, fname, space, latt, pos, atoms, box, namespace)
Write canonicalized xyz file with atom labels and positions in Angstroms.
subroutine, public zio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine, public dio_function_output_global(how, dir, fname, namespace, space, mesh, ff, unit, ierr)
subroutine, public io_function_read_what_how_when(namespace, space, what, how, output_interval, what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
subroutine, public zio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
subroutine, public write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
for format specification see: http:
subroutine zio_function_output_vector(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
subroutine write_xsf_geometry_supercell(iunit, space, latt, pos, atoms, mesh, centers, supercell, extra_atom)
for format specification see: http:
subroutine dio_function_output_supercell(how, dir, fname, mesh, space, latt, ff, centers, supercell, unit, ierr, namespace, pos, atoms, grp, root, is_global, extra_atom)
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine dio_function_output_vector(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
subroutine transpose3(in, out)
subroutine zio_function_output_supercell(how, dir, fname, mesh, space, latt, ff, centers, supercell, unit, ierr, namespace, pos, atoms, grp, root, is_global, extra_atom)
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
subroutine dio_function_output_global_bz(how, dir, fname, namespace, kpoints, ff, unit, ierr)
integer, parameter, private zoutput_kind
subroutine, public write_standard_xyz_file(namespace, fname, pos, species, header)
Write a standard xyz file with atom labels and positions (in Angstrom).
subroutine, public zio_function_output_global(how, dir, fname, namespace, space, mesh, ff, unit, ierr)
subroutine dio_function_output_vector_bz(how, dir, fname, namespace, space, kpt, kpoints, ff, unit, ierr, grp, root)
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
subroutine, public dio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
subroutine zio_function_output_global_bz(how, dir, fname, namespace, kpoints, ff, unit, ierr)
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_variable_is_block(namespace, name)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
type(mpi_grp_t), public mpi_world
Some general things and nomenclature:
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
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_t), public unit_angstrom
For XYZ files.
type(unit_system_t), public units_out
This module is intended to contain simple general-purpose utility functions and procedures.
class to tell whether a point is inside or outside
Describes mesh distribution to nodes.