field_M2E Subroutine

public subroutine field_M2E(eel, cpt, E)

This subroutine computes the potential generated by the static multipoles to a set of arbitrary coordinates, without applying any screening rules.

Arguments

Type IntentOptional Attributes Name
type(ommp_electrostatics_type), intent(in) :: eel

Electrostatics data structure

real(kind=rp), intent(in) :: cpt(:,:)

Coordinates at which the electric field is requested

real(kind=rp), intent(inout) :: E(:,:)

Electric field (results will be added)


Calls

proc~~field_m2e~~CallsGraph proc~field_m2e field_M2E proc~coulomb_kernel coulomb_kernel proc~field_m2e->proc~coulomb_kernel proc~q_elec_prop q_elec_prop proc~field_m2e->proc~q_elec_prop proc~mu_elec_prop mu_elec_prop proc~field_m2e->proc~mu_elec_prop proc~quad_elec_prop quad_elec_prop proc~field_m2e->proc~quad_elec_prop proc~fatal_error fatal_error proc~coulomb_kernel->proc~fatal_error proc~ommp_message ommp_message proc~fatal_error->proc~ommp_message proc~close_output close_output proc~fatal_error->proc~close_output proc~close_output->proc~ommp_message

Called by

proc~~field_m2e~~CalledByGraph proc~field_m2e field_M2E proc~ommp_field_mmpol2ext ommp_field_mmpol2ext proc~ommp_field_mmpol2ext->proc~field_m2e proc~electrostatic_for_grad electrostatic_for_grad proc~electrostatic_for_grad->proc~field_m2e proc~ommp_field_mm2ext ommp_field_mm2ext proc~ommp_field_mm2ext->proc~field_m2e proc~c_ommp_field_mm2ext C_ommp_field_mm2ext proc~c_ommp_field_mm2ext->proc~field_m2e proc~c_ommp_field_mmpol2ext C_ommp_field_mmpol2ext proc~c_ommp_field_mmpol2ext->proc~ommp_field_mmpol2ext proc~c_ommp_prepare_qm_ele_grd C_ommp_prepare_qm_ele_grd proc~c_ommp_prepare_qm_ele_grd->proc~electrostatic_for_grad

Contents

Source Code


Source Code

    subroutine field_M2E(eel, cpt, E)
        !! This subroutine computes the potential generated by the static
        !! multipoles to a set of arbitrary coordinates, without applying
        !! any screening rules.
        
        implicit none

        type(ommp_electrostatics_type), intent(in) :: eel
        !! Electrostatics data structure
        real(rp), intent(inout) :: E(:,:)
        !! Electric field (results will be added)
        real(rp), intent(in) :: cpt(:,:)
        !! Coordinates at which the electric field is requested

        integer(ip) :: i, j, n_cpt
        real(rp) :: kernel(5), dr(3), tmpV, tmpE(3), tmpEgr(6), tmpHE(10)

        n_cpt = size(cpt, 2)

        if(eel%amoeba) then
            !$omp parallel do default(shared) schedule(dynamic) &
            !$omp private(i,j,dr,kernel,tmpE,tmpV,tmpEgr,tmpHE) reduction(+:E)
            do i=1, eel%top%mm_atoms
                do j=1, n_cpt
                    dr = cpt(:,j) - eel%top%cmm(:,i)
                    call coulomb_kernel(dr, 4, kernel)
                    tmpE = 0.0_rp
                    
                    call q_elec_prop(eel%q(1,i), dr, kernel, .false., tmpV, &
                                     .true., tmpE, .false., tmpEgr, & 
                                     .false., tmpHE)
                    call mu_elec_prop(eel%q(2:4,i), dr, kernel, .false., tmpV, &
                                      .true., tmpE, .false., tmpEgr, & 
                                      .false., tmpHE)
                    call quad_elec_prop(eel%q(5:10,i), dr, kernel, .false., tmpV, &
                                        .true., tmpE, .false., tmpEgr, & 
                                        .false., tmpHE)

                    E(:,j) = E(:,j) + tmpE
                end do
            end do
        else
            !$omp parallel do default(shared) schedule(dynamic) &
            !$omp private(i,j,dr,kernel,tmpE,tmpV,tmpEgr,tmpHE) reduction(+:E)
            do i=1, eel%top%mm_atoms
                ! loop on sources
                do j=1, n_cpt
                    dr = cpt(:,j) - eel%top%cmm(:,i)
                    call coulomb_kernel(dr, 2, kernel)
                    tmpE = 0.0_rp
                    
                    call q_elec_prop(eel%q(1,i), dr, kernel, .false., tmpV, &
                                     .true., tmpE, .false., tmpEgr, & 
                                     .false., tmpHE)
                    
                    E(:,j) = E(:,j) + tmpE
                end do
            end do
        end if
    end subroutine field_M2E