ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 459
Committed: Fri Apr 4 19:57:01 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 6333 byte(s)
Log Message:
fixed a memory read bug in neighborlist

File Contents

# User Rev Content
1 mmeineke 377 module dipole_dipole
2    
3     use force_globals
4     use definitions
5     use atype_module
6     use vector_class
7     #ifdef IS_MPI
8     use mpiSimulation
9     #endif
10     implicit none
11    
12     real(kind=dp), save :: rrf
13     real(kind=dp), save :: rt
14     real(kind=dp), save :: rrfsq
15     real(kind=dp), save :: pre
16     logical, save :: dipole_initialized = .false.
17    
18     contains
19    
20     subroutine initialize_dipole(this_rrf, this_rt)
21     real(kind=dp), intent(in) :: this_rrf, this_rt
22     rrf = this_rrf
23     rt = this_rt
24     rrfsq = rrf * rrf
25     ! pre converts from mu in units of debye to kcal/mol
26 chuckv 388 pre = 14.38362_dp
27 gezelter 394
28 mmeineke 377 dipole_initialized = .true.
29    
30     return
31     end subroutine initialize_dipole
32    
33    
34 mmeineke 459 subroutine do_dipole_pair(natoms, atom1, atom2, d, rij, r2, pot, u_l, f, t, &
35 mmeineke 377 do_pot, do_stress)
36    
37     logical :: do_pot, do_stress
38    
39 mmeineke 459 integer :: natoms
40 mmeineke 377 integer atom1, atom2, me1, me2
41     real(kind=dp) :: rij, mu1, mu2
42     real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5, pre
43     real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
44     real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
45     real(kind=dp) :: taper, dtdr, vterm
46    
47     real( kind = dp ) :: pot
48     real( kind = dp ), dimension(3) :: d
49     real( kind = dp ), dimension(:,:) :: u_l
50     real( kind = dp ), dimension(:,:) :: f
51     real( kind = dp ), dimension(:,:) :: t
52    
53     real (kind = dp), dimension(3) :: ul1
54     real (kind = dp), dimension(3) :: ul2
55    
56     if (.not.dipole_initialized) then
57     write(default_error,*) 'Dipole-dipole not initialized!'
58     return
59     endif
60    
61     #ifdef IS_MPI
62     me1 = atid_Row(atom1)
63     ul1(1) = u_l_Row(1,atom1)
64     ul1(2) = u_l_Row(2,atom1)
65     ul1(3) = u_l_Row(3,atom1)
66    
67     me2 = atid_Col(atom2)
68     ul2(1) = u_l_Col(1,atom2)
69     ul2(2) = u_l_Col(2,atom2)
70     ul2(3) = u_l_Col(3,atom2)
71     #else
72     me1 = atid(atom1)
73     ul1(1) = u_l(1,atom1)
74     ul1(2) = u_l(2,atom1)
75     ul1(3) = u_l(3,atom1)
76    
77     me2 = atid(atom2)
78     ul2(1) = u_l(1,atom2)
79     ul2(2) = u_l(2,atom2)
80     ul2(3) = u_l(3,atom2)
81     #endif
82    
83 mmeineke 443 if( atom1 .eq. 2 )then
84     write (*,*) 'ul =', ul1(1), ul1(2), ul1(3)
85     endif
86    
87     if( atom2 .eq. 2 )then
88     write (*,*) 'ul =', ul2(1), ul2(2), ul2(3)
89     endif
90    
91    
92 mmeineke 377 call getElementProperty(atypes, me1, "dipole_moment", mu1)
93     call getElementProperty(atypes, me2, "dipole_moment", mu2)
94 chuckv 388
95 mmeineke 377 if (rij.le.rrf) then
96    
97     if (rij.lt.rt) then
98     taper = 1.0d0
99     dtdr = 0.0d0
100     else
101     taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
102     dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
103     endif
104    
105     r3 = r2*rij
106     r5 = r3*r2
107    
108     rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
109     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
110     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
111    
112     dip2 = pre * mu1 * mu2
113    
114     dfact1 = 3.0d0*dip2 / r2
115     dfact2 = 3.0d0*dip2 / r5
116    
117     vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
118    
119     if (do_pot) then
120     #ifdef IS_MPI
121     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
122     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
123     #else
124     pot = pot + vterm*taper
125     #endif
126     endif
127    
128     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
129     (5.0d0*(rdotu1*rdotu2)/r5)) - &
130     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
131     vterm*dtdr*d(1)/rij
132    
133     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
134     (5.0d0*(rdotu1*rdotu2)/r5)) - &
135     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
136     vterm*dtdr*d(2)/rij
137    
138     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
139     (5.0d0*(rdotu1*rdotu2)/r5)) - &
140     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
141     vterm*dtdr*d(3)/rij
142    
143     dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
144     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
145     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
146    
147     dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
148     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
149     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
150    
151    
152     #ifdef IS_MPI
153 gezelter 394 f_Row(1,atom1) = f_Row(1,atom1) - dudx
154     f_Row(2,atom1) = f_Row(2,atom1) - dudy
155     f_Row(3,atom1) = f_Row(3,atom1) - dudz
156 mmeineke 377
157 gezelter 394 f_Col(1,atom2) = f_Col(1,atom2) + dudx
158     f_Col(2,atom2) = f_Col(2,atom2) + dudy
159     f_Col(3,atom2) = f_Col(3,atom2) + dudz
160 mmeineke 377
161     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
162     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
163     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
164    
165     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
166     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
167     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
168     #else
169 gezelter 394 f(1,atom1) = f(1,atom1) - dudx
170     f(2,atom1) = f(2,atom1) - dudy
171     f(3,atom1) = f(3,atom1) - dudz
172 mmeineke 377
173 gezelter 394 f(1,atom2) = f(1,atom2) + dudx
174     f(2,atom2) = f(2,atom2) + dudy
175     f(3,atom2) = f(3,atom2) + dudz
176 mmeineke 377
177     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
178     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
179     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
180    
181     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
182     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
183     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
184     #endif
185    
186     if (do_stress) then
187     tau_Temp(1) = tau_Temp(1) + dudx * d(1)
188     tau_Temp(2) = tau_Temp(2) + dudx * d(2)
189     tau_Temp(3) = tau_Temp(3) + dudx * d(3)
190     tau_Temp(4) = tau_Temp(4) + dudy * d(1)
191     tau_Temp(5) = tau_Temp(5) + dudy * d(2)
192     tau_Temp(6) = tau_Temp(6) + dudy * d(3)
193     tau_Temp(7) = tau_Temp(7) + dudz * d(1)
194     tau_Temp(8) = tau_Temp(8) + dudz * d(2)
195     tau_Temp(9) = tau_Temp(9) + dudz * d(3)
196     virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
197     endif
198    
199     endif
200    
201     return
202     end subroutine do_dipole_pair
203    
204     end module dipole_dipole