elec_prop_D2D Subroutine

private subroutine elec_prop_D2D(eel, in_kind, do_V, do_E, do_Egrd, do_EHes)

Computes the electric field of a trial set of induced point dipoles at polarizable sites. This is intended to be used as matrix-vector routine in the solution of the linear system.

Arguments

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

Data structure for electrostatic part of the system

character, intent(in) :: in_kind
logical, intent(in) :: do_V

Flag to control which properties have to be computed.

logical, intent(in) :: do_E

Flag to control which properties have to be computed.

logical, intent(in) :: do_Egrd

Flag to control which properties have to be computed.

logical, intent(in) :: do_EHes

Flag to control which properties have to be computed.


Calls

proc~~elec_prop_d2d~~CallsGraph proc~elec_prop_d2d elec_prop_D2D proc~fatal_error fatal_error proc~elec_prop_d2d->proc~fatal_error proc~prepare_fmm_ipd prepare_fmm_ipd proc~elec_prop_d2d->proc~prepare_fmm_ipd proc~cart_propfar_at_ipart cart_propfar_at_ipart proc~elec_prop_d2d->proc~cart_propfar_at_ipart proc~mu_elec_prop mu_elec_prop proc~elec_prop_d2d->proc~mu_elec_prop proc~damped_coulomb_kernel damped_coulomb_kernel proc~elec_prop_d2d->proc~damped_coulomb_kernel proc~ommp_message ommp_message proc~fatal_error->proc~ommp_message proc~close_output close_output proc~fatal_error->proc~close_output proc~prepare_fmm_ipd->proc~fatal_error proc~prepare_fmm_ext_ipd prepare_fmm_ext_ipd proc~prepare_fmm_ipd->proc~prepare_fmm_ext_ipd proc~fmm_l2l fmm_l2l proc~cart_propfar_at_ipart->proc~fmm_l2l proc~ntot_sph_harm ntot_sph_harm proc~cart_propfar_at_ipart->proc~ntot_sph_harm proc~damped_coulomb_kernel->proc~fatal_error proc~coulomb_kernel coulomb_kernel proc~damped_coulomb_kernel->proc~coulomb_kernel interface~mallocate mallocate proc~prepare_fmm_ext_ipd->interface~mallocate proc~fmm_solve_for_multipoles fmm_solve_for_multipoles proc~prepare_fmm_ext_ipd->proc~fmm_solve_for_multipoles interface~mfree mfree proc~prepare_fmm_ext_ipd->interface~mfree proc~make_vfact make_vfact proc~fmm_l2l->proc~make_vfact proc~fmm_l2l_rotation_work fmm_l2l_rotation_work proc~fmm_l2l->proc~fmm_l2l_rotation_work proc~close_output->proc~ommp_message proc~coulomb_kernel->proc~fatal_error proc~i_alloc1 i_alloc1 interface~mallocate->proc~i_alloc1 proc~i_alloc3 i_alloc3 interface~mallocate->proc~i_alloc3 proc~l_alloc2 l_alloc2 interface~mallocate->proc~l_alloc2 proc~r_alloc1 r_alloc1 interface~mallocate->proc~r_alloc1 proc~r_alloc3 r_alloc3 interface~mallocate->proc~r_alloc3 proc~i_alloc2 i_alloc2 interface~mallocate->proc~i_alloc2 proc~r_alloc2 r_alloc2 interface~mallocate->proc~r_alloc2 proc~l_alloc1 l_alloc1 interface~mallocate->proc~l_alloc1 proc~fmm_solve_for_multipoles->proc~fatal_error proc~tree_p2m tree_p2m proc~fmm_solve_for_multipoles->proc~tree_p2m proc~tree_m2m tree_m2m proc~fmm_solve_for_multipoles->proc~tree_m2m proc~tree_l2l tree_l2l proc~fmm_solve_for_multipoles->proc~tree_l2l proc~tree_m2l tree_m2l proc~fmm_solve_for_multipoles->proc~tree_m2l proc~time_pull time_pull proc~fmm_solve_for_multipoles->proc~time_pull proc~time_push time_push proc~fmm_solve_for_multipoles->proc~time_push proc~i_free2 i_free2 interface~mfree->proc~i_free2 proc~r_free3 r_free3 interface~mfree->proc~r_free3 proc~r_free2 r_free2 interface~mfree->proc~r_free2 proc~l_free1 l_free1 interface~mfree->proc~l_free1 proc~l_free2 l_free2 interface~mfree->proc~l_free2 proc~i_free3 i_free3 interface~mfree->proc~i_free3 proc~i_free1 i_free1 interface~mfree->proc~i_free1 proc~r_free1 r_free1 interface~mfree->proc~r_free1 proc~fmm_sph_rotate_oz_adj_work fmm_sph_rotate_oz_adj_work proc~fmm_l2l_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~carttosph carttosph proc~fmm_l2l_rotation_work->proc~carttosph proc~fmm_sph_rotate_oxz_work fmm_sph_rotate_oxz_work proc~fmm_l2l_rotation_work->proc~fmm_sph_rotate_oxz_work proc~trgev trgev proc~fmm_l2l_rotation_work->proc~trgev proc~fmm_sph_rotate_oz_work fmm_sph_rotate_oz_work proc~fmm_l2l_rotation_work->proc~fmm_sph_rotate_oz_work proc~fmm_l2l_ztranslate_work fmm_l2l_ztranslate_work proc~fmm_l2l_rotation_work->proc~fmm_l2l_ztranslate_work proc~chk_free chk_free proc~i_free2->proc~chk_free proc~r_free3->proc~chk_free proc~chk_alloc chk_alloc proc~i_alloc1->proc~chk_alloc proc~memory_init memory_init proc~i_alloc1->proc~memory_init proc~i_alloc3->proc~chk_alloc proc~i_alloc3->proc~memory_init proc~l_alloc2->proc~chk_alloc proc~l_alloc2->proc~memory_init proc~tree_p2m->proc~ntot_sph_harm proc~fmm_error fmm_error proc~tree_p2m->proc~fmm_error proc~fmm_m2m fmm_m2m proc~tree_p2m->proc~fmm_m2m proc~r_free2->proc~chk_free proc~l_free1->proc~chk_free proc~l_free2->proc~chk_free proc~tree_m2m->proc~ntot_sph_harm proc~tree_m2m->proc~fmm_m2m proc~i_free3->proc~chk_free proc~i_free1->proc~chk_free proc~tree_l2l->proc~fmm_l2l proc~tree_l2l->proc~ntot_sph_harm proc~tree_m2l->proc~ntot_sph_harm proc~fmm_m2l fmm_m2l proc~tree_m2l->proc~fmm_m2l proc~time_pull->proc~fatal_error proc~time_pull->proc~ommp_message proc~mem_stat mem_stat proc~time_pull->proc~mem_stat proc~r_free1->proc~chk_free proc~r_alloc1->proc~chk_alloc proc~r_alloc1->proc~memory_init proc~r_alloc3->proc~chk_alloc proc~r_alloc3->proc~memory_init proc~i_alloc2->proc~chk_alloc proc~i_alloc2->proc~memory_init proc~r_alloc2->proc~chk_alloc proc~r_alloc2->proc~memory_init proc~l_alloc1->proc~chk_alloc proc~l_alloc1->proc~memory_init proc~time_push->proc~fatal_error proc~time_push->proc~mem_stat proc~chk_free->proc~fatal_error proc~chk_alloc->proc~fatal_error proc~fmm_m2m_rotation_work fmm_m2m_rotation_work proc~fmm_m2m->proc~fmm_m2m_rotation_work proc~fmm_m2l_rotation_work fmm_m2l_rotation_work proc~fmm_m2l->proc~fmm_m2l_rotation_work proc~mem_stat->proc~memory_init proc~fmm_m2m_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~fmm_m2m_rotation_work->proc~carttosph proc~fmm_m2m_rotation_work->proc~fmm_sph_rotate_oxz_work proc~fmm_m2m_rotation_work->proc~trgev proc~fmm_m2m_rotation_work->proc~fmm_sph_rotate_oz_work proc~fmm_m2m_ztranslate_work fmm_m2m_ztranslate_work proc~fmm_m2m_rotation_work->proc~fmm_m2m_ztranslate_work proc~fmm_m2l_rotation_work->proc~fmm_sph_rotate_oz_adj_work proc~fmm_m2l_rotation_work->proc~carttosph proc~fmm_m2l_rotation_work->proc~fmm_sph_rotate_oxz_work proc~fmm_m2l_rotation_work->proc~trgev proc~fmm_m2l_rotation_work->proc~fmm_sph_rotate_oz_work proc~fmm_m2l_ztranslate_work fmm_m2l_ztranslate_work proc~fmm_m2l_rotation_work->proc~fmm_m2l_ztranslate_work

