ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_reaction_field.F90
Revision: 621
Committed: Wed Jul 16 02:11:02 2003 UTC (20 years, 11 months ago) by gezelter
File size: 8077 byte(s)
Log Message:
more fixes for box changes

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