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