38  use, 
intrinsic :: iso_fortran_env
 
   77    real(real64),         
pointer :: density(:, :)
 
   78    type(states_elec_t),  
pointer :: st
 
   79    type(grid_t),         
pointer :: gr
 
   80    type(accel_mem_t)             :: buff_density
 
   93    type(density_calc_t),           
intent(out)   :: this
 
   94    type(states_elec_t),  
target,   
intent(in)    :: st
 
   95    type(grid_t),         
target,   
intent(in)    :: gr
 
   96    real(real64),         
target,   
intent(out)   :: density(:, :)
 
   99    logical :: correct_size
 
  106    this%density => density
 
  108    do ispin = 1, st%d%nspin
 
  111        this%density(ip, ispin) = 
m_zero 
  117    this%packed = .false.
 
  119    correct_size = ubound(this%density, dim = 1) == this%gr%np .or. &
 
  120      ubound(this%density, dim = 1) == this%gr%np_part
 
  130    type(density_calc_t),           
intent(inout)   :: this
 
  131    logical, 
optional,              
intent(in)      :: async
 
  151    type(density_calc_t),         
intent(inout) :: this
 
  152    type(wfs_elec_t),             
intent(in)    :: psib
 
  153    integer,                      
intent(in)    :: istin
 
  155    integer :: ist, ip, ispin, istin_, wgsize
 
  156    complex(real64)   :: term, psi1, psi2
 
  157    real(real64), 
allocatable :: weight(:)
 
  158    type(accel_mem_t) :: buff_weight
 
  159    type(accel_kernel_t), 
pointer :: kernel
 
  164    ispin = this%st%d%get_spin_index(psib%ik)
 
  168    safe_allocate(weight(1:psib%nst))
 
  170      weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
 
  171      if (psib%ist(ist) == istin) istin_ = ist
 
  183    select case (psib%status())
 
  185      select case (this%st%d%ispin)
 
  189          do ip = 1, this%gr%np
 
  190            this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
 
  194          do ip = 1, this%gr%np
 
  195            this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
 
  196              real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
 
  201        do ip = 1, this%gr%np
 
  202          psi1 = psib%zff(ip, 1, ist)
 
  203          psi2 = psib%zff(ip, 2, ist)
 
  204          this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
 
  205          this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
 
  206          term = weight(ist)*psi1*conjg(psi2)
 
  207          this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
 
  208          this%density(ip, 4) = this%density(ip, 4) + aimag(term)
 
  214      select case (this%st%d%ispin)
 
  218          do ip = 1, this%gr%np
 
  219            this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
 
  224          do ip = 1, this%gr%np
 
  225            this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
 
  226              real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
 
  230        assert(mod(psib%nst_linear, 2) == 0)
 
  232        do ip = 1, this%gr%np
 
  233          psi1 = psib%zff_pack(2*ist - 1, ip)
 
  234          psi2 = psib%zff_pack(2*ist,     ip)
 
  235          term = weight(ist)*psi1*conjg(psi2)
 
  237          this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
 
  238          this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
 
  239          this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
 
  240          this%density(ip, 4) = this%density(ip, 4) + aimag(term)
 
  250      select case (this%st%d%ispin)
 
  286    safe_deallocate_a(weight)
 
  299    logical, 
optional,            
intent(in)    :: async
 
  301    integer :: ist, ip, ispin, wgsize
 
  302    complex(real64)   :: term, psi1, psi2
 
  303    real(real64), 
allocatable :: weight(:)
 
  310    ispin = this%st%d%get_spin_index(psib%ik)
 
  312    safe_allocate(weight(1:psib%nst))
 
  314      weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
 
  317    select case (psib%status())
 
  319      select case (this%st%d%ispin)
 
  325            do ip = 1, this%gr%np
 
  326              this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
 
  333            do ip = 1, this%gr%np
 
  334              this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
 
  335                real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
 
  343          do ip = 1, this%gr%np
 
  344            psi1 = psib%zff(ip, 1, ist)
 
  345            psi2 = psib%zff(ip, 2, ist)
 
  346            this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
 
  347            this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
 
  348            term = weight(ist)*psi1*conjg(psi2)
 
  349            this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
 
  350            this%density(ip, 4) = this%density(ip, 4) + aimag(term)
 
  357      select case (this%st%d%ispin)
 
  361          do ip = 1, this%gr%np
 
  363              this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
 
  368          do ip = 1, this%gr%np
 
  370              this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
 
  371                real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
 
  376        assert(mod(psib%nst_linear, 2) == 0)
 
  378        do ip = 1, this%gr%np
 
  380            psi1 = psib%zff_pack(2*ist - 1, ip)
 
  381            psi2 = psib%zff_pack(2*ist,     ip)
 
  382            term = weight(ist)*psi1*conjg(psi2)
 
  384            this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
 
  385            this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
 
  386            this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
 
  387            this%density(ip, 4) = this%density(ip, 4) + aimag(term)
 
  398      select case (this%st%d%ispin)
 
  436    safe_deallocate_a(weight)
 
  439    if (this%st%d%ispin /= 
spinors) 
then 
  465    logical, 