Called by

proc~~elec_prop_d2d~~CalledByGraph proc~elec_prop_d2d elec_prop_D2D proc~prepare_polelec prepare_polelec proc~prepare_polelec->proc~elec_prop_d2d proc~polelec_geomgrad polelec_geomgrad proc~polelec_geomgrad->proc~prepare_polelec proc~ommp_get_polelec_energy ommp_get_polelec_energy proc~ommp_get_polelec_energy->proc~prepare_polelec proc~energy_mm_pol energy_MM_pol proc~ommp_get_polelec_energy->proc~energy_mm_pol proc~ommp_set_external_field ommp_set_external_field proc~ommp_set_external_field->proc~prepare_polelec proc~energy_mm_pol->proc~prepare_polelec proc~ommp_polelec_geomgrad ommp_polelec_geomgrad proc~ommp_polelec_geomgrad->proc~polelec_geomgrad proc~c_ommp_get_polelec_energy C_ommp_get_polelec_energy proc~c_ommp_get_polelec_energy->proc~ommp_get_polelec_energy proc~ommp_get_full_ele_energy ommp_get_full_ele_energy proc~ommp_get_full_ele_energy->proc~ommp_get_polelec_energy proc~ommp_set_external_field_nomm ommp_set_external_field_nomm proc~ommp_set_external_field_nomm->proc~ommp_set_external_field proc~ommp_full_geomgrad ommp_full_geomgrad proc~ommp_full_geomgrad->proc~polelec_geomgrad proc~c_ommp_set_external_field_nomm C_ommp_set_external_field_nomm proc~c_ommp_set_external_field_nomm->proc~ommp_set_external_field proc~c_ommp_set_external_field C_ommp_set_external_field proc~c_ommp_set_external_field->proc~ommp_set_external_field proc~c_ommp_polelec_geomgrad C_ommp_polelec_geomgrad proc~c_ommp_polelec_geomgrad->proc~ommp_polelec_geomgrad proc~c_ommp_full_geomgrad C_ommp_full_geomgrad proc~c_ommp_full_geomgrad->proc~ommp_full_geomgrad proc~c_ommp_get_full_ele_energy C_ommp_get_full_ele_energy proc~c_ommp_get_full_ele_energy->proc~ommp_get_full_ele_energy proc~ommp_get_full_energy ommp_get_full_energy proc~ommp_get_full_energy->proc~ommp_get_full_ele_energy proc~c_ommp_get_full_energy C_ommp_get_full_energy proc~c_ommp_get_full_energy->proc~ommp_get_full_energy

