36#if defined(HAVE_ETSF_IO) 
   52  use, 
intrinsic :: iso_fortran_env
 
  112  subroutine output_init(outp, namespace, space, st, gr, nst, ks)
 
  113    type(output_t),            
intent(out)   :: outp
 
  114    type(namespace_t),         
intent(in)    :: namespace
 
  115    class(space_t),            
intent(in)    :: space
 
  116    type(states_elec_t),       
intent(in)    :: st
 
  117    type(grid_t),              
intent(in)    :: gr
 
  118    integer,                   
intent(in)    :: nst
 
  119    type(v_ks_t),              
intent(inout) :: ks
 
  123    character(len=80) :: nst_string, default
 
  130    if (outp%what(option__output__wfs_fourier)) 
then 
  132        message(1) = 
"Wave functions in Fourier space not supported on GPUs." 
  139    if (outp%what(option__output__elf) .or. outp%what(option__output__elf_basins)) 
then 
  140      if (space%dim /= 2 .and. space%dim /= 3) 
then 
  141        outp%what(option__output__elf) = .false.
 
  142        outp%what(option__output__elf_basins) = .false.
 
  143        write(
message(1), 
'(a)') 
'Cannot calculate ELF except in 2D and 3D.' 
  149    if (outp%what(option__output__mmb_wfs)) 
then 
  153    if (outp%what(option__output__xc_torque)) 
then 
  154      if (st%d%ispin /= 
spinors) 
then 
  155        write(
message(1), 
'(a)') 
'The output xc_torque can only be computed for spinors.' 
  158      if (space%dim /= 3) 
then 
  159        write(
message(1), 
'(a)') 
'The output xc_torque can only be computed in the 3D case.' 
  163    if (outp%what(option__output__mmb_den)) 
then 
  172    if (outp%what(option__output__energy_density)) 
call messages_experimental(
"'Output = energy_density'", namespace=namespace)
 
  173    if (outp%what(option__output__heat_current)) 
call messages_experimental(
"'Output = heat_current'", namespace=namespace)
 
  175    if (outp%what(option__output__wfs) .or. outp%what(option__output__wfs_sqmod)) 
then 
  187      write(nst_string,
'(i6)') nst
 
  188      write(default,
'(a,a)') 
"1-", trim(adjustl(nst_string))
 
  189      call parse_variable(namespace, 
'OutputWfsNumber', default, outp%wfs_list)
 
  192    if (
parse_block(namespace, 
'CurrentThroughPlane', blk) == 0) 
then 
  193      if (.not. outp%what(option__output__j_flow)) 
then 
  194        outp%what(option__output__j_flow) = .
true.
 
  195        call parse_variable(namespace, 
'OutputInterval', 50, outp%output_interval(option__output__j_flow))
 
  246      select case (space%dim)
 
  264        norm = norm2(outp%plane%u(1:3))
 
  266          write(
message(1), 
'(a)') 
'u-vector for CurrentThroughPlane cannot have norm zero.' 
  269        outp%plane%u(1:3) = outp%plane%u(1:3) / norm
 
  271        norm = norm2(outp%plane%v(1:3))
 
  273          write(
message(1), 
'(a)') 
'v-vector for CurrentThroughPlane cannot have norm zero.' 
  276        outp%plane%v(1:3) = outp%plane%v(1:3) / norm
 
  278        outp%plane%n(1) = outp%plane%u(2)*outp%plane%v(3) - outp%plane%u(3)*outp%plane%v(2)
 
  279        outp%plane%n(2) = outp%plane%u(3)*outp%plane%v(1) - outp%plane%u(1)*outp%plane%v(3)
 
  280        outp%plane%n(3) = outp%plane%u(1)*outp%plane%v(2) - outp%plane%u(2)*outp%plane%v(1)
 
  292        norm = norm2(outp%line%u(1:2))
 
  294          write(
message(1), 
'(a)') 
'u-vector for CurrentThroughPlane cannot have norm zero.' 
  297        outp%line%u(1:2) = outp%line%u(1:2) / norm
 
  299        outp%line%n(1) = -outp%line%u(2)
 
  300        outp%line%n(2) =  outp%line%u(1)
 
  314    if (outp%what(option__output__matrix_elements)) 