optional,    
intent(in)    :: allreduce
 
  466    logical, 
optional,    
intent(in)    :: symmetrize
 
  469    real(real64), 
allocatable :: tmpdensity(:)
 
  473    if (this%packed) 
then 
  474      safe_allocate(tmpdensity(1:this%gr%np))
 
  477      do ispin = 1, this%st%d%nspin
 
  478        call accel_read_buffer(this%buff_density, int(this%gr%np, int64), tmpdensity, offset = (ispin - 1)*this%pnp)
 
  481        do ip = 1, this%gr%np
 
  482          this%density(ip, ispin) = this%density(ip, ispin) + tmpdensity(ip)
 
  486      this%packed = .false.
 
  488      safe_deallocate_a(tmpdensity)
 
  492    if ((this%st%parallel_in_states .or. this%st%d%kpt%parallel) .and. 
optional_default(allreduce, .
true.)) 
then 
  494      call comm_allreduce(this%st%st_kpt_mpi_grp, this%density, dim = (/this%gr%np, this%st%d%nspin/))
 
  499      do ispin = 1, this%st%d%nspin
 
  516    type(
grid_t),            
intent(in)  :: gr
 
  517    real(real64),            
intent(out) :: density(:, :)
 
  518    integer, 
optional,       
intent(in)  :: istin
 
  525    assert(ubound(density, dim = 1) == gr%np .or. ubound(density, dim = 1) == gr%np_part)
 
  529    if (.not. 
present(istin)) 
then 
  531      do ik = st%d%kpt%start, st%d%kpt%end
 
  532        do ib = st%group%block_start, st%group%block_end
 
  539      do ik = st%d%kpt%start, st%d%kpt%end
 
  540        do ib = st%group%block_start, st%group%block_end
 
  541          if (st%group%psib(ib, ik)%ist(1) <= istin .and. st%group%psib(ib, ik)%ist(st%group%psib(ib, ik)%nst) >= istin) 
then 
  560    class(
space_t),      
intent(in)    :: space
 
  561    type(
grid_t),        
intent(in)    :: gr
 
  564    integer,             
intent(in)    :: n
 
  565    logical,             
intent(in)    :: family_is_mgga
 
  567    integer :: ist, ik, ib, nblock
 
  568    integer :: nodeto, nodefr
 
  570    complex(real64), 
allocatable :: psi(:, :, :), rec_buffer(:,:)
 
  576    if (n >= st%nst) 
then 
  577      write(
message(1),
'(a)') 
'Attempting to freeze a number of orbitals which is larger or equal to' 
  578      write(
message(2),
'(a)') 
'the total number. The program has to stop.' 
  584    if (.not. 
allocated(st%frozen_rho)) 
then 
  585      safe_allocate(st%frozen_rho(1:gr%np, 1:st%d%nspin))
 
  590    do ik = st%d%kpt%start, st%d%kpt%end
 
  591      if (n < st%st_start) cycle
 
  593      do ib =  st%group%block_start, st%group%block_end
 
  604          safe_allocate(psi(1:gr%np, 1:st%d%dim, 1:nblock))
 
  613          safe_deallocate_a(psi)
 
  624    if (family_is_mgga) 
then 
  625      if (.not. 
allocated(st%frozen_tau)) 
then 
  626        safe_allocate(st%frozen_tau(1:gr%np, 1:st%d%nspin))
 
  628      if (.not. 
allocated(st%frozen_gdens)) 
then 
  629        safe_allocate(st%frozen_gdens(1:gr%np, 1:space%dim, 1:st%d%nspin))
 
  631      if (.not. 
allocated(st%frozen_ldens)) 
then 
  632        safe_allocate(st%frozen_ldens(1:gr%np, 1:st%d%nspin))
 
  636        density_gradient = st%frozen_gdens, density_laplacian = st%frozen_ldens, st_end = n)
 
  644    safe_allocate(psi(1:gr%np, 1:st%d%dim, 1))
 
  645    safe_allocate(rec_buffer(1:gr%np, 1:st%d%dim))
 
  647    do ik = st%d%kpt%start, st%d%kpt%end
 
  651        nodeto = st%node(ist)
 
  653        nodefr = staux%node(ist+n)
 
  655        if (st%mpi_grp%rank == nodeto .and. st%mpi_grp%rank == nodefr) 
then 
  661          if (st%mpi_grp%rank == nodefr) 
