ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_reaction_field.F90
Revision: 394
Committed: Mon Mar 24 21:55:34 2003 UTC (21 years, 3 months ago) by gezelter
File size: 7541 byte(s)
Log Message:
electrostatic changes for dipole / RF separation

File Contents

# User Rev Content
1 mmeineke 377 module reaction_field
2     use force_globals
3     use definitions
4     use atype_module
5     use vector_class
6     #ifdef IS_MPI
7     use mpiSimulation
8     #endif
9     implicit none
10    
11     real(kind=dp), save :: rrf
12     real(kind=dp), save :: rt
13     real(kind=dp), save :: dielect
14     real(kind=dp), save :: rrfsq
15     real(kind=dp), save :: pre
16     logical, save :: rf_initialized = .false.
17    
18     contains
19    
20     subroutine initialize_rf(this_rrf, this_rt, this_dielect)
21     real(kind=dp), intent(in) :: this_rrf, this_rt, this_dielect
22    
23     rrf = this_rrf
24     rt = this_rt
25     dielect = this_dielect
26    
27     rrfsq = rrf * rrf
28     pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf)
29    
30 gezelter 394
31     write(*,*) 'rrf = ', rrf
32     write(*,*) 'rt = ', rt
33     write(*,*) 'dielect = ', dielect
34     write(*,*) 'pre = ', pre
35 mmeineke 377 rf_initialized = .true.
36    
37     return
38     end subroutine initialize_rf
39    
40     subroutine accumulate_rf(atom1, atom2, rij, u_l)
41    
42     integer, intent(in) :: atom1, atom2
43     real (kind = dp), intent(in) :: rij
44     real (kind = dp), dimension(:,:) :: u_l
45    
46     integer :: me1, me2
47     real (kind = dp) :: taper, mu1, mu2
48     real (kind = dp), dimension(3) :: ul1
49     real (kind = dp), dimension(3) :: ul2
50    
51     if (.not.rf_initialized) then
52     write(default_error,*) 'Reaction field not initialized!'
53     return
54     endif
55    
56     if (rij.le.rrf) then
57    
58     if (rij.lt.rt) then
59     taper = 1.0d0
60     else
61     taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
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     call getElementProperty(atypes, me1, "dipole_moment", mu1)
87     call getElementProperty(atypes, me2, "dipole_moment", mu2)
88    
89    
90     #ifdef IS_MPI
91     rf_Row(1,atom1) = rf_Row(1,atom1) + ul2(1)*mu2*taper
92     rf_Row(2,atom1) = rf_Row(2,atom1) + ul2(2)*mu2*taper
93     rf_Row(3,atom1) = rf_Row(3,atom1) + ul2(3)*mu2*taper
94    
95 gezelter 394 rf_Col(1,atom2) = rf_Col(1,atom2) + ul1(1)*mu1*taper
96     rf_Col(2,atom2) = rf_Col(2,atom2) + ul1(2)*mu1*taper
97     rf_Col(3,atom2) = rf_Col(3,atom2) + ul1(3)*mu1*taper
98 mmeineke 377 #else
99     rf(1,atom1) = rf(1,atom1) + ul2(1)*mu2*taper
100     rf(2,atom1) = rf(2,atom1) + ul2(2)*mu2*taper
101     rf(3,atom1) = rf(3,atom1) + ul2(3)*mu2*taper
102    
103 gezelter 394 rf(1,atom2) = rf(1,atom2) + ul1(1)*mu1*taper
104     rf(2,atom2) = rf(2,atom2) + ul1(2)*mu1*taper
105     rf(3,atom2) = rf(3,atom2) + ul1(3)*mu1*taper
106 mmeineke 377 #endif
107    
108     endif
109     return
110     end subroutine accumulate_rf
111    
112     subroutine accumulate_self_rf(atom1, mu1, u_l)
113    
114     integer, intent(in) :: atom1
115     real(kind=dp), intent(in) :: mu1
116     real(kind=dp), dimension(:,:) :: u_l
117    
118     !! should work for both MPI and non-MPI version since this is not pairwise.
119     rf(1,atom1) = rf(1,atom1) + u_l(1,atom1)*mu1
120     rf(2,atom1) = rf(2,atom1) + u_l(2,atom1)*mu1
121     rf(3,atom1) = rf(3,atom1) + u_l(3,atom1)*mu1
122    
123     return
124     end subroutine accumulate_self_rf
125    
126     subroutine reaction_field_final(a1, mu1, u_l, rfpot, t, do_pot)
127    
128     integer, intent(in) :: a1
129     real (kind=dp), intent(in) :: mu1
130     real (kind=dp), intent(inout) :: rfpot
131     logical, intent(in) :: do_pot
132     real (kind = dp), dimension(:,:) :: u_l
133     real (kind = dp), dimension(:,:) :: t
134    
135     if (.not.rf_initialized) then
136     write(default_error,*) 'Reaction field not initialized!'
137     return
138     endif
139    
140     ! compute torques on dipoles:
141     ! pre converts from mu in units of debye to kcal/mol
142    
143 gezelter 394 ! The torque contribution is dipole cross reaction_field
144    
145 mmeineke 377 t(1,a1) = t(1,a1) + pre*mu1*(u_l(2,a1)*rf(3,a1) - u_l(3,a1)*rf(2,a1))
146     t(2,a1) = t(2,a1) + pre*mu1*(u_l(3,a1)*rf(1,a1) - u_l(1,a1)*rf(3,a1))
147     t(3,a1) = t(3,a1) + pre*mu1*(u_l(1,a1)*rf(2,a1) - u_l(2,a1)*rf(1,a1))
148    
149     ! the potential contribution is -1/2 dipole dot reaction_field
150    
151     if (do_pot) then
152     rfpot = rfpot - 0.5d0 * pre * mu1 * &
153     (rf(1,a1)*u_l(1,a1) + rf(2,a1)*u_l(2,a1) + rf(3,a1)*u_l(3,a1))
154     endif
155    
156     return
157     end subroutine reaction_field_final
158    
159     subroutine rf_correct_forces(atom1, atom2, d, rij, u_l, f, do_stress)
160    
161     integer, intent(in) :: atom1, atom2
162     real(kind=dp), dimension(3), intent(in) :: d
163     real(kind=dp), intent(in) :: rij
164     real( kind = dp ), dimension(:,:) :: u_l
165     real( kind = dp ), dimension(:,:) :: f
166     logical, intent(in) :: do_stress
167    
168     real (kind = dp), dimension(3) :: ul1
169     real (kind = dp), dimension(3) :: ul2
170     real (kind = dp) :: dtdr
171     real (kind = dp) :: dudx, dudy, dudz, u1dotu2
172     integer :: me1, me2
173     real (kind = dp) :: mu1, mu2
174    
175     if (.not.rf_initialized) then
176     write(default_error,*) 'Reaction field not initialized!'
177     return
178     endif
179    
180     if (rij.le.rrf) then
181    
182     if (rij.lt.rt) then
183     dtdr = 0.0d0
184     else
185     dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
186     endif
187    
188     #ifdef IS_MPI
189     me1 = atid_Row(atom1)
190     ul1(1) = u_l_Row(1,atom1)
191     ul1(2) = u_l_Row(2,atom1)
192     ul1(3) = u_l_Row(3,atom1)
193    
194     me2 = atid_Col(atom2)
195     ul2(1) = u_l_Col(1,atom2)
196     ul2(2) = u_l_Col(2,atom2)
197     ul2(3) = u_l_Col(3,atom2)
198     #else
199     me1 = atid(atom1)
200     ul1(1) = u_l(1,atom1)
201     ul1(2) = u_l(2,atom1)
202     ul1(3) = u_l(3,atom1)
203    
204     me2 = atid(atom2)
205     ul2(1) = u_l(1,atom2)
206     ul2(2) = u_l(2,atom2)
207     ul2(3) = u_l(3,atom2)
208     #endif
209    
210     call getElementProperty(atypes, me1, "dipole_moment", mu1)
211     call getElementProperty(atypes, me2, "dipole_moment", mu2)
212    
213     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
214    
215     dudx = - pre*mu1*mu2*u1dotu2*dtdr*d(1)/rij
216     dudy = - pre*mu1*mu2*u1dotu2*dtdr*d(2)/rij
217     dudz = - pre*mu1*mu2*u1dotu2*dtdr*d(3)/rij
218    
219     #ifdef IS_MPI
220     f_Row(1,atom1) = f_Row(1,atom1) + dudx
221     f_Row(2,atom1) = f_Row(2,atom1) + dudy
222     f_Row(3,atom1) = f_Row(3,atom1) + dudz
223    
224     f_Col(1,atom2) = f_Col(1,atom2) - dudx
225     f_Col(2,atom2) = f_Col(2,atom2) - dudy
226     f_Col(3,atom2) = f_Col(3,atom2) - dudz
227     #else
228     f(1,atom1) = f(1,atom1) + dudx
229     f(2,atom1) = f(2,atom1) + dudy
230     f(3,atom1) = f(3,atom1) + dudz
231    
232     f(1,atom2) = f(1,atom2) - dudx
233     f(2,atom2) = f(2,atom2) - dudy
234     f(3,atom2) = f(3,atom2) - dudz
235     #endif
236    
237     if (do_stress) then
238     tau_Temp(1) = tau_Temp(1) + dudx * d(1)
239     tau_Temp(2) = tau_Temp(2) + dudx * d(2)
240     tau_Temp(3) = tau_Temp(3) + dudx * d(3)
241     tau_Temp(4) = tau_Temp(4) + dudy * d(1)
242     tau_Temp(5) = tau_Temp(5) + dudy * d(2)
243     tau_Temp(6) = tau_Temp(6) + dudy * d(3)
244     tau_Temp(7) = tau_Temp(7) + dudz * d(1)
245     tau_Temp(8) = tau_Temp(8) + dudz * d(2)
246     tau_Temp(9) = tau_Temp(9) + dudz * d(3)
247     virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
248     endif
249     endif
250    
251     return
252     end subroutine rf_correct_forces
253     end module reaction_field