--- trunk/OOPSE/libmdtools/calc_sticky_pair.F90 2003/03/21 17:42:12 378 +++ trunk/OOPSE/libmdtools/calc_sticky_pair.F90 2003/04/09 04:06:43 483 @@ -9,12 +9,13 @@ !! @author Matthew Meineke !! @author Christopher Fennel !! @author J. Daniel Gezelter -!! @version $Id: calc_sticky_pair.F90,v 1.1.1.1 2003-03-21 17:42:12 mmeineke Exp $, $Date: 2003-03-21 17:42:12 $, $Name: not supported by cvs2svn $, $Revision: 1.1.1.1 $ +!! @version $Id: calc_sticky_pair.F90,v 1.8 2003-04-09 04:06:43 gezelter Exp $, $Date: 2003-04-09 04:06:43 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $ module sticky_pair use force_globals use definitions + use simulation #ifdef IS_MPI use mpiSimulation #endif @@ -24,11 +25,11 @@ module sticky_pair PRIVATE logical, save :: sticky_initialized = .false. - real( kind = dp ), save :: SSD_w0 - real( kind = dp ), save :: SSD_v0 - real( kind = dp ), save :: SSD_rl - real( kind = dp ), save :: SSD_ru - real( kind = dp ), save :: SSD_rup + real( kind = dp ), save :: SSD_w0 = 0.0_dp + real( kind = dp ), save :: SSD_v0 = 0.0_dp + real( kind = dp ), save :: SSD_rl = 0.0_dp + real( kind = dp ), save :: SSD_ru = 0.0_dp + real( kind = dp ), save :: SSD_rup = 0.0_dp public :: check_sticky_FF public :: set_sticky_params @@ -74,9 +75,9 @@ contains real (kind=dp), intent(inout) :: rij, r2 real (kind=dp), dimension(3), intent(in) :: d real (kind=dp) :: pot - real (kind=dp), dimension(:,:) :: A - real (kind=dp), dimension(:,:) :: f - real (kind=dp), dimension(:,:) :: t + real (kind=dp), dimension(9,getNlocal()) :: A + real (kind=dp), dimension(3,getNlocal()) :: f + real (kind=dp), dimension(3,getNlocal()) :: t logical, intent(in) :: do_pot, do_stress real (kind=dp) :: xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2 @@ -92,6 +93,7 @@ contains real (kind=dp) :: fxii, fyii, fzii, fxjj, fyjj, fzjj real (kind=dp) :: fxij, fyij, fzij, fxji, fyji, fzji real (kind=dp) :: fxradial, fyradial, fzradial + real (kind=dp) :: rijtest, rjitest if (.not.sticky_initialized) then write(*,*) 'Sticky forces not initialized!' @@ -140,7 +142,7 @@ contains xj2 = xj*xj yj2 = yj*yj zj2 = zj*zj - + call calc_sw_fnc(rij, s, sp, dsdr, dspdr) wi = 2.0d0*(xi2-yi2)*zi / r3 @@ -157,7 +159,6 @@ contains wjp = zjf*zjf*zjs*zjs - SSD_w0 wp = wip + wjp - if (do_pot) then #ifdef IS_MPI pot_row(atom1) = pot_row(atom1) + 0.25d0*SSD_v0*(s*w + sp*wp) @@ -303,36 +304,32 @@ contains f_Row(2,atom1) = f_Row(2,atom1) + fyradial f_Row(3,atom1) = f_Row(3,atom1) + fzradial - f_Col(1,atom2) = f_Col(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - & - dspdr*drdx*wp + fxjj + fxij) - f_Col(2,atom2) = f_Col(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - & - dspdr*drdy*wp + fyjj + fyij) - f_Col(3,atom2) = f_Col(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - & - dspdr*drdz*wp + fzjj + fzij) + f_Col(1,atom2) = f_Col(1,atom2) - fxradial + f_Col(2,atom2) = f_Col(2,atom2) - fyradial + f_Col(3,atom2) = f_Col(3,atom2) - fzradial #else f(1,atom1) = f(1,atom1) + fxradial f(2,atom1) = f(2,atom1) + fyradial f(3,atom1) = f(3,atom1) + fzradial - f(1,atom2) = f(1,atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + & - fxjj + fxij) - f(2,atom2) = f(2,atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + & - fyjj + fyij) - f(3,atom2) = f(3,atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + & - fzjj + fzij) + f(1,atom2) = f(1,atom2) - fxradial + f(2,atom2) = f(2,atom2) - fyradial + f(3,atom2) = f(3,atom2) - fzradial #endif if (do_stress) then - tau_Temp(1) = tau_Temp(1) + fxradial * d(1) - tau_Temp(2) = tau_Temp(2) + fxradial * d(2) - tau_Temp(3) = tau_Temp(3) + fxradial * d(3) - tau_Temp(4) = tau_Temp(4) + fyradial * d(1) - tau_Temp(5) = tau_Temp(5) + fyradial * d(2) - tau_Temp(6) = tau_Temp(6) + fyradial * d(3) - tau_Temp(7) = tau_Temp(7) + fzradial * d(1) - tau_Temp(8) = tau_Temp(8) + fzradial * d(2) - tau_Temp(9) = tau_Temp(9) + fzradial * d(3) - virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) + if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then + tau_Temp(1) = tau_Temp(1) + fxradial * d(1) + tau_Temp(2) = tau_Temp(2) + fxradial * d(2) + tau_Temp(3) = tau_Temp(3) + fxradial * d(3) + tau_Temp(4) = tau_Temp(4) + fyradial * d(1) + tau_Temp(5) = tau_Temp(5) + fyradial * d(2) + tau_Temp(6) = tau_Temp(6) + fyradial * d(3) + tau_Temp(7) = tau_Temp(7) + fzradial * d(1) + tau_Temp(8) = tau_Temp(8) + fzradial * d(2) + tau_Temp(9) = tau_Temp(9) + fzradial * d(3) + virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) + endif endif end subroutine do_sticky_pair