ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_dipole_dipole.F90
Revision: 312
Committed: Tue Mar 11 17:46:18 2003 UTC (21 years, 3 months ago) by gezelter
File size: 5387 byte(s)
Log Message:
Bunch o' stuff, particularly the vector_class.F90 module

File Contents

# Content
1 module dipole_dipole
2
3 use simulation
4 use definitions
5 use forceGlobals
6 use atype_typedefs
7 #ifdef IS_MPI
8 use mpiSimulation
9 #endif
10 implicit none
11
12 contains
13
14 subroutine do_dipole_pair(atom1, atom2, atype1, atype2, dx, dy, dz, rij, &
15 pot, u_l, f, t)
16
17 integer atom1, atom2
18 double precision dx, dy, dz, rij
19 double precision dfact1, dfact2, dip2, r2, r3, r5, pre
20 double precision dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
21 double precision dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
22 double precision taper, dtdr, vterm
23
24 real( kind = dp ) :: pot, rt, rrf
25 real( kind = dp ), dimension(3,getNlocal()) :: u_l
26 real( kind = dp ), dimension(3,getNlocal()) :: f
27 real( kind = dp ), dimension(3,getNlocal()) :: t
28
29 real (kind = dp), dimension(3) :: ul1
30 real (kind = dp), dimension(3) :: ul2
31 type (atype), pointer :: atype1
32 type (atype), pointer :: atype2
33
34 #ifdef IS_MPI
35 ul1(1) = u_l_Row(1,atom1)
36 ul1(2) = u_l_Row(2,atom1)
37 ul1(3) = u_l_Row(3,atom1)
38
39 ul2(1) = u_l_Col(1,atom2)
40 ul2(2) = u_l_Col(2,atom2)
41 ul2(3) = u_l_Col(3,atom2)
42 #else
43 ul1(1) = u_l(1,atom1)
44 ul1(2) = u_l(2,atom1)
45 ul1(3) = u_l(3,atom1)
46
47 ul1(1) = u_l(1,atom2)
48 ul1(2) = u_l(2,atom2)
49 ul1(3) = u_l(3,atom2)
50 #endif
51
52 rrf = getRrf()
53 rt = getRt()
54
55 ! pre converts from mu in units of debye to kcal/mol
56 pre = 14.38362d0
57
58 if (rij.le.rrf) then
59
60 if (rij.lt.rt) then
61 taper = 1.0d0
62 dtdr = 0.0d0
63 else
64 taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
65 dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
66 endif
67
68 r2 = rij*rij
69 r3 = r2*rij
70 r5 = r3*r2
71
72 rdotu1 = dx*ul1(1) + dy*ul1(2) + dz*ul1(3)
73 rdotu2 = dx*ul2(1) + dy*ul2(2) + dz*ul2(3)
74 u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
75
76 dip2 = pre * atype1%dipole_moment * atype2%dipole_moment
77
78 dfact1 = 3.0d0*dip2 / r2
79 dfact2 = 3.0d0*dip2 / r5
80
81 vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
82
83 #ifdef IS_MPI
84 pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
85 pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
86 #else
87 pot = pot + vterm*taper
88 #endif
89
90 dudx = (-dfact1 * dx * ((u1dotu2/r3) - &
91 (5.0d0*(rdotu1*rdotu2)/r5)) - &
92 dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
93 vterm*dtdr*dx/rij
94
95 dudy = (-dfact1 * dy * ((u1dotu2/r3) - &
96 (5.0d0*(rdotu1*rdotu2)/r5)) - &
97 dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
98 vterm*dtdr*dy/rij
99
100 dudz = (-dfact1 * dz * ((u1dotu2/r3) - &
101 (5.0d0*(rdotu1*rdotu2)/r5)) - &
102 dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
103 vterm*dtdr*dz/rij
104
105 dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*dx*rdotu2/r5)))*taper
106 dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*dy*rdotu2/r5)))*taper
107 dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*dz*rdotu2/r5)))*taper
108
109 dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*dx*rdotu1/r5)))*taper
110 dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*dy*rdotu1/r5)))*taper
111 dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*dz*rdotu1/r5)))*taper
112
113
114 #ifdef IS_MPI
115 f_Row(1,atom1) = f_Row(1,atom1) + dudx
116 f_Row(2,atom1) = f_Row(2,atom1) + dudy
117 f_Row(3,atom1) = f_Row(3,atom1) + dudz
118
119 f_Col(1,atom2) = f_Col(1,atom2) - dudx
120 f_Col(2,atom2) = f_Col(2,atom2) - dudy
121 f_Col(3,atom2) = f_Col(3,atom2) - dudz
122
123 t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
124 t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
125 t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
126
127 t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
128 t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
129 t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
130 #else
131 f(1,atom1) = f(1,atom1) + dudx
132 f(2,atom1) = f(2,atom1) + dudy
133 f(3,atom1) = f(3,atom1) + dudz
134
135 f(1,atom2) = f(1,atom2) - dudx
136 f(2,atom2) = f(2,atom2) - dudy
137 f(3,atom2) = f(3,atom2) - dudz
138
139 t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
140 t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
141 t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
142
143 t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
144 t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
145 t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
146 #endif
147
148 if (doStress()) then
149 tau_Temp(1) = tau_Temp(1) + dudx * dx
150 tau_Temp(2) = tau_Temp(2) + dudx * dy
151 tau_Temp(3) = tau_Temp(3) + dudx * dz
152 tau_Temp(4) = tau_Temp(4) + dudy * dx
153 tau_Temp(5) = tau_Temp(5) + dudy * dy
154 tau_Temp(6) = tau_Temp(6) + dudy * dz
155 tau_Temp(7) = tau_Temp(7) + dudz * dx
156 tau_Temp(8) = tau_Temp(8) + dudz * dy
157 tau_Temp(9) = tau_Temp(9) + dudz * dz
158 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
159 endif
160
161 endif
162
163 return
164 end subroutine do_dipole_pair
165
166 end module dipole_dipole