compute forces from electrostatic interactions
fix for multiple subsystems
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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(:) | :: | gr_sr_atom_start | index of the first short-range atom for this process group |
|
integer, | intent(in), | dimension(:) | :: | gr_sr_atom_end | index of the last short-range atom for this process group |
|
integer, | intent(in), | dimension(:) | :: | gr_lr_atom_start | index of the first long-ange atom for this process group |
|
integer, | intent(in), | dimension(:) | :: | gr_lr_atom_end | index of the last long-range atom for this process group |
|
integer, | intent(in), | dimension(:) | :: | pr_sr_atom_start | index of the first short-range atom for this process |
|
integer, | intent(in), | dimension(:) | :: | pr_sr_atom_end | index of the last short-range atom for this process |
|
integer, | intent(in), | dimension(:) | :: | pr_lr_atom_start | index of the first long-ange atom for this process |
|
integer, | intent(in), | dimension(:) | :: | pr_lr_atom_end | index of the last long-range atom for this process |
|
logical, | intent(in) | :: | root | logical for root process |
subroutine mimic_compute_forces(subsystems, &
quantum_fragment, &
x_start, &
x_end, &
gr_sr_atom_start, &
gr_sr_atom_end, &
gr_lr_atom_start, &
gr_lr_atom_end, &
pr_sr_atom_start, &
pr_sr_atom_end, &
pr_lr_atom_start, &
pr_lr_atom_end, &
root)
!> 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 group
integer, dimension(:), intent(in) :: gr_sr_atom_start
!> index of the last short-range atom for this process group
integer, dimension(:), intent(in) :: gr_sr_atom_end
!> index of the first long-ange atom for this process group
integer, dimension(:), intent(in) :: gr_lr_atom_start
!> index of the last long-range atom for this process group
integer, dimension(:), intent(in) :: gr_lr_atom_end
!> index of the first short-range atom for this process
integer, dimension(:), intent(in) :: pr_sr_atom_start
!> index of the last short-range atom for this process
integer, dimension(:), intent(in) :: pr_sr_atom_end
!> index of the first long-ange atom for this process
integer, dimension(:), intent(in) :: pr_lr_atom_start
!> index of the last long-range atom for this process
integer, dimension(:), intent(in) :: pr_lr_atom_end
!> logical for root process
logical, intent(in) :: root
integer :: i
call timer_start("mimic_compute_forces")
do i = 1, size(subsystems)
call compute_sr_forces_electrons(subsystems(i)%sr_atoms, &
quantum_fragment%density%field, &
quantum_fragment%cell, &
x_start, x_end, &
gr_sr_atom_start(i), &
gr_sr_atom_end(i))
call compute_sr_forces_ions(subsystems(i)%sr_atoms, &
quantum_fragment%nuclei, &
pr_sr_atom_start(i), &
pr_sr_atom_end(i))
if (subsystems(i)%num_lr_atoms > 0) then
call compute_lr_forces(subsystems(i)%lr_atoms, &
quantum_fragment, &
pr_lr_atom_start(i), &
pr_lr_atom_end(i))
end if
end do
call timer_stop
end subroutine mimic_compute_forces