Contents

Source Code


Source Code

    subroutine elec_prop_D2D(eel, in_kind, do_V, do_E, do_Egrd, do_EHes)
        !! Computes the electric field of a trial set of induced point dipoles
        !! at polarizable sites. This is intended to be used as matrix-vector
        !! routine in the solution of the linear system.
        implicit none

        type(ommp_electrostatics_type), intent(inout) :: eel
        !! Data structure for electrostatic part of the system
        logical, intent(in) :: do_V, do_E, do_Egrd, do_EHes
        !! Flag to control which properties have to be computed.
        character, intent(in) :: in_kind

        integer(ip) :: i, j, jpol, ipol, ij, idx, ikernel, knd
        logical :: to_scale, to_do
        real(rp) :: kernel(5), dr(3), tmpV, tmpE(3), tmpEgr(6), tmpHE(10), scalf

        knd = 1 ! Default
        if(in_kind == 'P') then
            knd = _amoeba_P_
        elseif(in_kind == 'D') then 
            knd = _amoeba_D_
        elseif(eel%amoeba) then
            call fatal_error("Unrecognized interaction '"//in_kind//"' in elec&
                             &_prop_D2D.")
        else
            knd = 1
        end if

        if(do_EHes) then
            ikernel = 5 
        elseif(do_Egrd) then
            ikernel = 4
        elseif(do_E) then
            ikernel = 3
        elseif(do_V) then
            ikernel = 2
        else
            return
        end if

        if(eel%use_fmm) then
            call prepare_fmm_ipd(eel, knd)


            !$omp parallel do default(shared) schedule(dynamic) &
            !$omp private(i,j,ij,ipol,jpol,idx,dr,kernel,to_do,to_scale,scalf,tmpV,tmpE,tmpEgr,tmpHE) 
            do ipol=1, eel%pol_atoms 
                i = eel%polar_mm(ipol)
                
                call cart_propfar_at_ipart(eel%fmm_ipd(knd), i, &
                                           do_V, eel%V_D2D(ipol,knd), &
                                           do_E, eel%E_D2D(:,ipol,knd), &
                                           do_Egrd, eel%Egrd_D2D(:,ipol,knd), &
                                           do_EHes, eel%EHes_D2D(:,ipol,knd))

                ! Near field is computed internally because dumped kernel is required
                do ij=eel%fmm_near_field_list%ri(i), eel%fmm_near_field_list%ri(i+1)-1
                    j = eel%fmm_near_field_list%ci(ij)
                    jpol = eel%polar_mm(j)
                    ! If the atom is not polarizable, skip
                    if(jpol < 1) cycle 

                    !loop on target
                    to_do = .true.
                    to_scale = .false.
                    scalf = 1.0

                    ! Check if the element should be scaled
                    do idx=eel%list_P_P%ri(ipol), eel%list_P_P%ri(ipol+1)-1
                        if(eel%list_P_P%ci(idx) == jpol) then
                            to_scale = .true.
                            exit
                        end if
                    end do

                    !If it should set the correct variables
                    if(to_scale) then
                        to_do = eel%todo_P_P(idx)
                        scalf = eel%scalef_P_P(idx)
                    end if
                    
                    if(to_do) then
                        call damped_coulomb_kernel(eel, j, i,& 
                                                   ikernel, kernel, dr)
                        

                        if(do_V) tmpV = 0.0_rp
                        if(do_E) tmpE = 0.0_rp
                        if(do_Egrd) tmpEgr = 0.0_rp
                        if(do_EHes) tmpHE = 0.0_rp

                        call mu_elec_prop(eel%ipd(:,jpol,knd), dr, kernel, &
                                        do_V, tmpV, &
                                        do_E, tmpE, &
                                        do_Egrd, tmpEgr, & 
                                        do_EHes, tmpHE)
                        if(to_scale) then
                            if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV * scalf
                            if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE * scalf
                            if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr * scalf
                            if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE * scalf
                        else
                            if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV
                            if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE
                            if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr
                            if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE
                        end if
                    end if
                end do
            end do

            if(allocated(eel%list_P_P_fmm_far)) then
                ! Now remove screened interactions from far-field
                do ipol=1, eel%pol_atoms
                    i = eel%polar_mm(ipol)
                    
                    do idx=eel%list_P_P_fmm_far%ri(ipol), eel%list_P_P_fmm_far%ri(ipol+1)-1
                        jpol = eel%list_P_P_fmm_far%ci(idx)
                        j = eel%polar_mm(jpol)

                        scalf = 1.0 - eel%scalef_P_P_fmm_far(idx)
                        
                        call damped_coulomb_kernel(eel, j, i,& 
                                                    2, kernel(1:3), dr)
                        
                        if(do_V) tmpV = 0.0_rp
                        if(do_E) tmpE = 0.0_rp
                        if(do_Egrd) tmpEgr = 0.0_rp
                        if(do_EHes) tmpHE = 0.0_rp

                        call mu_elec_prop(eel%ipd(:,jpol,knd), dr, kernel, &
                                        do_V, tmpV, &
                                        do_E, tmpE, &
                                        do_Egrd, tmpEgr, & 
                                        do_EHes, tmpHE)

                        if(do_V) eel%V_D2D(ipol,knd) = eel%V_D2D(ipol,knd) + tmpV * scalf
                        if(do_E) eel%E_D2D(:, ipol,knd) = eel%E_D2D(:, ipol,knd) + tmpE * scalf
                        if(do_Egrd) eel%Egrd_D2D(:, ipol,knd) = eel%Egrd_D2D(:, ipol,knd) + tmpEgr * scalf
                        if(do_EHes) eel%EHes_D2D(:, ipol,knd) = eel%EHes_D2D(:, ipol,knd) + tmpHE * scalf
                    end do
                end do
            end if
        else
        !$omp parallel do default(shared) schedule(dynamic) &
        !$omp private(i,j,idx,dr,kernel,to_do,to_scale,scalf,tmpV,tmpE,tmpEgr,tmpHE) 
        do j=1, eel%pol_atoms
            do i=1, eel%pol_atoms
                if(j == i) cycle
                !loop on target
                to_do = .true.
                to_scale = .false.
                scalf = 1.0

                ! Check if the element should be scaled
                do idx=eel%list_P_P%ri(i), eel%list_P_P%ri(i+1)-1
                    if(eel%list_P_P%ci(idx) == j) then
                        to_scale = .true.
                        exit
                    end if
                end do

                !If it should set the correct variables
                if(to_scale) then
                    to_do = eel%todo_P_P(idx)
                    scalf = eel%scalef_P_P(idx)
                end if
                
                if(to_do) then
                    call damped_coulomb_kernel(eel, eel%polar_mm(i), &
                                               eel%polar_mm(j),& 
                                               ikernel, kernel, dr)
                    
                    if(do_V) tmpV = 0.0_rp
                    if(do_E) tmpE = 0.0_rp
                    if(do_Egrd) tmpEgr = 0.0_rp
                    if(do_EHes) tmpHE = 0.0_rp

                    call mu_elec_prop(eel%ipd(:,i,knd), dr, kernel, &
                                      do_V, tmpV, &
                                      do_E, tmpE, &
                                      do_Egrd, tmpEgr, & 
                                      do_EHes, tmpHE)
                    if(to_scale) then
                        if(do_V) eel%V_D2D(j,knd) = eel%V_D2D(j,knd) + tmpV * scalf
                        if(do_E) eel%E_D2D(:, j,knd) = eel%E_D2D(:, j,knd) + tmpE * scalf
                        if(do_Egrd) eel%Egrd_D2D(:, j,knd) = eel%Egrd_D2D(:, j,knd) + tmpEgr * scalf
                        if(do_EHes) eel%EHes_D2D(:, j,knd) = eel%EHes_D2D(:, j,knd) + tmpHE * scalf
                    else
                        if(do_V) eel%V_D2D(j,knd) = eel%V_D2D(j,knd) + tmpV
                        if(do_E) eel%E_D2D(:, j,knd) = eel%E_D2D(:, j,knd) + tmpE
                        if(do_Egrd) eel%Egrd_D2D(:, j,knd) = eel%Egrd_D2D(:, j,knd) + tmpEgr
                        if(do_EHes) eel%EHes_D2D(:, j,knd) = eel%EHes_D2D(:, j,knd) + tmpHE
                    end if
                end if
            end do
        end do
        end if
    end subroutine elec_prop_D2D