ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 378
Committed: Fri Mar 21 17:42:12 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 6133 byte(s)
Log Message:
This commit was generated by cvs2svn to compensate for changes in r377,
which included commits to RCS files with non-trunk default branches.

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     pre = 14.38362d0
27    
28     dipole_initialized = .true.
29    
30     return
31     end subroutine initialize_dipole
32    
33    
34     subroutine do_dipole_pair(atom1, atom2, d, rij, r2, pot, u_l, f, t, &
35     do_pot, do_stress)
36    
37     logical :: do_pot, do_stress
38    
39     integer atom1, atom2, me1, me2
40     real(kind=dp) :: rij, mu1, mu2
41     real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5, pre
42     real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
43     real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
44     real(kind=dp) :: taper, dtdr, vterm
45    
46     real( kind = dp ) :: pot
47     real( kind = dp ), dimension(3) :: d
48     real( kind = dp ), dimension(:,:) :: u_l
49     real( kind = dp ), dimension(:,:) :: f
50     real( kind = dp ), dimension(:,:) :: t
51    
52     real (kind = dp), dimension(3) :: ul1
53     real (kind = dp), dimension(3) :: ul2
54    
55     if (.not.dipole_initialized) then
56     write(default_error,*) 'Dipole-dipole not initialized!'
57     return
58     endif
59    
60     #ifdef IS_MPI
61     me1 = atid_Row(atom1)
62     ul1(1) = u_l_Row(1,atom1)
63     ul1(2) = u_l_Row(2,atom1)
64     ul1(3) = u_l_Row(3,atom1)
65    
66     me2 = atid_Col(atom2)
67     ul2(1) = u_l_Col(1,atom2)
68     ul2(2) = u_l_Col(2,atom2)
69     ul2(3) = u_l_Col(3,atom2)
70     #else
71     me1 = atid(atom1)
72     ul1(1) = u_l(1,atom1)
73     ul1(2) = u_l(2,atom1)
74     ul1(3) = u_l(3,atom1)
75    
76     me2 = atid(atom2)
77     ul2(1) = u_l(1,atom2)
78     ul2(2) = u_l(2,atom2)
79     ul2(3) = u_l(3,atom2)
80     #endif
81    
82     call getElementProperty(atypes, me1, "dipole_moment", mu1)
83     call getElementProperty(atypes, me2, "dipole_moment", mu2)
84    
85     if (rij.le.rrf) then
86    
87     if (rij.lt.rt) then
88     taper = 1.0d0
89     dtdr = 0.0d0
90     else
91     taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
92     dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
93     endif
94    
95     r3 = r2*rij
96     r5 = r3*r2
97    
98     rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
99     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
100     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
101    
102     dip2 = pre * mu1 * mu2
103    
104     dfact1 = 3.0d0*dip2 / r2
105     dfact2 = 3.0d0*dip2 / r5
106    
107     vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
108    
109     if (do_pot) then
110     #ifdef IS_MPI
111     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
112     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
113     #else
114     pot = pot + vterm*taper
115     #endif
116     endif
117    
118     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
119     (5.0d0*(rdotu1*rdotu2)/r5)) - &
120     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
121     vterm*dtdr*d(1)/rij
122    
123     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
124     (5.0d0*(rdotu1*rdotu2)/r5)) - &
125     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
126     vterm*dtdr*d(2)/rij
127    
128     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
129     (5.0d0*(rdotu1*rdotu2)/r5)) - &
130     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
131     vterm*dtdr*d(3)/rij
132    
133     dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
134     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
135     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
136    
137     dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
138     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
139     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
140    
141    
142     #ifdef IS_MPI
143     f_Row(1,atom1) = f_Row(1,atom1) + dudx
144     f_Row(2,atom1) = f_Row(2,atom1) + dudy
145     f_Row(3,atom1) = f_Row(3,atom1) + dudz
146    
147     f_Col(1,atom2) = f_Col(1,atom2) - dudx
148     f_Col(2,atom2) = f_Col(2,atom2) - dudy
149     f_Col(3,atom2) = f_Col(3,atom2) - dudz
150    
151     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
152     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
153     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
154    
155     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
156     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
157     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
158     #else
159     f(1,atom1) = f(1,atom1) + dudx
160     f(2,atom1) = f(2,atom1) + dudy
161     f(3,atom1) = f(3,atom1) + dudz
162    
163     f(1,atom2) = f(1,atom2) - dudx
164     f(2,atom2) = f(2,atom2) - dudy
165     f(3,atom2) = f(3,atom2) - dudz
166    
167     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
168     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
169     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
170    
171     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
172     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
173     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
174     #endif
175    
176     if (do_stress) then
177     tau_Temp(1) = tau_Temp(1) + dudx * d(1)
178     tau_Temp(2) = tau_Temp(2) + dudx * d(2)
179     tau_Temp(3) = tau_Temp(3) + dudx * d(3)
180     tau_Temp(4) = tau_Temp(4) + dudy * d(1)
181     tau_Temp(5) = tau_Temp(5) + dudy * d(2)
182     tau_Temp(6) = tau_Temp(6) + dudy * d(3)
183     tau_Temp(7) = tau_Temp(7) + dudz * d(1)
184     tau_Temp(8) = tau_Temp(8) + dudz * d(2)
185     tau_Temp(9) = tau_Temp(9) + dudz * d(3)
186     virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
187     endif
188    
189     endif
190    
191     return
192     end subroutine do_dipole_pair
193    
194     end module dipole_dipole