51    real(real64), 
allocatable, 
public :: v_old(:, :, :)
 
   52    real(real64), 
allocatable, 
public :: vtau_old(:, :, :)
 
   53    logical                    :: mgga_with_exc
 
   56  integer :: interpolation_steps
 
   62    type(potential_interpolation_t), 
intent(inout) :: vkso
 
   63    type(potential_interpolation_t), 
intent(in)    :: vksi
 
   67    safe_allocate_source_a(vkso%v_old, vksi%v_old)
 
   68    safe_allocate_source_a(vkso%vtau_old, vksi%vtau_old)
 
   69    vkso%mgga_with_exc = vksi%mgga_with_exc
 
   77    type(potential_interpolation_t), 
intent(inout) :: potential_interpolation
 
   78    integer, 
intent(in) :: np, nspin
 
   79    logical, 
intent(in) :: mgga_with_exc
 
   80    integer, 
optional, 
intent(in) :: order
 
   86    safe_allocate(potential_interpolation%v_old(1:np, 1:nspin, 0:interpolation_steps))
 
   87    potential_interpolation%v_old(:, :, :) = 
m_zero 
   89    potential_interpolation%mgga_with_exc = mgga_with_exc
 
   91    if (potential_interpolation%mgga_with_exc) 
then 
   93      safe_allocate(potential_interpolation%vtau_old(1:np, 1:nspin, 0:interpolation_steps))
 
   94      potential_interpolation%vtau_old(:, :, :) = 
m_zero 
  104    type(potential_interpolation_t), 
intent(inout) :: potential_interpolation
 
  108    assert(
allocated(potential_interpolation%v_old))
 
  109    safe_deallocate_a(potential_interpolation%v_old)
 
  110    safe_deallocate_a(potential_interpolation%vtau_old)
 
  118    type(potential_interpolation_t), 
intent(inout) :: potential_interpolation
 
  119    integer,           
intent(in)    :: np, nspin
 
  120    real(real64),      
intent(in)    :: vhxc(:, :)
 
  121    real(real64), 
optional,   
intent(in)    :: vtau(:, :)
 
  123    integer :: i, ispin, ip
 
  127    do i = 1, interpolation_steps
 
  130          potential_interpolation%v_old(ip, ispin, i) = vhxc(ip, ispin)
 
  135    if (potential_interpolation%mgga_with_exc) 
then 
  136      assert(
present(vtau))
 
  137      do i = 1, interpolation_steps
 
  140            potential_interpolation%vtau_old(ip, ispin, i) = vtau(ip, ispin)
 
  153    integer,           
intent(in)    :: np, nspin
 
  154    real(real64),      
intent(in)    :: time, dt
 
  155    real(real64), 
contiguous, 
intent(in)    :: vhxc(:, :)
 
  156    real(real64), 
optional,   
intent(in)    :: vtau(:, :)
 
  158    real(real64), 
allocatable :: times(:)
 
  163    safe_allocate(times(1:interpolation_steps))
 
  164    do j = 1, interpolation_steps
 
  165      times(j) = time - j*dt
 
  168    do j = interpolation_steps, 2, -1
 
  169      call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, j-1), potential_interpolation%v_old(:, :, j))
 
  171    call lalg_copy(np, nspin, vhxc(:, :),     potential_interpolation%v_old(:, :, 1))
 
  172    call interpolate( times, potential_interpolation%v_old(:, :, 1:interpolation_steps), &
 
  173      time, potential_interpolation%v_old(:, :, 0))
 
  175    if (potential_interpolation%mgga_with_exc) 
then 
  176      assert(
present(vtau))
 
  177      do j = interpolation_steps, 2, -1
 
  178        call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, j-1), potential_interpolation%vtau_old(:, :, j))
 
  180      call lalg_copy(np, nspin, vtau(:, :),     potential_interpolation%vtau_old(:, :, 1))
 
  181      call interpolate( times, potential_interpolation%vtau_old(:, :, 1:interpolation_steps), &
 
  182        time, potential_interpolation%vtau_old(:, :, 0))
 
  185    safe_deallocate_a(times)
 
  193    integer,           
intent(in)    :: np, nspin
 
  194    integer,           
intent(in)    :: i
 
  195    real(real64),      
intent(inout) :: vhxc(:, :)
 
  196    real(real64), 
optional,   
intent(in)    :: vtau(:, :)
 
  200    call lalg_copy(np, nspin, vhxc, potential_interpolation%v_old(:, :, i))
 
  202    if (potential_interpolation%mgga_with_exc) 