then 
  317      outp%me%what = .false.
 
  320    if (outp%what(option__output__berkeleygw)) 
then 
  322        message(1) = 
"BerkeleyGW is not compatible with GPUs." 
  329    if (outp%what(option__output__potential_gradient) .and. .not. outp%what(option__output__potential)) 
then 
  330      outp%what(option__output__potential) = .
true.
 
  331      outp%output_interval(option__output__potential) = outp%output_interval(option__output__potential_gradient)
 
  343    call parse_variable(namespace, 
'OutputDuringSCF', .false., outp%duringscf)
 
  356    call parse_variable(namespace, 
'RestartWriteInterval', 50, outp%restart_write_interval)
 
  357    if (outp%restart_write_interval <= 0) 
then 
  358      message(1) = 
"RestartWriteInterval must be > 0." 
  373    call parse_variable(namespace, 
'OutputIterDir', 
"output_iter", outp%iter_dir)
 
  374    if (any(outp%what) .and. any(outp%output_interval > 0)) 
then 
  375      call io_mkdir(outp%iter_dir, namespace)
 
  388    if (outp%what(option__output__current_dia)) 
then 
  389      message(1) = 
"The diamagnetic current will be calculated only if CalculateDiamagneticCurrent = yes." 
  397  subroutine output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
 
  400    class(
space_t),           
intent(in)    :: space
 
  401    character(len=*),         
intent(in)    :: dir
 
  402    type(
grid_t),             
intent(in)    :: gr
 
  403    type(
ions_t),             
intent(in)    :: ions
 
  404    integer,                  
intent(in)    :: iter
 
  407    type(
v_ks_t),             
intent(inout) :: ks
 
  409    integer :: idir, ierr, iout, iunit
 
  410    character(len=80) :: fname
 
  415    if (any(outp%what)) 
then 
  416      message(1) = 
"Info: Writing output to " // trim(dir)
 
  421    if (outp%what_now(option__output__mesh_r, iter)) 
then 
  422      do idir = 1, space%dim
 
  423        write(fname, 
'(a,a)') 
'mesh_r-', 
index2axis(idir)
 
  425          gr, gr%x(:,idir), 
units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
 
  429    call output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
 
  430    call output_hamiltonian(outp, namespace, space, dir, hm, st, gr%der, ions, gr, iter, st%st_kpt_mpi_grp)
 
  433    if (outp%what_now(option__output__el_pressure, iter)) 
then 
  437      if (st%d%spin_channels > 1) 
then 
  444    if (outp%what_now(option__output__j_flow, iter)) 
then 
  448    if (outp%what_now(option__output__geometry, iter)) 
then 
  449      if (
bitand(outp%how(option__output__geometry), option__outputformat__xcrysden) /= 0) 
then 
  452      if (
bitand(outp%how(option__output__geometry), option__outputformat__xyz) /= 0) 
then 
  453        call ions%write_xyz(trim(dir)//
'/geometry')
 
  454        if (ions%space%is_periodic()) 
then 
  455          call ions%write_crystal(dir)
 
  458      if (
bitand(outp%how(option__output__geometry), option__outputformat__vtk) /= 0) 
then 
  459        call ions%write_vtk_geometry(trim(dir)//
'/geometry')
 
  463    if (outp%what_now(option__output__forces, iter)) 
then 
  464      if (
bitand(outp%how(option__output__forces), option__outputformat__bild) /= 0) 
then 
  465        call ions%write_bild_forces_file(dir, 
"forces")
 
  468          gr, namespace, total_forces = ions%tot_force)
 
  472    if (outp%what_now(option__output__matrix_elements, iter)) 
