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: 317
Committed: Tue Mar 11 23:13:06 2003 UTC (21 years, 4 months ago) by gezelter
File size: 5593 byte(s)
Log Message:
Bug fixes

File Contents

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