ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/mpole.F90
Revision: 1464
Committed: Fri Jul 9 19:29:05 2010 UTC (15 years ago) by gezelter
File size: 15258 byte(s)
Log Message:
removing cruft (atom numbers, do_pot, do_stress) from many modules and
force managers

File Contents

# Content
1 subroutine doMultipolePair(atom1, atom2, d, rij, r2, rcut, sw, &
2 vpair, fpair, pot, eFrame, f, t)
3
4 !***************************************************************
5 ! doMultipolePair evaluates the potential, forces, and torques
6 ! between point multipoles. It is based on the ewald2 real
7 ! space routine in the MDMULP code written by W. Smith and
8 ! available from the CCP5 website.
9 !
10 ! We're using the damped real space portion of the Ewald sum
11 ! as the entire interaction. Details on this portion of the
12 ! sum can be found in:
13 !
14 ! "Point Multipoles in the Ewald Summation (Revisited)," by
15 ! W. Smith, CCP5 Newsletter, 46, pp. 18-30 (1998).
16 !
17 !**************************************************************
18
19 integer, intent(in) :: atom1, atom2
20 real(kind=dp), intent(in) :: rij, r2, sw, rcut
21 real(kind=dp), intent(in), dimension(3) :: d
22 real(kind=dp), intent(inout) :: vpair
23 real(kind=dp), intent(inout), dimension(3) :: fpair
24
25 real( kind = dp ) :: pot
26 real( kind = dp ), dimension(9,nLocal) :: eFrame
27 real( kind = dp ), dimension(3,nLocal) :: f
28 real( kind = dp ), dimension(3,nLocal) :: t
29 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
30 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
31 integer :: me1, me2, id1, id2
32 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
33 real (kind=dp) :: qxx_i, qyy_i, qzz_i
34 real (kind=dp) :: qxx_j, qyy_j, qzz_j
35
36 real (kind=dp) :: epot, qii, alsq2, alsq2n, exp2a, ralpha
37 real (kind=dp), dimension(3) :: ftm, ttm_i, ttm_j
38 real (kind=dp), dimension(3) :: qir, qjr, qidj, qjdi, qiqjr, qjqir
39 real (kind=dp), dimension(3) :: dixdj, qixqj, dixr, djxr, rxqir, rxqjr
40 real (kind=dp), dimension(3) :: dixqjr, djxqir, rxqidj, rxqjdi, rxqijr
41 real (kind=dp), dimension(3) :: rxqjir, qjrxqir
42 real (kind=dp), dimension(6) :: bn
43 real (kind=dp), dimension(10) :: sc
44 real (kind=dp), dimension(9) :: gl
45
46 real(kind=dp) :: a1, a2, a3, a4, a5, p
47
48 DATA A1,A2,A3/0.254829592,-0.284496736,1.421413741/
49 DATA A4,A5,P/-1.453152027,1.061405429,0.3275911/
50
51 if (.not.summationMethodChecked) then
52 call checkSummationMethod()
53 endif
54
55 #ifdef IS_MPI
56 me1 = atid_Row(atom1)
57 me2 = atid_Col(atom2)
58 #else
59 me1 = atid(atom1)
60 me2 = atid(atom2)
61 #endif
62
63 ! logicals
64 i_is_Charge = ElectrostaticMap(me1)%is_Charge
65 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
66 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
67
68 j_is_Charge = ElectrostaticMap(me2)%is_Charge
69 j_is_Dipole = ElectrostaticMap(me2)%is_Dipole
70 j_is_Quadrupole = ElectrostaticMap(me2)%is_Quadrupole
71
72 if (i_is_Charge) then
73 q_i = ElectrostaticMap(me1)%charge
74 else
75 q_i = 0.0_dp
76 endif
77
78 if (j_is_Charge) then
79 q_j = ElectrostaticMap(me2)%charge
80 else
81 q_j = 0.0_dp
82 endif
83
84 if (i_is_Dipole) then
85 mu_i = ElectrostaticMap(me1)%dipole_moment
86 #ifdef IS_MPI
87 d_i(1) = eFrame_Row(3,atom1)
88 d_i(2) = eFrame_Row(6,atom1)
89 d_i(3) = eFrame_Row(9,atom1)
90 #else
91 d_i(1) = eFrame(3,atom1)
92 d_i(2) = eFrame(6,atom1)
93 d_i(3) = eFrame(9,atom1)
94 #endif
95 d_i = d_i * mu_i
96 else
97 d_i = 0.0_dp
98 endif
99
100 if (j_is_Dipole) then
101 mu_j = ElectrostaticMap(me2)%dipole_moment
102 #ifdef IS_MPI
103 d_j(1) = eFrame_Col(3,atom2)
104 d_j(2) = eFrame_Col(6,atom2)
105 d_j(3) = eFrame_Col(9,atom2)
106 #else
107 d_j(1) = eFrame(3,atom2)
108 d_j(2) = eFrame(6,atom2)
109 d_j(3) = eFrame(9,atom2)
110 #endif
111 d_j = d_j * mu_j
112 else
113 d_j = 0.0_dp
114 endif
115
116 if (i_is_Quadrupole) then
117 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
118 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
119 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
120 #ifdef IS_MPI
121 qpole_i(1) = qxx_i * eFrame_Row(1,atom1)
122 qpole_i(2) = qxx_i * eFrame_Row(4,atom1)
123 qpole_i(3) = qxx_i * eFrame_Row(7,atom1)
124 qpole_i(4) = qyy_i * eFrame_Row(2,atom1)
125 qpole_i(5) = qyy_i * eFrame_Row(5,atom1)
126 qpole_i(6) = qyy_i * eFrame_Row(8,atom1)
127 qpole_i(7) = qzz_i * eFrame_Row(3,atom1)
128 qpole_i(8) = qzz_i * eFrame_Row(6,atom1)
129 qpole_i(9) = qzz_i * eFrame_Row(9,atom1)
130 #else
131 qpole_i(1) = qxx_i * eFrame(1,atom1)
132 qpole_i(2) = qxx_i * eFrame(4,atom1)
133 qpole_i(3) = qxx_i * eFrame(7,atom1)
134 qpole_i(4) = qyy_i * eFrame(2,atom1)
135 qpole_i(5) = qyy_i * eFrame(5,atom1)
136 qpole_i(6) = qyy_i * eFrame(8,atom1)
137 qpole_i(7) = qzz_i * eFrame(3,atom1)
138 qpole_i(8) = qzz_i * eFrame(6,atom1)
139 qpole_i(9) = qzz_i * eFrame(9,atom1)
140 #endif
141 else
142 qpole_i = 0.0_dp
143 endif
144
145 if (j_is_Quadrupole) then
146 qxx_j = ElectrostaticMap(me2)%quadrupole_moments(1)
147 qyy_j = ElectrostaticMap(me2)%quadrupole_moments(2)
148 qzz_j = ElectrostaticMap(me2)%quadrupole_moments(3)
149 #ifdef IS_MPI
150 qpole_j(1) = qxx_j * eFrame_Col(1,atom2)
151 qpole_j(2) = qxx_j * eFrame_Col(4,atom2)
152 qpole_j(3) = qxx_j * eFrame_Col(7,atom2)
153 qpole_j(4) = qyy_j * eFrame_Col(2,atom2)
154 qpole_j(5) = qyy_j * eFrame_Col(5,atom2)
155 qpole_j(6) = qyy_j * eFrame_Col(8,atom2)
156 qpole_j(7) = qzz_j * eFrame_Col(3,atom2)
157 qpole_j(8) = qzz_j * eFrame_Col(6,atom2)
158 qpole_j(9) = qzz_j * eFrame_Col(9,atom2)
159 #else
160 qpole_j(1) = qxx_j * eFrame(1,atom2)
161 qpole_j(2) = qxx_j * eFrame(4,atom2)
162 qpole_j(3) = qxx_j * eFrame(7,atom2)
163 qpole_j(4) = qyy_j * eFrame(2,atom2)
164 qpole_j(5) = qyy_j * eFrame(5,atom2)
165 qpole_j(6) = qyy_j * eFrame(8,atom2)
166 qpole_j(7) = qzz_j * eFrame(3,atom2)
167 qpole_j(8) = qzz_j * eFrame(6,atom2)
168 qpole_j(9) = qzz_j * eFrame(9,atom2)
169 #endif
170 else
171 qpole_j = 0.0_dp
172 endif
173
174 !
175 ! INITIALISE TEMPORARY FORCE AND TORQUE ACCUMULATORS
176 !
177
178 ftm = 0.0
179 ttm_i = 0.0
180 ttm_j = 0.0
181
182 ri = 1.0_dp / rij
183 ri2 = ri * ri
184
185 !
186 ! CALCULATE BN FUNCTIONS
187 !
188 if (screeningMethod .eq. DAMPED) then
189 ! assemble the damping variables
190 call lookupUniformSpline1d(erfcSpline, rij, erfcVal, derfcVal)
191 else
192 dampingAlpha = 0.0_dp
193 erfcVal = 1.0_dp
194 endif
195
196 rrtpi = 1.0 / sqrt(pi)
197 alsq2 = 2.0 * dampingAlpha**2
198 ralpi = 0.0
199 if(dampingAlpha .gt. 0.0) ralpi = rrtpi / dampingAlpha
200
201 alphar = dampingAlpha * rij
202 t = 1.0 / (1.0 + p*alphar)
203 exp2a = exp(-alphar**2)
204 bn(1) = ((((a5*t+a4)*t+a3)*t+a2)*t+a1) * t * exp2a * ri
205 alsq2n = ralpi
206
207 do n = 1, 5
208 bfac = real(n+n-1, kind=dp)
209 alsq2n = alsq2*alsq2n
210 bn(n+1) = ri2 * (bfac*bn(n) + alsq2n*exp2a)
211 enddo
212
213 ralpha = dampingAlpha * rij
214 bn(0) = erfcVal * ri
215 alsq2 = 2.0d0 * dampingAlpha**2
216 alsq2n = 0.0d0
217 if (dampingAlpha .gt. 0.0_dp) alsq2n = 1.0_dp / (sqrtpi*dampingAlpha)
218 exp2a = exp(-ralpha**2)
219 do m = 1, 5
220 bfac = real(m+m-1, kind=dp)
221 alsq2n = alsq2 * alsq2n
222 bn(m) = (bfac*bn(m-1)+alsq2n*exp2a) / r2
223 end do
224
225 !
226 ! CONSTRUCT AUXILIARY VECTORS
227 !
228
229 dixdj(1) = d_i(2)*d_j(3) - d_i(3)*d_j(2)
230 dixdj(2) = d_i(3)*d_j(1) - d_i(1)*d_j(3)
231 dixdj(3) = d_i(1)*d_j(2) - d_i(2)*d_j(1)
232
233 dixr(1) = d_i(2)*d(3) - d_i(3)*d(2)
234 dixr(2) = d_i(3)*d(1) - d_i(1)*d(3)
235 dixr(3) = d_i(1)*d(2) - d_i(2)*d(1)
236
237 djxr(1) = d_j(2)*d(3) - d_j(3)*d(2)
238 djxr(2) = d_j(3)*d(1) - d_j(1)*d(3)
239 djxr(3) = d_j(1)*d(2) - d_j(2)*d(1)
240
241 qir(1) = qpole_i(1)*d(1) + qpole_i(4)*d(2) + qpole_i(7)*d(3)
242 qir(2) = qpole_i(2)*d(1) + qpole_i(5)*d(2) + qpole_i(8)*d(3)
243 qir(3) = qpole_i(3)*d(1) + qpole_i(6)*d(2) + qpole_i(9)*d(3)
244
245 qjr(1) = qpole_j(1)*d(1) + qpole_j(4)*d(2) + qpole_j(7)*d(3)
246 qjr(2) = qpole_j(2)*d(1) + qpole_j(5)*d(2) + qpole_j(8)*d(3)
247 qjr(3) = qpole_j(3)*d(1) + qpole_j(6)*d(2) + qpole_j(9)*d(3)
248
249 qiqjr(1) = qpole_i(1)*qjr(1) + qpole_i(4)*qjr(2) + qpole_i(7)*qjr(3)
250 qiqjr(2) = qpole_i(2)*qjr(1) + qpole_i(5)*qjr(2) + qpole_i(8)*qjr(3)
251 qiqjr(3) = qpole_i(3)*qjr(1) + qpole_i(6)*qjr(2) + qpole_i(9)*qjr(3)
252
253
254 qjqir(1) = qpole_j(1)*QIR(1) + qpole_j(4)*QIR(2) + qpole_j(7)*QIR(3)
255 qjqir(2) = qpole_j(2)*QIR(1) + qpole_j(5)*QIR(2) + qpole_j(8)*QIR(3)
256 qjqir(3) = qpole_j(3)*QIR(1) + qpole_j(6)*QIR(2) + qpole_j(9)*QIR(3)
257
258 qixqj(1) = qpole_i(2)*qpole_j(3) + qpole_i(5)*qpole_j(6) + &
259 qpole_i(8)*qpole_j(9) - qpole_i(3)*qpole_j(2) - &
260 qpole_i(6)*qpole_j(5) - qpole_i(9)*qpole_j(8)
261
262 qixqj(2) = qpole_i(3)*qpole_j(1) + qpole_i(6)*qpole_j(4) + &
263 qpole_i(9)*qpole_j(7) - qpole_i(1)*qpole_j(3) - &
264 qpole_i(4)*qpole_j(6) - qpole_i(7)*qpole_j(9)
265
266 qixqj(3) = qpole_i(1)*qpole_j(2) + qpole_i(4)*qpole_j(5) + &
267 qpole_i(7)*qpole_j(8) - qpole_i(2)*qpole_j(1) - &
268 qpole_i(5)*qpole_j(4) - qpole_i(8)*qpole_j(7)
269
270 rxqir(1) = d(2)*qir(3) - d(3)*qir(2)
271 rxqir(2) = d(3)*qir(1) - d(1)*qir(3)
272 rxqir(3) = d(1)*qir(2) - d(2)*qir(1)
273
274 rxqjr(1) = d(2)*qjr(3) - d(3)*qjr(2)
275 rxqjr(2) = d(3)*qjr(1) - d(1)*qjr(3)
276 rxqjr(3) = d(1)*qjr(2) - d(2)*qjr(1)
277
278 rxqijr(1) = d(2)*qiqjr(3) - d(3)*qiqjr(2)
279 rxqijr(2) = d(3)*qiqjr(1) - d(1)*qiqjr(3)
280 rxqijr(3) = d(1)*qiqjr(2) - d(2)*qiqjr(1)
281
282 rxqjir(1) = d(2)*qjqir(3) - d(3)*qjqir(2)
283 rxqjir(2) = d(3)*qjqir(1) - d(1)*qjqir(3)
284 rxqjir(3) = d(1)*qjqir(2) - d(2)*qjqir(1)
285
286 qjrxqir(1) = qjr(2)*qir(3) - qjr(3)*qir(2)
287 qjrxqir(2) = qjr(3)*qir(1) - qjr(1)*qir(3)
288 qjrxqir(3) = qjr(1)*qir(2) - qjr(2)*qir(1)
289
290 qidj(1) = qpole_i(1)*d_j(1) + qpole_i(4)*d_j(2) + qpole_i(7)*d_j(3)
291 qidj(2) = qpole_i(2)*d_j(1) + qpole_i(5)*d_j(2) + qpole_i(8)*d_j(3)
292 qidj(3) = qpole_i(3)*d_j(1) + qpole_i(6)*d_j(2) + qpole_i(9)*d_j(3)
293
294 qjdi(1) = qpole_j(1)*d_i(1) + qpole_j(4)*d_i(2) + qpole_j(7)*d_i(3)
295 qjdi(2) = qpole_j(2)*d_i(1) + qpole_j(5)*d_i(2) + qpole_j(8)*d_i(3)
296 qjdi(3) = qpole_j(3)*d_i(1) + qpole_j(6)*d_i(2) + qpole_j(9)*d_i(3)
297
298 dixqjr(1) = d_i(2)*qjr(3) - d_i(3)*qjr(2)
299 dixqjr(2) = d_i(3)*qjr(1) - d_i(1)*qjr(3)
300 dixqjr(3) = d_i(1)*qjr(2) - d_i(2)*qjr(1)
301
302 djxqir(1) = d_j(2)*qir(3) - d_j(3)*qir(2)
303 djxqir(2) = d_j(3)*qir(1) - d_j(1)*qir(3)
304 djxqir(3) = d_j(1)*qir(2) - d_j(2)*qir(1)
305
306 rxqidj(1) = d(2)*qidj(3) - d(3)*qidj(2)
307 rxqidj(2) = d(3)*qidj(1) - d(1)*qidj(3)
308 rxqidj(3) = d(1)*qidj(2) - d(2)*qidj(1)
309
310 rxqjdi(1) = d(2)*qjdi(3) - d(3)*qjdi(2)
311 rxqjdi(2) = d(3)*qjdi(1) - d(1)*qjdi(3)
312 rxqjdi(3) = d(1)*qjdi(2) - d(2)*qjdi(1)
313
314 !
315 ! CALCULATE SCALAR PRODUCTS
316 !
317
318 qii = qpole_i(1) + qpole_i(5) + qpole_i(9)
319
320 sc(1) = qpole_j(1) + qpole_j(5) + qpole_j(9)
321 sc(2) = d_i(1)*d_j(1) + d_i(2)*d_j(2) + d_i(3)*d_j(3)
322 sc(3) = d_i(1)*d(1) + d_i(2)*d(2) + d_i(3)*d(3)
323 sc(4) = d_j(1)*d(1) + d_j(2)*d(2) + d_j(3)*d(3)
324 sc(5) = qir(1)*d(1) + qir(2)*d(2) + qir(3)*d(3)
325 sc(6) = qjr(1)*d(1) + qjr(2)*d(2) + qjr(3)*d(3)
326 sc(7) = qir(1)*d_j(1) + qir(2)*d_j(2) + qir(3)*d_j(3)
327 sc(8) = qjr(1)*d_i(1) + qjr(2)*d_i(2) + qjr(3)*d_i(3)
328 sc(9) = qir(1)*qjr(1) + qir(2)*qjr(2) + qir(3)*qjr(3)
329 sc(10) = qpole_i(1)*qpole_j(1) + qpole_i(2)*qpole_j(2) + &
330 qpole_i(3)*qpole_j(3) + qpole_i(4)*qpole_j(4) + &
331 qpole_i(5)*qpole_j(5) + qpole_i(6)*qpole_j(6) + &
332 qpole_i(7)*qpole_j(7) + qpole_i(8)*qpole_j(8) + &
333 qpole_i(9)*qpole_j(9)
334
335 !
336 ! CALCULATE GL FUNCTIONS FOR POTENTIAL
337 !
338
339 gl(1) = q_i*q_j
340 gl(2) = q_j*sc(3) - q_i*sc(4)
341 gl(3) = q_i*sc(6) + q_j*sc(5) - sc(3)*sc(4)
342 gl(4) = sc(3)*sc(6) - sc(4)*sc(5)
343 gl(5) = sc(5)*sc(6)
344 gl(7) = sc(2) - q_j*qii - q_i*sc(1)
345 gl(8) = sc(4)*qii - sc(1)*sc(3) + 2.0*(sc(7) - sc(8))
346 gl(9) = 2.0*sc(10) + sc(1)*qii
347 gl(6) = -(sc(1)*sc(5) + sc(6)*qii + 4.0*sc(9))
348
349 !
350 ! CALCULATE POTENTIAL AND VIRIAL
351 !
352
353 epot = bn(1)*gl(1) + bn(2)*(gl(2) + gl(7)) + bn(3)*(gl(3) &
354 + gl(8) + gl(9)) + bn(4)*(gl(4) + gl(6)) + bn(5)*gl(5)
355
356 vpair = vpair + epot
357
358 #ifdef IS_MPI
359 pot_row(ELECTROSTATIC_POT,atom1) = pot_row(ELECTROSTATIC_POT,atom1) + &
360 0.5_dp*epot*sw
361 pot_col(ELECTROSTATIC_POT,atom2) = pot_col(ELECTROSTATIC_POT,atom2) + &
362 0.5_dp*epot*sw
363 #else
364 pot = pot + epot*sw
365 #endif
366
367 !
368 ! CALCULATE FORCE AND TORQUE COEFFICIENTS
369 !
370
371 gl(1) = bn(2)*gl(1) + bn(3)*(gl(2) + gl(7)) + bn(4)*(gl(3) &
372 + gl(8) + gl(9)) + bn(5)*(gl(4) + gl(6)) + bn(6)*gl(5)
373 gl(2) = -q_j*bn(2) + (sc(4) + sc(1))*bn(3) - sc(6)*bn(4)
374 gl(3) = q_i*bn(2) + (sc(3) - qii)*bn(3) + sc(5)*bn(4)
375 gl(4) = 2.0*bn(3)
376 gl(5) = 2.0*( - q_j*bn(3) + (sc(1) + sc(4))*bn(4) - sc(6)* bn(5))
377 gl(6) = 2.0*( - q_i*bn(3) + (qii - sc(3))*bn(4) - sc(5)*bn(5))
378 gl(7) = 4.0*bn(4)
379
380 !
381 ! CALCULATE FORCES BETWEEN MULTIPOLES I AND J
382 !
383
384 ftm(1) = gl(1)*d(1) + gl(2)*d_i(1) + gl(3)*d_j(1) + &
385 gl(4)*(qjdi(1) - qidj(1)) + gl(5)*qir(1) + &
386 gl(6)*qjr(1) + gl(7)*(qiqjr(1) + qjqir(1))
387 ftm(2) = gl(1)*d(2) + gl(2)*d_i(2) + gl(3)*d_j(2) + &
388 gl(4)*(qjdi(2) - qidj(2)) + gl(5)*qir(2) + &
389 gl(6)*qjr(2) + gl(7)*(qiqjr(2) + qjqir(2))
390 ftm(3) = gl(1)*d(3) + gl(2)*d_i(3) + gl(3)*d_j(3) + &
391 gl(4)*(qjdi(3) - qidj(3)) + gl(5)*qir(3) + &
392 gl(6)*qjr(3) + gl(7)*(qiqjr(3) + qjqir(3))
393
394 !
395 ! CALCULATE TORQUES BETWEEN MULTIPOLES I AND J
396 !
397
398 ttm_i(1) = - bn(2)*dixdj(1) + gl(2)*dixr(1) + gl(4)* &
399 (dixqjr(1) + djxqir(1) + rxqidj(1) - 2.0*qixqj(1)) - &
400 gl(5)*rxqir(1) - gl(7)*(rxqijr(1) + qjrxqir(1))
401 ttm_i(2) = - bn(2)*dixdj(2) + gl(2)*dixr(2) + gl(4)* &
402 (dixqjr(2) + djxqir(2) + rxqidj(2) - 2.0*qixqj(2)) - &
403 gl(5)*rxqir(2) - gl(7)*(rxqijr(2) + qjrxqir(2))
404 ttm_i(3) = - bn(2)*dixdj(3) + gl(2)*dixr(3) + gl(4)* &
405 (dixqjr(3) + djxqir(3) + rxqidj(3) - 2.0*qixqj(3)) - &
406 gl(5)*rxqir(3) - gl(7)*(rxqijr(3) + qjrxqir(3))
407 ttm_j(1) = bn(2)*dixdj(1) + gl(3)*djxr(1) - gl(4)* &
408 (dixqjr(1) + djxqir(1) + rxqjdi(1) - 2.0*qixqj(1)) - &
409 gl(6)*rxqjr(1) - gl(7)*(rxqjir(1) - qjrxqir(1))
410 ttm_j(2) = bn(2)*dixdj(2) + gl(3)*djxr(2) - gl(4)* &
411 (dixqjr(2) + djxqir(2) + rxqjdi(2) - 2.0*qixqj(2)) - &
412 gl(6)*rxqjr(2) - gl(7)*(rxqjir(2) - qjrxqir(2))
413 ttm_j(3) = bn(2)*dixdj(3) + gl(3)*djxr(3) - gl(4)* &
414 (dixqjr(3) + djxqir(3) + rxqjdi(3) - 2.0*qixqj(3)) - &
415 gl(6)*rxqjr(3) - gl(7)*(rxqjir(3) - qjrxqir(3))
416
417 ftm = ftm*sw
418 ttm_i = ttm_i*sw
419 ttm_j = ttm_j*sw
420
421 #ifdef IS_MPI
422 f_Row(1,atom1) = f_Row(1,atom1) + ftm(1)
423 f_Row(2,atom1) = f_Row(2,atom1) + ftm(2)
424 f_Row(3,atom1) = f_Row(3,atom1) + ftm(3)
425
426 f_Col(1,atom2) = f_Col(1,atom2) - ftm(1)
427 f_Col(2,atom2) = f_Col(2,atom2) - ftm(2)
428 f_Col(3,atom2) = f_Col(3,atom2) - ftm(3)
429
430 t_Row(1,atom1) = t_Row(1,atom1) + ttm_i(1)
431 t_Row(2,atom1) = t_Row(2,atom1) + ttm_i(2)
432 t_Row(3,atom1) = t_Row(3,atom1) + ttm_i(3)
433
434 t_Col(1,atom2) = t_Col(1,atom2) + ttm_j(1)
435 t_Col(2,atom2) = t_Col(2,atom2) + ttm_j(2)
436 t_Col(3,atom2) = t_Col(3,atom2) + ttm_j(3)
437 #else
438 f(1,atom1) = f(1,atom1) + ftm(1)
439 f(2,atom1) = f(2,atom1) + ftm(2)
440 f(3,atom1) = f(3,atom1) + ftm(3)
441
442 f(1,atom2) = f(1,atom2) - ftm(1)
443 f(2,atom2) = f(2,atom2) - ftm(2)
444 f(3,atom2) = f(3,atom2) - ftm(3)
445
446 t(1,atom1) = t(1,atom1) + ttm_i(1)
447 t(2,atom1) = t(2,atom1) + ttm_i(2)
448 t(3,atom1) = t(3,atom1) + ttm_i(3)
449
450 t(1,atom2) = t(1,atom2) + ttm_j(1)
451 t(2,atom2) = t(2,atom2) + ttm_j(2)
452 t(3,atom2) = t(3,atom2) + ttm_j(3)
453 #endif
454
455 #ifdef IS_MPI
456 id1 = AtomRowToGlobal(atom1)
457 id2 = AtomColToGlobal(atom2)
458 #else
459 id1 = atom1
460 id2 = atom2
461 #endif
462
463 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
464
465 fpair(1) = fpair(1) + ftm(1)
466 fpair(2) = fpair(2) + ftm(2)
467 fpair(3) = fpair(3) + ftm(3)
468
469 endif
470
471 return
472
473 end subroutine domultipolepair

Properties

Name Value
svn:keywords Author Id Revision Date