mimic_compute_potential Subroutine

public subroutine mimic_compute_potential(subsystems, quantum_fragment, x_start, x_end, sr_atom_start, sr_atom_end, tensor_sums)

compute electrostatic potential

Arguments

TypeIntentOptionalAttributesName
type(subsystem_type), intent(inout), dimension(:):: subsystems

array of subsystems associated with client codes

type(quantum_fragment_type), intent(inout) :: quantum_fragment

quantum subsystem

integer, intent(in) :: x_start

index of the first yz-plane for this process

integer, intent(in) :: x_end

index of the last yz-plane for this process

integer, intent(in), dimension(:):: sr_atom_start

index of the first short-range atom for this process

integer, intent(in), dimension(:):: sr_atom_end

index of the last short-range atom for this process

real(kind=dp), intent(inout), dimension(:,:), allocatable:: tensor_sums

Array stroing tensor sums


Contents


Source Code

subroutine mimic_compute_potential(subsystems, &
                                 quantum_fragment, &
                                 x_start, x_end, &
                                 sr_atom_start, sr_atom_end, &
                                 tensor_sums)

    !> array of subsystems associated with client codes
    type(subsystem_type), dimension(:), intent(inout) :: subsystems
    !> quantum subsystem
    type(quantum_fragment_type), intent(inout) :: quantum_fragment
    !> index of the first yz-plane for this process
    integer, intent(in) :: x_start
    !> index of the last yz-plane for this process
    integer, intent(in) :: x_end
    !> index of the first short-range atom for this process
    integer, dimension(:), intent(in) :: sr_atom_start
    !> index of the last short-range atom for this process
    integer, dimension(:), intent(in) :: sr_atom_end
    !> Array stroing tensor sums
    real(dp), dimension(:,:), intent(inout), allocatable :: tensor_sums

    integer :: i

    call timer_start("mimic_compute_potential")

    quantum_fragment%potential%field = 0.0_dp

    do i = 1, size(subsystems)
        call compute_sr_potential_electrons(subsystems(i)%sr_atoms, &
                                    quantum_fragment%potential%field, &
                                    quantum_fragment%cell, &
                                    x_start, x_end, &
                                    sr_atom_start(i), &
                                    sr_atom_end(i))
        if (subsystems(i)%num_lr_atoms > 0) then
            call compute_lr_potential(quantum_fragment, &
                                      x_start, x_end, &
                                      tensor_sums(:,i))
        end if
    end do

    call timer_stop

end subroutine mimic_compute_potential