then 
  203      assert(
present(vtau))
 
  204      call lalg_copy(np, nspin, vtau, potential_interpolation%vtau_old(:, :, i))
 
  214    integer,                              
intent(in)    :: np, nspin
 
  215    integer,                              
intent(in)    :: i
 
  216    real(real64), 
contiguous,             
intent(inout) :: vhxc(:, :)
 
  217    real(real64), 
optional, 
contiguous,   
intent(inout) :: vtau(:, :)
 
  221    call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, i), vhxc)
 
  223    if (potential_interpolation%mgga_with_exc) 
then 
  224      assert(
present(vtau))
 
  225      call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, i), vtau)
 
  236    integer,           
intent(in)    :: order
 
  237    real(real64),      
intent(in)    :: time, dt, t
 
  238    real(real64),      
intent(inout) :: vhxc(:, :)
 
  239    real(real64), 
optional,   
intent(inout) :: vtau(:, :)
 
  242    real(real64), 
allocatable :: times(:)
 
  246    safe_allocate(times(1:interpolation_steps))
 
  247    do j = 1, interpolation_steps
 
  248      times(j) = time - (j-1)*dt
 
  251    call interpolate( times(1:order), potential_interpolation%v_old(:, :, 0:order-1), t, vhxc(:, :))
 
  253    if (potential_interpolation%mgga_with_exc) 
then 
  254      assert(
present(vtau))
 
  255      call interpolate( times(1:order), potential_interpolation%vtau_old(:, :, 0:order-1), t, vtau(:, :))
 
  258    safe_deallocate_a(times)
 
  267    class(
space_t),    
intent(in)    :: space
 
  268    class(
mesh_t),     
intent(in)    :: mesh
 
  269    integer,           
intent(in)    :: nspin
 
  270    integer,           
intent(out)   :: err2
 
  272    integer :: ii, is, err
 
  273    character(len=256) :: filename
 
  279    do ii = 1, interpolation_steps - 1
 
  281        write(filename,
'(a6,i2.2,i3.3)') 
'vprev_', ii, is
 
  283          potential_interpolation%v_old(1:mesh%np, is, ii), err)
 
  286        if (err /= 0) err2 = err2 + 1
 
  291    if (potential_interpolation%mgga_with_exc) 
then 
  293      do ii = 1, interpolation_steps - 1
 
  295          write(filename,
'(a9,i2.2,i3.3)') 
'vtauprev_', ii, is
 
  297            potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
 
  300          if (err /= 0) err2 = err2 + 1
 
  313    class(
space_t),    
intent(in)    :: space
 
  315    class(
mesh_t),     
intent(in)    :: mesh
 
  316    integer,           
intent(in)    :: nspin
 
  317    integer,           
intent(out)   :: err2
 
  319    integer :: ii, is, err
 
  320    character(len=256) :: filename
 
  325    do ii = 1, interpolation_steps - 1
 
  327        write(filename,
'(a,i2.2,i3.3)') 
'vprev_', ii, is
 
  329          potential_interpolation%v_old(1:mesh%np, is, ii), err)
 
  332          message(1) = 
"Unable to read VKS restart file '" // trim(filename) // 
"'" 
  338    if (potential_interpolation%mgga_with_exc) 
then 
  340      do ii = 1, interpolation_steps - 1
 
  342          write(filename,
'(a,i2.2,i3.3)') 
'vtauprev_', ii, is
 
  344            potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
 
  347            message(1) = 
"Unable to read VKS restart file '" // trim(filename) // 
"'" 
Copies a vector x, to a vector y.
 
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
 
real(real64), parameter, public m_zero
 
This module is intended to contain "only mathematical" functions and procedures.
 
This module defines the meshes, which are used in Octopus.
 
subroutine, public messages_warning(no_lines, all_nodes, namespace)
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public potential_interpolation_copy(vkso, vksi)
 
subroutine, public potential_interpolation_end(potential_interpolation)
 
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
 
subroutine, public potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
 
subroutine, public potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
 
subroutine, public potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
 
subroutine, public potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
 
subroutine, public potential_interpolation_load(potential_interpolation, namespace, space, restart, mesh, nspin, err2)
 
subroutine, public potential_interpolation_dump(potential_interpolation, space, restart, mesh, nspin, err2)
 
subroutine, public potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
 
Describes mesh distribution to nodes.