then 
  473      call output_me(outp%me, namespace, space, dir, st, gr, ions, hm)
 
  476    do iout = lbound(outp%how, 1), ubound(outp%how, 1)
 
  477      if (
bitand(outp%how(iout), option__outputformat__etsf) /= 0) 
then 
  478        call output_etsf(outp, namespace, space, dir, st, gr, hm%kpoints, ions, iter)
 
  483    if (outp%what_now(option__output__berkeleygw, iter)) 
then 
  487    if (outp%what_now(option__output__energy_density, iter)) 
then 
  491    if (outp%what_now(option__output__stress, iter)) 
then 
  493      iunit = 
io_open(trim(dir)//
'/stress', namespace, action=
'write')
 
  494      call output_stress(iunit, space%periodic_dim, st%stress_tensors)
 
  499      if (outp%what_now(option__output__occ_matrices, iter))&
 
  502      if (outp%what_now(option__output__effectiveu, iter))&
 
  505      if (outp%what_now(option__output__magnetization, iter))&
 
  508      if (outp%what_now(option__output__local_orbitals, iter))&
 
  509        call output_dftu_orbitals(outp, dir, namespace, space, hm%lda_u, st, gr, ions, hm%phase%is_allocated())
 
  511      if (outp%what_now(option__output__kanamoriu, iter))&
 
  515        if (outp%what_now(option__output__photon_correlator, iter)) 
