ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 482
Committed: Tue Apr 8 22:38:43 2003 UTC (21 years, 3 months ago) by chuckv
File size: 6311 byte(s)
Log Message:
It works (kinda)...

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 chuckv 460 use simulation
8 mmeineke 377 #ifdef IS_MPI
9     use mpiSimulation
10     #endif
11     implicit none
12    
13 mmeineke 469 real(kind=dp), save :: rrf = 0.0
14     real(kind=dp), save :: rt = 0.0
15     real(kind=dp), save :: rrfsq = 0.0
16     real(kind=dp), save :: pre = 0.0
17 mmeineke 377 logical, save :: dipole_initialized = .false.
18    
19     contains
20    
21     subroutine initialize_dipole(this_rrf, this_rt)
22     real(kind=dp), intent(in) :: this_rrf, this_rt
23     rrf = this_rrf
24     rt = this_rt
25     rrfsq = rrf * rrf
26     ! pre converts from mu in units of debye to kcal/mol
27 chuckv 388 pre = 14.38362_dp
28 gezelter 394
29 mmeineke 377 dipole_initialized = .true.
30    
31     return
32     end subroutine initialize_dipole
33    
34    
35 chuckv 460 subroutine do_dipole_pair(atom1, atom2, d, rij, r2, pot, u_l, f, t, &
36 mmeineke 377 do_pot, do_stress)
37    
38     logical :: do_pot, do_stress
39    
40     integer atom1, atom2, me1, me2
41     real(kind=dp) :: rij, mu1, mu2
42 mmeineke 469 real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5
43 mmeineke 377 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 chuckv 460 real( kind = dp ), dimension(3,getNlocal()) :: u_l
50     real( kind = dp ), dimension(3,getNlocal()) :: f
51     real( kind = dp ), dimension(3,getNlocal()) :: t
52 mmeineke 377
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 mmeineke 469 write (0,*) 'ul =', ul1(1), ul1(2), ul1(3)
85 mmeineke 443 endif
86    
87     if( atom2 .eq. 2 )then
88 mmeineke 469 write (0,*) 'ul =', ul2(1), ul2(2), ul2(3)
89 mmeineke 443 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 mmeineke 469
105 mmeineke 377 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 mmeineke 469
112 mmeineke 377 dip2 = pre * mu1 * mu2
113     dfact1 = 3.0d0*dip2 / r2
114     dfact2 = 3.0d0*dip2 / r5
115 mmeineke 469
116 mmeineke 377 vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
117 mmeineke 469
118 mmeineke 377 if (do_pot) then
119     #ifdef IS_MPI
120     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
121     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
122     #else
123     pot = pot + vterm*taper
124     #endif
125     endif
126    
127     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
128     (5.0d0*(rdotu1*rdotu2)/r5)) - &
129     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
130     vterm*dtdr*d(1)/rij
131    
132     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
133     (5.0d0*(rdotu1*rdotu2)/r5)) - &
134     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
135     vterm*dtdr*d(2)/rij
136    
137     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
138     (5.0d0*(rdotu1*rdotu2)/r5)) - &
139     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
140     vterm*dtdr*d(3)/rij
141 mmeineke 469
142 mmeineke 377 dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
143     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
144     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
145 mmeineke 469
146 mmeineke 377 dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
147     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
148     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
149    
150 mmeineke 469
151 mmeineke 377 #ifdef IS_MPI
152 chuckv 482 f_Row(1,atom1) = f_Row(1,atom1) + dudx
153     f_Row(2,atom1) = f_Row(2,atom1) + dudy
154     f_Row(3,atom1) = f_Row(3,atom1) + dudz
155 mmeineke 377
156 chuckv 482 f_Col(1,atom2) = f_Col(1,atom2) - dudx
157     f_Col(2,atom2) = f_Col(2,atom2) - dudy
158     f_Col(3,atom2) = f_Col(3,atom2) - dudz
159 mmeineke 377
160     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
161     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
162     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
163    
164     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
165     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
166     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
167     #else
168 chuckv 482 f(1,atom1) = f(1,atom1) + dudx
169     f(2,atom1) = f(2,atom1) + dudy
170     f(3,atom1) = f(3,atom1) + dudz
171 mmeineke 377
172 chuckv 482 f(1,atom2) = f(1,atom2) - dudx
173     f(2,atom2) = f(2,atom2) - dudy
174     f(3,atom2) = f(3,atom2) - dudz
175 mmeineke 377
176     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
177     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
178     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
179    
180     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
181     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
182     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
183     #endif
184    
185     if (do_stress) then
186     tau_Temp(1) = tau_Temp(1) + dudx * d(1)
187     tau_Temp(2) = tau_Temp(2) + dudx * d(2)
188     tau_Temp(3) = tau_Temp(3) + dudx * d(3)
189     tau_Temp(4) = tau_Temp(4) + dudy * d(1)
190     tau_Temp(5) = tau_Temp(5) + dudy * d(2)
191     tau_Temp(6) = tau_Temp(6) + dudy * d(3)
192     tau_Temp(7) = tau_Temp(7) + dudz * d(1)
193     tau_Temp(8) = tau_Temp(8) + dudz * d(2)
194     tau_Temp(9) = tau_Temp(9) + dudz * d(3)
195     virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
196     endif
197    
198     endif
199    
200     return
201     end subroutine do_dipole_pair
202    
203     end module dipole_dipole