ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 843
Committed: Wed Oct 29 20:41:39 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 6652 byte(s)
Log Message:
fixed a stdlib.h include error in bass.l

fixed a little bug in the first time step, regarding the setting of ecr and est in fortran

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