then 
  516          write(fname, 
'(a)') 
'photon_correlator' 
  517          call dio_function_output(outp%how(option__output__photon_correlator), dir, trim(fname), namespace, space, &
 
  518            gr, ks%oep_photon%pt%correlator(:,1), 
units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
 
  524    if (outp%what_now(option__output__xc_torque, iter)) 
then 
  526        write(
message(1), 
'(a)') 
'The output xc_torque can only be computed when there is a xc potential.' 
  542    class(
space_t),           
intent(in)    :: space
 
  543    character(len=*),         
intent(in)    :: dir
 
  546    type(
grid_t),             
intent(in)    :: gr
 
  547    type(
ions_t),             
intent(in)    :: ions
 
  548    integer,                  
intent(in)    :: iter
 
  550    real(real64), 
allocatable :: f_loc(:,:)
 
  551    character(len=256) :: fname
 
  552    integer :: is, ierr, imax
 
  557    mpi_grp = st%dom_st_kpt_mpi_grp
 
  563    safe_allocate(f_loc(1:gr%np, 1:imax))
 
  566    if (outp%what_now(option__output__elf, iter) .or. outp%what_now(option__output__elf_basins, iter)) 
then 
  567      assert(space%dim /= 1)
 
  569      call elf_calc(space, st, gr, hm%kpoints, f_loc)
 
  572      if (outp%what_now(option__output__elf, iter)) 
then 
  573        write(fname, 
'(a)') 
'elf_rs' 
  574        call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
 
  575          f_loc(:,imax), 
unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
 
  580            write(fname, 
'(a,i1)') 
'elf_rs-sp', is
 
  581            call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
 
  582              f_loc(:, is), 
unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
 
  588      if (outp%what_now(option__output__elf_basins, iter)) 
then 
  589        call out_basins(f_loc(:,1), 
"elf_rs_basins", outp%how(option__output__elf_basins))
 
  594    if (outp%what_now(option__output__bader, iter)) 
then 
  595      do is = 1, st%d%nspin
 
  597        if (st%d%nspin == 1) 
then 
  598          write(fname, 
'(a)') 
'bader' 
  600          write(fname, 
'(a,i1)') 
'bader-sp', is
 
  602        call dio_function_output(outp%how(option__output__bader), dir, trim(fname), namespace, space, gr, &
 
  603          f_loc(:,is), 
units_out%length**(-2 - space%dim), ierr, &
 
  604          pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
 
  606        if (st%d%nspin == 1) 
then 
  607          write(fname, 
'(a)') 
'bader_basins' 
  609          write(fname, 
'(a,i1)') 
'bader_basins-sp', is
 
  611        call out_basins(f_loc(:,1), fname, outp%how(option__output__bader))
 
  616    if (outp%what_now(option__output__el_pressure, iter)) 
then 
  618      call dio_function_output(outp%how(option__output__el_pressure), dir, 
"el_pressure", namespace, space, gr, &
 
  619        f_loc(:,1), 
unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
 
  623    safe_deallocate_a(f_loc)
 
  629    subroutine out_basins(ff, filename, output_how)
 
  630      real(real64),     
intent(in)    :: ff(:)
 
  631      character(len=*), 
intent(in)    :: filename
 
  632      integer(int64),      
intent(in)    :: output_how
 
  634      character(len=256) :: fname
 
  641      call basins_analyze(basins, namespace, gr, ff(:), st%rho, 0.01_real64)
 
  644        real(basins%map, real64) , unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
 
  647      write(fname,
'(4a)') trim(dir), 
'/', trim(filename), 
'.info' 
  648      iunit = 
io_open(trim(fname), namespace, action = 
'write')
 
  664    type(
grid_t),             
intent(in)    :: gr
 
  665    real(real64),             
intent(out)   :: pressure(:)
 
  667    real(real64), 
allocatable :: rho(:,:), lrho(:), tau(:,:)
 
  668    real(real64)   :: p_tf, dens
 
  673    safe_allocate( rho(1:gr%np_part, 1:st%d%nspin))
 
  674    safe_allocate(lrho(1:gr%np))
 
  675    safe_allocate( tau(1:gr%np, 1:st%d%nspin))
 
  682    do is = 1, st%d%spin_channels
 
  686      pressure(:) = pressure(:) + &
 
  691      dens = sum(rho(ii,1:st%d%spin_channels))
 
  698      pressure(ii) = pressure(ii) + (dens*hm%vxc(ii,1) - hm%energy%exchange - hm%energy%correlation)
 
  700      pressure(ii) = pressure(ii)/p_tf
 
  712    class(
space_t),            
intent(in) :: space
 
  713    character(len=*),          
intent(in) :: dir
 
  715    type(
v_ks_t),              
intent(inout) :: ks
 
  717    type(
ions_t),              
intent(in) :: ions
 
  718    type(
grid_t),              
intent(in) :: gr
 
  720    integer :: is, ierr, ip
 
  721    character(len=MAX_PATH_LEN) :: fname
 
  723    real(real64), 
allocatable :: energy_density(:, :)
 
  724    real(real64), 
allocatable :: ex_density(:)
 
  725    real(real64), 
allocatable :: ec_density(:)
 
  730    safe_allocate(energy_density(1:gr%np, 1:st%d%nspin))
 
  736    do is = 1, st%d%nspin
 
  738        energy_density(ip, is) = energy_density(ip, is) + st%rho(ip, is)*hm%ep%vpsl(ip)
 
  743    do is = 1, st%d%nspin
 
  745        energy_density(ip, is) = energy_density(ip, is) + 
m_half*st%rho(ip, is)*hm%vhartree(ip)
 
  750    safe_allocate(ex_density(1:gr%np))
 
  751    safe_allocate(ec_density(1:gr%np))
 
  753    call xc_get_vxc(gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, st%rho, st%d%ispin, &
 
  754      hm%ions%latt%rcell_volume, ex_density = ex_density, ec_density = ec_density)
 
  755    do is = 1, st%d%nspin
 
  757        energy_density(ip, is) = energy_density(ip, is) + ex_density(ip) + ec_density(ip)
 
  761    safe_deallocate_a(ex_density)
 
  762    safe_deallocate_a(ec_density)
 
  764    select case (st%d%ispin)
 
  766      write(fname, 
'(a)') 
'energy_density' 
  767      call dio_function_output(outp%how(option__output__energy_density), dir, trim(fname), namespace, space, gr, &
 
  768        energy_density(:,1), 
unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
 
  771        write(fname, 
'(a,i1)') 
'energy_density-sp', is
 
  772        call dio_function_output(outp%how(option__output__energy_density), dir, trim(fname), namespace, space, gr, &
 
  773          energy_density(:, is), 
unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
 
  776    safe_deallocate_a(energy_density)
 
  786    need_exx =(outp%what(option__output__berkeleygw) &
 
  787      .or. outp%me%what(option__outputmatrixelements__two_body) &
 
  788      .or. outp%me%what(option__outputmatrixelements__two_body_exc_k))
 
  795    character(len=*),    
intent(in) :: dir
 
  797    class(
space_t),      
intent(in) :: space
 
  798    type(
lda_u_t),       
intent(in) :: this
 
  800    class(
mesh_t),       
intent(in) :: mesh
 
  801    type(
ions_t),        
intent(in) :: ions
 
  802    logical,             
intent(in) :: has_phase
 
  804    integer :: ios, im, ik, idim, ierr
 
  805    complex(real64), 
allocatable :: tmp(:)
 
  806    real(real64), 
allocatable :: dtmp(:)
 
  809    character(len=32) :: fname
 
  815    if (this%basis%use_submesh) 
then 
  817        safe_allocate(dtmp(1:mesh%np))
 
  819        safe_allocate(tmp(1:mesh%np))
 
  823    do ios = 1, this%norbsets
 
  824      os => this%orbsets(ios)
 
  825      do ik = st%d%kpt%start, st%d%kpt%end
 
  826        do im = 1, this%orbsets(ios)%norbs
 
  827          do idim = 1, min(os%ndim, st%d%dim)
 
  829              if (min(os%ndim, st%d%dim) > 1) 
then 
  830                write(fname, 
'(a,i1,a,i3.3,a,i8.8,a,i1)') 
'orb', im, 
'-os', ios, 
'-k', ik, 
'-sp', idim
 
  832                write(fname, 
'(a,i1,a,i3.3,a,i8.8)') 
'orb', im, 
'-os', ios, 
'-k', ik
 
  835              if (min(os%ndim, st%d%dim) > 1) 
then 
  836                write(fname, 
'(a,i1,a,i3.3,a,i1)') 
'orb', im, 
'-os', ios, 
'-sp', idim
 
  838                write(fname, 
'(a,i1,a,i3.3)') 
'orb', im, 
'-os', ios
 
  842              if (.not. this%basis%use_submesh) 
then 
  843                call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
 
  844                  mesh, os%eorb_mesh(1:mesh%np,im,idim,ik), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  848                call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
 
  849                  mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  852              if (.not. this%basis%use_submesh) 
then 
  854                  call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
 
  855                    os%dorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  857                  call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
 
  858                    os%zorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  864                  call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
 
  865                    mesh, dtmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  869                  call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
 
  870                    mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
 
  879    safe_deallocate_a(tmp)
 
  880    safe_deallocate_a(dtmp)
 
  888    logical,             
intent(in) :: states_are_real
 
  892    if (outp%what(option__output__current) &
 
  893      .or. outp%what(option__output__current_dia) &
 
  894      .or. outp%what(option__output__heat_current) &
 
  895      .or. outp%what(option__output__current_kpt)) 
then 
  896      if (.not. states_are_real) 
then 
  899        message(1) = 
'No current density output for real states since it is identically zero.' 
  907#include "output_etsf_inc.F90" 
  909#include "output_states_inc.F90" 
  911#include "output_h_inc.F90" 
  914#include "complex.F90" 
  915#include "output_linear_response_inc.F90" 
  919#include "output_linear_response_inc.F90" 
pure logical function, public accel_is_enabled()
 
subroutine, public basins_write(this, mesh, iunit)
 
subroutine, public basins_init(this, namespace, mesh)
 
subroutine, public basins_analyze(this, namespace, mesh, f, rho, threshold)
 
subroutine, public basins_end(this)
 
This module implements a calculator for the density and defines related functions.
 
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
 
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
 
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
 
integer, parameter, public unpolarized
Parameters...
 
integer, parameter, public spinors
 
integer, parameter, public spin_polarized
 
subroutine, public elf_calc(space, st, gr, kpoints, elf, de)
(time-dependent) electron localization function, (TD)ELF.
 
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_four
 
real(real64), parameter, public m_pi
some mathematical constants
 
complex(real64), parameter, public m_z0
 
real(real64), parameter, public m_epsilon
 
real(real64), parameter, public m_half
 
real(real64), parameter, public m_one
 
real(real64), parameter, public m_three
 
real(real64), parameter, public m_five
 
This module implements the underlying real-space grid.
 
integer, parameter, public generalized_kohn_sham_dft
 
integer, parameter, public kohn_sham_dft
 
This module defines classes and functions for interaction partners.
 
subroutine, public zio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
 
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 dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
 
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
 
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)
 
subroutine, public lda_u_write_occupation_matrices(dir, this, st, namespace)
Prints the occupation matrices at the end of the scf calculation.
 
subroutine, public lda_u_write_kanamoriu(dir, st, this, namespace)
 
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
 
subroutine, public lda_u_write_effectiveu(dir, this, namespace)
 
integer, parameter, public dft_u_none
 
This module defines various routines, operating on mesh functions.
 
This module defines the meshes, which are used in Octopus.
 
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)
 
subroutine, public messages_experimental(name, namespace)
 
subroutine, public output_berkeleygw(bgw, namespace, space, dir, st, gr, ks, hm, ions)
 
subroutine, public output_berkeleygw_init(nst, namespace, bgw, periodic_dim)
 
this module contains the output system
 
subroutine, public output_me_init(this, namespace, space, st, gr, nst)
 
subroutine, public output_me(this, namespace, space, dir, st, gr, ions, hm)
 
this module contains the output system
 
subroutine calc_electronic_pressure(st, hm, gr, pressure)
 
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
 
logical function, public output_needs_current(outp, states_are_real)
 
subroutine, public output_hamiltonian(outp, namespace, space, dir, hm, st, der, ions, gr, iter, grp)
 
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
 
logical function, public output_need_exchange(outp)
 
subroutine, public output_init(outp, namespace, space, st, gr, nst, ks)
 
subroutine output_xc_torque(outp, namespace, dir, mesh, hm, st, ions, space)
 
subroutine, public output_current_flow(outp, namespace, space, dir, gr, st, kpoints)
 
subroutine, public zoutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
 
subroutine, public doutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
 
subroutine, public output_scalar_pot(outp, namespace, space, dir, mesh, ions, ext_partners, time)
 
subroutine output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
 
subroutine output_dftu_orbitals(outp, dir, namespace, space, this, st, mesh, ions, has_phase)
 
subroutine output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
 
subroutine output_etsf(outp, namespace, space, dir, st, gr, kpoints, ions, iter)
To create an etsf file one has to do the following:
 
integer function, public parse_block(namespace, name, blk, check_varinfo_)
 
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
 
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
 
pure logical function, public states_are_real(st)
 
This module defines routines to write information about states.
 
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
 
This module implements the calculation of the stress tensor.
 
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
 
subroutine, public add_last_slash(str)
Adds a '/' in the end of the string, only if it missing. Useful for directories.
 
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_out
 
type(unit_system_t), public units_inp
the units systems for reading and writing
 
type(unit_t), public unit_one
some special units required for particular quantities
 
This module is intended to contain simple general-purpose utility functions and procedures.
 
character pure function, public index2axis(idir)
 
subroutine, public v_ks_calculate_current(this, calc_cur)
 
subroutine, public xc_get_vxc(gr, xcs, st, kpoints, psolver, namespace, space, rho, ispin, rcell_volume, vxc, ex, ec, deltaxc, vtau, ex_density, ec_density, stress_xc, force_orbitalfree)
 
integer, parameter, public oep_level_full
 
subroutine out_basins(ff, filename, output_how)
 
Description of the grid, containing information on derivatives, stencil, and symmetries.
 
Class to describe DFT+U parameters.
 
Describes mesh distribution to nodes.
 
This is defined even when running serial.
 
The states_elec_t class contains all electronic wave functions.