46#if defined(HAVE_NETCDF)
81#if defined(HAVE_NETCDF)
88 integer,
parameter,
private :: &
93 character(len=3),
parameter :: &
94 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')
473 what_no_how = (/ option__output__matrix_elements, option__output__berkeleygw, option__output__dos, &
474 option__output__tpa, option__output__mmb_den, option__output__j_flow, option__output__occ_matrices, &
475 option__output__effectiveu, option__output__magnetization, option__output__kanamoriu, option__output__stress, option__output__jdos /)
477 if (
parse_block(namespace, what_tag, blk) == 0)
then
479 do iout = 0, nrows - 1
495 what(what_i) = .
true.
496 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
497 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
498 .or. (what_tag /=
'Output'))
then
506 else if (ncols == 2)
then
519 what(what_i) = .
true.
520 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
521 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
522 .or. (what_tag /=
'Output'))
then
524 if (how(what_i) == 0)
call parse_variable(namespace, how_tag, 0, how(what_i))
551 what(what_i) = .
true.
552 do column_index = 0, 1
554 if (output_column_marker ==
'output_interval')
then
556 else if (output_column_marker ==
'output_format')
then
557 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
558 .or. (what_tag /=
'Output'))
then
564 else if (len_trim(output_column_marker) /= 0)
then
569 if (output_interval(what_i) == 0)
then
570 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
572 if (how(what_i) == 0)
then
573 if (((what_tag ==
'Output') .and. (.not. any(what_no_how == what_i)))&
574 .or. (what_tag /=
'Output'))
then
593 call parse_variable(namespace, output_interval_tag, 50, output_interval(0))
597 do what_i = lbound(what, 1), ubound(what, 1)
598 if (what_tag ==
'Output')
then
599 if (what(what_i) .and. (.not. any(what_no_how == what_i)))
then
604 if (how(what_i) == 0 .and. .not.
optional_default(ignore_error, .false.))
then
605 write(
message(1),
'(a)')
'Must specify output method with variable OutputFormat.'
610 if (space%dim == 1)
then
611 if (
bitand(how(what_i), option__outputformat__axis_y) /= 0)
then
612 message(1) =
"OutputFormat = axis_y not available with Dimensions = 1."
615 if (
bitand(how(what_i), option__outputformat__plane_z) /= 0)
then
616 message(1) =
"OutputFormat = plane_z not available with Dimensions = 1."
619 if (
bitand(how(what_i), option__outputformat__xcrysden) /= 0)
then
620 message(1) =
"OutputFormat = xcrysden not available with Dimensions = 1."
625 if (space%dim <= 2)
then
626 if (
bitand(how(what_i), option__outputformat__axis_z) /= 0)
then
627 message(1) =
"OutputFormat = axis_z not available with Dimensions <= 2."
630 if (
bitand(how(what_i), option__outputformat__plane_x) /= 0)
then
631 message(1) =
"OutputFormat = plane_x not available with Dimensions <= 2."
634 if (
bitand(how(what_i), option__outputformat__plane_y) /= 0)
then
635 message(1) =
"OutputFormat = plane_y not available with Dimensions <= 2."
638 if (
bitand(how(what_i), option__outputformat__integrate_xy) /= 0)
then
639 message(1) =
"OutputFormat = integrate_xy not available with Dimensions <= 2."
642 if (
bitand(how(what_i), option__outputformat__integrate_xz) /= 0)
then
643 message(1) =
"OutputFormat = integrate_xz not available with Dimensions <= 2."
646 if (
bitand(how(what_i), option__outputformat__integrate_yz) /= 0)
then
647 message(1) =
"OutputFormat = integrate_yz not available with Dimensions <= 2."
650 if (
bitand(how(what_i), option__outputformat__dx) /= 0)
then
651 message(1) =
"OutputFormat = dx not available with Dimensions <= 2."
654 if (
bitand(how(what_i), option__outputformat__cube) /= 0)
then
655 message(1) =
"OutputFormat = cube not available with Dimensions <= 2."
660#if !defined(HAVE_NETCDF)
661 if (
bitand(how(what_i), option__outputformat__netcdf) /= 0)
then
662 message(1) =
'Octopus was compiled without NetCDF support.'
663 message(2) =
'It is not possible to write output in NetCDF format.'
667#if !defined(HAVE_ETSF_IO)
668 if (
bitand(how(what_i), option__outputformat__etsf) /= 0)
then
669 message(1) =
'Octopus was compiled without ETSF_IO support.'
670 message(2) =
'It is not possible to write output in ETSF format.'
679 if (output_interval(what_i) < 0)
then
680 message(1) =
"OutputInterval must be >= 0."
694 character(len=*),
intent(in) :: where
699 if (index(
where,
"AxisX") /= 0) how = ior(how, option__outputformat__axis_x)
700 if (index(
where,
"AxisY") /= 0) how = ior(how, option__outputformat__axis_y)
701 if (index(
where,
"AxisZ") /= 0) how = ior(how, option__outputformat__axis_z)
702 if (index(
where,
"PlaneX") /= 0) how = ior(how, option__outputformat__plane_x)
703 if (index(
where,
"PlaneY") /= 0) how = ior(how, option__outputformat__plane_y)
704 if (index(
where,
"PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
705 if (index(
where,
"IntegrateXY") /= 0) how = ior(how, option__outputformat__integrate_xy)
706 if (index(
where,
"IntegrateXZ") /= 0) how = ior(how, option__outputformat__integrate_xz)
707 if (index(
where,
"IntegrateYZ") /= 0) how = ior(how, option__outputformat__integrate_yz)
708 if (index(
where,
"PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
709 if (index(
where,
"DX") /= 0) how = ior(how, option__outputformat__dx)
710 if (index(
where,
"XCrySDen") /= 0) how = ior(how, option__outputformat__xcrysden)
711 if (index(
where,
"Binary") /= 0) how = ior(how, option__outputformat__binary)
712 if (index(
where,
"MeshIndex") /= 0) how = ior(how, option__outputformat__mesh_index)
713 if (index(
where,
"XYZ") /= 0) how = ior(how, option__outputformat__xyz)
714#if defined(HAVE_NETCDF)
715 if (index(
where,
"NETCDF") /= 0) how = ior(how, option__outputformat__netcdf)
717 if (index(
where,
"Cube") /= 0) how = ior(how, option__outputformat__cube)
718 if (index(
where,
"VTK") /= 0) how = ior(how, option__outputformat__vtk)
729 character(len=*),
intent(in) :: dir
730 character(len=*),
intent(in) :: fname
731 class(
space_t),
intent(in) :: space
733 real(real64),
intent(in) :: pos(:,:)
734 type(
atom_t),
intent(in) :: atoms(:)
735 class(
box_t),
intent(in) :: box
740 real(real64) :: position
746 iunit =
io_open(trim(dir)//
'/'//trim(fname)//
'.xyz', namespace, action=
'write', position=
'asis')
748 write(iunit,
'(i6)')
size(atoms)
749 write(iunit,
'(a,a,a)', advance=
'no') trim(space%short_info()),
'; ', trim(box%short_info(
unit_angstrom))
750 if (space%is_periodic())
then
751 write(iunit,
'(a,a)')
'; ', trim(latt%short_info(
unit_angstrom))
757 do iatom = 1,
size(atoms)
758 write(iunit,
'(10a)', advance=
'no') atoms(iatom)%label
760 if (idir <= space%dim)
then
761 position = pos(idir, iatom)
777 character(len=*),
intent(in) :: dir, fname
778 class(
space_t),
intent(in) :: space
780 real(real64),
intent(in) :: pos(:,:)
781 type(
atom_t),
intent(in) :: atoms(:)
782 class(
mesh_t),
intent(in) :: mesh
784 real(real64),
optional,
intent(in) :: total_forces(:,:)
787 real(real64),
allocatable :: forces(:,:)
794 iunit =
io_open(trim(dir)//
'/'//trim(fname)//
'.xsf', namespace, action=
'write', position=
'asis')
796 if (
present(total_forces))
then
797 safe_allocate(forces(1:space%dim, 1:
size(atoms)))
800 safe_deallocate_a(forces)
814 integer,
intent(in) :: iunit
815 class(
space_t),
intent(in) :: space
817 real(real64),
intent(in) :: pos(:,:)
818 type(
atom_t),
intent(in) :: atoms(:)
819 class(
mesh_t),
intent(in) :: mesh
820 real(real64),
optional,
intent(in) :: forces(:, :)
821 integer,
optional,
intent(in) :: index
823 integer :: idir, idir2, iatom, index_, natoms
824 character(len=7) :: index_str
825 real(real64) :: offset(space%dim)
826 real(real64) :: rlattice(space%dim, space%dim)
830 if (
present(index))
then
831 write(index_str,
'(a,i6)')
' ', index
834 write(index_str,
'(a)')
''
837 natoms =
size(pos, dim=2)
843 offset(1:space%dim) = latt%red_to_cart(spread(-
m_half, 1, space%dim))
845 do idir = space%periodic_dim + 1, space%dim
846 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
849 if (space%is_periodic())
then
850 if (index_ == 1)
then
851 select case (space%periodic_dim)
853 write(iunit,
'(a)')
'CRYSTAL'
855 write(iunit,
'(a)')
'SLAB'
857 write(iunit,
'(a)')
'POLYMER'
861 write(iunit,
'(a)')
'PRIMVEC'//trim(index_str)
864 rlattice = latt%rlattice
865 do idir = space%periodic_dim+1, space%dim
866 rlattice(:,idir) = rlattice(:,idir)*
m_two*mesh%box%bounding_box_l(idir)
869 do idir = 1, space%dim
873 write(iunit,
'(a)')
'PRIMCOORD'//trim(index_str)
874 write(iunit,
'(i10, a)') natoms,
' 1'
876 write(iunit,
'(a)')
'ATOMS'//trim(index_str)
881 write(iunit,
'(a10, 3f12.6)', advance=
'no') trim(atoms(iatom)%label), &
883 if (
present(forces))
then
884 write(iunit,
'(5x, 3f12.6)', advance=
'no') forces(:, iatom)
896 integer,
intent(in) :: iunit
897 class(
space_t),
intent(in) :: space
899 real(real64),
intent(in) :: pos(:,:)
900 type(
atom_t),
intent(in) :: atoms(:)
901 class(
mesh_t),
intent(in) :: mesh
902 real(real64),
intent(in) :: centers(:, :)
903 integer,
intent(in) :: supercell(:)
904 real(real64),
optional,
intent(in) :: extra_atom(:)
906 integer :: idir, idir2, iatom, index_
907 character(len=7) :: index_str
908 real(real64) :: offset(3)
909 integer :: irep, nreplica
913 write(index_str,
'(a)')
''
916 nreplica = product(supercell(1:space%dim))
921 offset(1:space%dim) = latt%red_to_cart(spread(-
m_half, 1, space%dim))
922 offset(1:3) = offset(1:3) + centers(1:3,1)
924 do idir = space%periodic_dim + 1, 3
925 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
928 if(space%is_periodic())
then
930 select case(space%periodic_dim)
932 write(iunit,
'(a)')
'CRYSTAL'
934 write(iunit,
'(a)')
'SLAB'
936 write(iunit,
'(a)')
'POLYMER'
940 write(iunit,
'(a)')
'PRIMVEC'//trim(index_str)
942 do idir = 1, space%dim
944 latt%rlattice(idir2, idir)*supercell(idir)), idir2 = 1, space%dim)
947 write(iunit,
'(a)')
'PRIMCOORD'//trim(index_str)
948 if(.not.
present(extra_atom))
then
949 write(iunit,
'(i10, a)')
size(atoms)*nreplica,
' 1'
951 write(iunit,
'(i10, a)')
size(atoms)*nreplica+1,
' 1'
954 write(iunit,
'(a)')
'ATOMS'//trim(index_str)
958 do irep = 1, nreplica
960 do iatom = 1,
size(atoms)
961 write(iunit,
'(a10, 3f12.6)', advance=
'no') trim(atoms(iatom)%label), &
963 - offset(idir)), idir = 1, space%dim)
967 write(iunit,
'(a10, 3f12.6)', advance=
'no')
'X', &
975#if defined(HAVE_NETCDF)
977 subroutine ncdf_error(func, status, filename, namespace, ierr)
978 character(len=*),
intent(in) :: func
979 integer,
intent(in) :: status
980 character(len=*),
intent(in) :: filename
982 integer,
intent(inout) :: ierr
986 if (status == nf90_noerr)
then
991 write(
message(1),
'(3a)')
"NETCDF error in function '" , trim(func) ,
"'"
992 write(
message(2),
'(3a)')
"(reading/writing ", trim(filename) ,
")"
993 write(
message(3),
'(6x,a,a)')
'Error code = ', trim(nf90_strerror(status))
1003 real(real64),
intent(in) :: in(:, :, :)
1004 real(real64),
intent(out) :: out(:, :, :)
1005 integer :: ix, iy, iz
1009 do ix = lbound(in, 1), ubound(in, 1)
1010 do iy = lbound(in, 2), ubound(in, 2)
1011 do iz = lbound(in, 3), ubound(in, 3)
1012 out(iz, iy, ix) = in(ix, iy, iz)
1022#include "io_function_inc.F90"
1025#include "complex.F90"
1026#include "io_function_inc.F90"
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_zero
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. Includes information about ...
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, public dout_cf_netcdf(filename, ierr, cf, cube, space, spacing, transpose, unit, namespace)
Writes a cube_function in netcdf format.
subroutine dio_function_output_global_bz(how, dir, fname, namespace, kpoints, ff, unit, ierr)
subroutine ncdf_error(func, status, filename, namespace, ierr)
integer, parameter, private zoutput_kind
subroutine, public zio_function_output_global(how, dir, fname, namespace, space, mesh, ff, unit, ierr)
subroutine, public zout_cf_netcdf(filename, ierr, cf, cube, space, spacing, transpose, unit, namespace)
Writes a cube_function in netcdf format.
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)
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
Some general things and nomenclature:
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.