ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_reaction_field.F90
Revision: 377
Committed: Fri Mar 21 17:42:12 2003 UTC (21 years, 3 months ago) by mmeineke
Original Path: branches/mmeineke/OOPSE/libmdtools/calc_reaction_field.F90
File size: 7419 byte(s)
Log Message:
New OOPSE Tree

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