then 
  663            call st%mpi_grp%send(psi(1, 1, 1), gr%np*st%d%dim, mpi_double_complex, nodeto, ist)
 
  665          if (st%mpi_grp%rank == nodeto) 
then 
  666            call st%mpi_grp%recv(rec_buffer(1, 1), gr%np*st%d%dim, mpi_double_complex, nodefr, ist)
 
  673    safe_deallocate_a(psi)
 
  674    safe_deallocate_a(rec_buffer)
 
  678        st%occ(ist, ik) = staux%occ(n+ist, ik)
 
  679        st%eigenval(ist, ik) = staux%eigenval(n+ist, ik)
 
  691    class(
mesh_t),       
intent(in)    :: mesh
 
  693    integer,             
intent(in)    :: nn
 
  703    safe_deallocate_a(st%eigenval)
 
  704    safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
 
  705    st%eigenval = huge(st%eigenval)
 
  707    safe_deallocate_a(st%occ)
 
  708    safe_allocate(st%occ     (1:st%nst, 1:st%nik))
 
  730    do ik = st%d%kpt%start, st%d%kpt%end
 
  731      do ist = st%st_start, st%st_end
 
  732        st%qtot = st%qtot + st%occ(ist, ik) * st%kweights(ik)
 
  736    if (st%parallel_in_states .or. st%d%kpt%parallel) 
then 
  757    class(
mesh_t),        
intent(in)  :: mesh
 
  758    real(real64),         
intent(out) :: total_rho(:,:)
 
  764    do is = 1, st%d%nspin
 
  766        total_rho(ip, is) = st%rho(ip, is)
 
  770    if (
allocated(st%rho_core)) 
then 
  772      do is = 1, st%d%spin_channels
 
  775          total_rho(ip, is) = total_rho(ip, is) + st%rho_core(ip)/st%d%spin_channels
 
  783    if (
allocated(st%frozen_rho)) 
then 
  785      do is = 1, st%d%nspin
 
  788          total_rho(ip, is) = total_rho(ip, is) + st%frozen_rho(ip, is)
 
  800#include "density_inc.F90" 
  803#include "complex.F90" 
  804#include "density_inc.F90" 
type(accel_kernel_t), target, save, public kernel_density_real
 
subroutine, public accel_finish()
 
type(accel_kernel_t), target, save, public kernel_density_spinors
 
integer, parameter, public accel_mem_read_write
 
subroutine, public accel_release_buffer(this)
 
type(accel_kernel_t), target, save, public kernel_density_complex
 
integer function, public accel_kernel_workgroup_size(kernel)
 
integer, parameter, public accel_mem_read_only
 
This module implements batches of mesh functions.
 
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
 
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
 
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
 
This module implements common operations on batches of mesh functions.
 
This module implements a calculator for the density and defines related functions.
 
subroutine, public density_calc_end(this, allreduce, symmetrize)
Finalize the density calculation.
 
subroutine density_calc_state(this, psib, istin)
Calculate contribution to the density from state istin.
 
subroutine, public density_calc_accumulate(this, psib, async)
Accumulate weighted orbital densities for a batch psib.
 
subroutine, public states_elec_freeze_adjust_qtot(st)
 
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
 
subroutine, public density_calc_init(this, st, gr, density)
initialize the density calculator
 
subroutine density_calc_pack(this, async)
prepare the density calculator for GPU use
 
subroutine, public ddensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate  within a batch.
 
subroutine, public zdensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate  within a batch.
 
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
 
subroutine, public states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
Calculate partial density for frozen orbitals.
 
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.
 
integer, parameter, public unpolarized
Parameters...
 
integer, parameter, public spinors
 
integer, parameter, public spin_polarized
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public m_epsilon
 
This module implements the underlying real-space grid.
 
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
 
This module is intended to contain "only mathematical" functions and procedures.
 
This module defines various routines, operating on mesh functions.
 
This module defines the meshes, which are used in Octopus.
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
 
This module handles the communicators for the various parallelization strategies.
 
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.
 
integer, parameter, public smear_fixed_occ
 
pure logical function, public states_are_complex(st)
 
pure logical function, public states_are_real(st)
 
This module handles spin dimensions of the states and the k-point distribution.
 
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
 
subroutine, public states_elec_distribute_nodes(st, namespace, mc)
@Brief. Distribute states over the processes for states parallelization
 
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
 
subroutine, public states_elec_end(st)
finalize the states_elec_t object
 
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
 
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
 
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
 
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
 
type(type_t), public type_float
 
type(type_t), public type_cmplx
 
Description of the grid, containing information on derivatives, stencil, and symmetries.
 
Describes mesh distribution to nodes.
 
Stores all communicators and groups.
 
The states_elec_t class contains all electronic wave functions.
 
batches of electronic states