dipole_T Subroutine

private subroutine dipole_T(eel, i, j, tens)

Uses

  • proc~~dipole_t~~UsesGraph proc~dipole_t dipole_T module~mod_electrostatics mod_electrostatics proc~dipole_t->module~mod_electrostatics module~mod_constants mod_constants proc~dipole_t->module~mod_constants module~mod_electrostatics->module~mod_constants module~mod_memory mod_memory module~mod_electrostatics->module~mod_memory module~mod_profiling mod_profiling module~mod_electrostatics->module~mod_profiling module~mod_topology mod_topology module~mod_electrostatics->module~mod_topology module~mod_io mod_io module~mod_electrostatics->module~mod_io module~mod_adjacency_mat mod_adjacency_mat module~mod_electrostatics->module~mod_adjacency_mat iso_c_binding iso_c_binding module~mod_constants->iso_c_binding module~mod_memory->module~mod_constants module~mod_memory->module~mod_io module~mod_memory->iso_c_binding module~mod_profiling->module~mod_constants module~mod_profiling->module~mod_memory module~mod_profiling->module~mod_io module~mod_topology->module~mod_memory module~mod_topology->module~mod_adjacency_mat module~mod_io->module~mod_constants module~mod_adjacency_mat->module~mod_memory

This subroutine compute the interaction tensor (rank 3) between two polarizable sites i and j. This tensor is built according to the following rules: ... TODO

Arguments

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

The electostatic data structure for which the interaction tensor should be computed

integer(kind=ip), intent(in) :: i

Index (in the list of polarizable sites) of the source site

integer(kind=ip), intent(in) :: j

Index (in the list of polarizable sites) of the target site

real(kind=rp), intent(out), dimension(3, 3) :: tens

Interaction tensor between sites i and j


Calls

proc~~dipole_t~~CallsGraph proc~dipole_t dipole_T proc~screening_rules screening_rules proc~dipole_t->proc~screening_rules proc~damped_coulomb_kernel damped_coulomb_kernel proc~dipole_t->proc~damped_coulomb_kernel proc~fatal_error fatal_error proc~screening_rules->proc~fatal_error proc~damped_coulomb_kernel->proc~fatal_error proc~coulomb_kernel coulomb_kernel proc~damped_coulomb_kernel->proc~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~coulomb_kernel->proc~fatal_error proc~close_output->proc~ommp_message

Called by

proc~~dipole_t~~CalledByGraph proc~dipole_t dipole_T proc~create_tmat create_tmat proc~create_tmat->proc~dipole_t proc~polarization polarization proc~polarization->proc~create_tmat proc~ommp_get_polelec_energy ommp_get_polelec_energy proc~ommp_get_polelec_energy->proc~polarization proc~ommp_set_external_field ommp_set_external_field proc~ommp_set_external_field->proc~polarization proc~polelec_geomgrad polelec_geomgrad proc~polelec_geomgrad->proc~polarization proc~ommp_get_full_ele_energy ommp_get_full_ele_energy proc~ommp_get_full_ele_energy->proc~ommp_get_polelec_energy program~test_si_potential test_SI_potential program~test_si_potential->proc~ommp_get_polelec_energy program~test_si_potential->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~ommp_set_external_field_nomm ommp_set_external_field_nomm proc~ommp_set_external_field_nomm->proc~ommp_set_external_field proc~c_ommp_get_polelec_energy C_ommp_get_polelec_energy proc~c_ommp_get_polelec_energy->proc~ommp_get_polelec_energy 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~ommp_polelec_geomgrad ommp_polelec_geomgrad proc~ommp_polelec_geomgrad->proc~polelec_geomgrad proc~ommp_full_geomgrad ommp_full_geomgrad proc~ommp_full_geomgrad->proc~polelec_geomgrad proc~ommp_get_full_energy ommp_get_full_energy proc~ommp_get_full_energy->proc~ommp_get_full_ele_energy proc~c_ommp_polelec_geomgrad C_ommp_polelec_geomgrad proc~c_ommp_polelec_geomgrad->proc~ommp_polelec_geomgrad proc~ommptest_fakeqm_internal_geomgrad ommptest_fakeqm_internal_geomgrad proc~ommptest_fakeqm_internal_geomgrad->proc~ommp_full_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~ommptest_totalqmmm_geomgrad ommptest_totalqmmm_geomgrad proc~ommptest_totalqmmm_geomgrad->proc~ommp_full_geomgrad proc~ommptest_fakeqm_linkatom_geomgrad ommptest_fakeqm_linkatom_geomgrad proc~ommptest_fakeqm_linkatom_geomgrad->proc~ommp_full_geomgrad 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 dipole_T(eel, i, j, tens)
        !! This subroutine compute the interaction tensor (rank 3) between
        !! two polarizable sites i and j.
        !! This tensor is built according to the following rules: ... TODO
        use mod_electrostatics, only: screening_rules, damped_coulomb_kernel
        use mod_constants, only: eps_rp

        implicit none
        !                      
        ! Compute element of the polarization tensor TTens between
        ! polarizable cpol atom I and polarizable cpol atom J. On the
        ! TTens diagonal (I=J) are inverse polarizabilities and on the 
        ! off-diagonal dipole field.
        !
        ! Polarizabilities pol are defined for polarizable atoms only while 
        ! Thole factors are defined for all of them
        
        type(ommp_electrostatics_type), intent(inout) :: eel
        !! The electostatic data structure  for which the 
        !! interaction tensor should be computed
        integer(ip), intent(in) :: i 
        !! Index (in the list of polarizable sites) of the source site
        integer(ip), intent(in) :: j
        !! Index (in the list of polarizable sites) of the target site
        real(rp), dimension(3, 3), intent(out) :: tens
        !! Interaction tensor between sites i and j
        
        real(rp) :: dr(3)
        real(rp) ::  kernel(3), scalf
        logical :: to_do, to_scale

        integer(ip) :: ii, jj
         
        tens = 0.0_rp
        
        if(i == j) then
            tens(1, 1) = 1.0_rp / eel%pol(i)
            tens(2, 2) = 1.0_rp / eel%pol(i)
            tens(3, 3) = 1.0_rp / eel%pol(i)
        else
            scalf = screening_rules(eel, i, 'P', j, 'P', '-')

            if(abs(scalf) > eps_rp) then
                call damped_coulomb_kernel(eel, eel%polar_mm(i), &
                                           eel%polar_mm(j), 2, kernel, dr)
                ! Fill the matrix elemets
                do ii=1, 3
                    do jj=1, 3
                        if(ii == jj) then
                            tens(ii, ii) = kernel(2) - 3.0_rp * kernel(3) * dr(ii) ** 2
                        else
                            tens(jj, ii) = -3.0_rp * kernel(3) * dr(ii) * dr(jj)
                        end if
                    end do
                end do
                ! Scale if needed
                if(abs(scalf-1.0) > eps_rp) tens = tens * scalf
            
            end if
        end if
    end subroutine dipole_T