ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 945
Committed: Thu Jan 15 01:14:55 2004 UTC (20 years, 5 months ago) by gezelter
File size: 7794 byte(s)
Log Message:
More work for adding charges

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 chuckv 901 use status
9 mmeineke 377 #ifdef IS_MPI
10     use mpiSimulation
11     #endif
12     implicit none
13    
14 mmeineke 626 PRIVATE
15 gezelter 945 real(kind=dp), save :: ecr = 0.0
16 mmeineke 469 real(kind=dp), save :: rt = 0.0
17 mmeineke 626 real(kind=dp), save :: pre = 0.0
18 chuckv 901 logical, save :: haveCutoffs = .false.
19     logical, save :: haveMomentMap = .false.
20 mmeineke 377
21 mmeineke 626 public::setCutoffsDipole
22     public::do_dipole_pair
23    
24 chuckv 901 type :: MomentList
25     real(kind=DP) :: dipole_moment = 0.0_DP
26     end type MomentList
27    
28     type(MomentList), dimension(:),allocatable :: MomentMap
29    
30 mmeineke 377 contains
31    
32 gezelter 945 subroutine setCutoffsDipole(this_ecr, this_rt)
33     real(kind=dp), intent(in) :: this_ecr, this_rt
34     ecr = this_ecr
35 mmeineke 377 rt = this_rt
36 mmeineke 843
37     ! pre converts from mu in units of debye to kcal/mol
38 chuckv 388 pre = 14.38362_dp
39 gezelter 394
40 chuckv 901 haveCutoffs = .true.
41 mmeineke 377
42     return
43 mmeineke 626 end subroutine setCutoffsDipole
44 mmeineke 377
45 chuckv 901 subroutine createMomentMap(status)
46     integer :: nAtypes
47     integer :: status
48     integer :: i
49     real (kind=DP) :: thisDP
50     logical :: thisProperty
51    
52     status = 0
53    
54     nAtypes = getSize(atypes)
55    
56     if (nAtypes == 0) then
57     status = -1
58     return
59     end if
60    
61     if (.not. allocated(MomentMap)) then
62     allocate(MomentMap(nAtypes))
63     endif
64    
65     do i = 1, nAtypes
66    
67     call getElementProperty(atypes, i, "is_DP", thisProperty)
68    
69     if (thisProperty) then
70     call getElementProperty(atypes, i, "dipole_moment", thisDP)
71     MomentMap(i)%dipole_moment = thisDP
72     endif
73    
74     end do
75    
76     haveMomentMap = .true.
77    
78     end subroutine createMomentMap
79    
80    
81 chuckv 460 subroutine do_dipole_pair(atom1, atom2, d, rij, r2, pot, u_l, f, t, &
82 mmeineke 377 do_pot, do_stress)
83    
84     logical :: do_pot, do_stress
85    
86 tim 727 integer atom1, atom2, me1, me2, id1, id2
87 chuckv 901 integer :: localError
88 mmeineke 377 real(kind=dp) :: rij, mu1, mu2
89 mmeineke 469 real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5
90 mmeineke 377 real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
91     real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
92     real(kind=dp) :: taper, dtdr, vterm
93    
94     real( kind = dp ) :: pot
95     real( kind = dp ), dimension(3) :: d
96 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
97     real( kind = dp ), dimension(3,nLocal) :: f
98     real( kind = dp ), dimension(3,nLocal) :: t
99 mmeineke 377
100     real (kind = dp), dimension(3) :: ul1
101     real (kind = dp), dimension(3) :: ul2
102    
103 chuckv 901 if (.not. haveCutoffs) then
104     write(default_error,*) 'Dipole-dipole does not have cutoffs set!'
105 mmeineke 377 return
106     endif
107 chuckv 901
108     if (.not.haveMomentMap) then
109     localError = 0
110     call createMomentMap(localError)
111     if ( localError .ne. 0 ) then
112     call handleError("dipole-dipole", "MomentMap creation failed!")
113     return
114     end if
115     endif
116    
117 mmeineke 377 #ifdef IS_MPI
118     me1 = atid_Row(atom1)
119     ul1(1) = u_l_Row(1,atom1)
120     ul1(2) = u_l_Row(2,atom1)
121     ul1(3) = u_l_Row(3,atom1)
122    
123     me2 = atid_Col(atom2)
124     ul2(1) = u_l_Col(1,atom2)
125     ul2(2) = u_l_Col(2,atom2)
126     ul2(3) = u_l_Col(3,atom2)
127     #else
128     me1 = atid(atom1)
129     ul1(1) = u_l(1,atom1)
130     ul1(2) = u_l(2,atom1)
131     ul1(3) = u_l(3,atom1)
132    
133     me2 = atid(atom2)
134     ul2(1) = u_l(1,atom2)
135     ul2(2) = u_l(2,atom2)
136     ul2(3) = u_l(3,atom2)
137     #endif
138    
139 chuckv 901 mu1 = MomentMap(me1)%dipole_moment
140     mu2 = MomentMap(me2)%dipole_moment
141 mmeineke 443
142 gezelter 945 if (rij.le.ecr) then
143 mmeineke 377
144     if (rij.lt.rt) then
145     taper = 1.0d0
146     dtdr = 0.0d0
147     else
148 gezelter 945 taper = (ecr + 2.0d0*rij - 3.0d0*rt)*(ecr-rij)**2/ ((ecr-rt)**3)
149     dtdr = 6.0d0*(rij*rij - rij*rt - rij*ecr +ecr*rt)/((ecr-rt)**3)
150 mmeineke 377 endif
151 chuckv 901
152 mmeineke 377 r3 = r2*rij
153     r5 = r3*r2
154    
155     rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
156     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
157     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
158 chuckv 901
159 mmeineke 377 dip2 = pre * mu1 * mu2
160     dfact1 = 3.0d0*dip2 / r2
161     dfact2 = 3.0d0*dip2 / r5
162 chuckv 901
163 mmeineke 377 vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
164 chuckv 901
165 mmeineke 377 if (do_pot) then
166     #ifdef IS_MPI
167     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
168     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
169     #else
170     pot = pot + vterm*taper
171     #endif
172     endif
173    
174     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
175     (5.0d0*(rdotu1*rdotu2)/r5)) - &
176     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
177     vterm*dtdr*d(1)/rij
178    
179     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
180     (5.0d0*(rdotu1*rdotu2)/r5)) - &
181     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
182     vterm*dtdr*d(2)/rij
183    
184     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
185     (5.0d0*(rdotu1*rdotu2)/r5)) - &
186     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
187     vterm*dtdr*d(3)/rij
188 mmeineke 469
189 mmeineke 377 dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
190     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
191     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
192 mmeineke 469
193 mmeineke 377 dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
194     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
195     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
196    
197 mmeineke 469
198 mmeineke 377 #ifdef IS_MPI
199 chuckv 482 f_Row(1,atom1) = f_Row(1,atom1) + dudx
200     f_Row(2,atom1) = f_Row(2,atom1) + dudy
201     f_Row(3,atom1) = f_Row(3,atom1) + dudz
202 mmeineke 377
203 chuckv 482 f_Col(1,atom2) = f_Col(1,atom2) - dudx
204     f_Col(2,atom2) = f_Col(2,atom2) - dudy
205     f_Col(3,atom2) = f_Col(3,atom2) - dudz
206 mmeineke 377
207     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
208     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
209     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
210    
211     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
212     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
213     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
214     #else
215 chuckv 482 f(1,atom1) = f(1,atom1) + dudx
216     f(2,atom1) = f(2,atom1) + dudy
217     f(3,atom1) = f(3,atom1) + dudz
218 mmeineke 377
219 chuckv 482 f(1,atom2) = f(1,atom2) - dudx
220     f(2,atom2) = f(2,atom2) - dudy
221     f(3,atom2) = f(3,atom2) - dudz
222 mmeineke 377
223     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
224     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
225     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
226    
227     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
228     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
229     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
230     #endif
231 tim 727
232 gezelter 730 if (do_stress) then
233    
234 tim 727 #ifdef IS_MPI
235     id1 = tagRow(atom1)
236     id2 = tagColumn(atom2)
237     #else
238     id1 = atom1
239     id2 = atom2
240 gezelter 730 #endif
241 tim 727 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
242 gezelter 483
243 gezelter 611 ! because the d vector is the rj - ri vector, and
244     ! because dudx, dudy, dudz are the (positive) force on
245     ! atom i (negative on atom j) we need a negative sign here:
246    
247     tau_Temp(1) = tau_Temp(1) - d(1) * dudx
248     tau_Temp(2) = tau_Temp(2) - d(1) * dudy
249     tau_Temp(3) = tau_Temp(3) - d(1) * dudz
250     tau_Temp(4) = tau_Temp(4) - d(2) * dudx
251     tau_Temp(5) = tau_Temp(5) - d(2) * dudy
252     tau_Temp(6) = tau_Temp(6) - d(2) * dudz
253     tau_Temp(7) = tau_Temp(7) - d(3) * dudx
254     tau_Temp(8) = tau_Temp(8) - d(3) * dudy
255     tau_Temp(9) = tau_Temp(9) - d(3) * dudz
256    
257 gezelter 483 virial_Temp = virial_Temp + &
258     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
259    
260     endif
261 mmeineke 377 endif
262    
263     endif
264    
265     return
266     end subroutine do_dipole_pair
267    
268     end module dipole_dipole