1 |
!! do_Forces.F90 |
2 |
!! module do_Forces |
3 |
!! Calculates Long Range forces. |
4 |
|
5 |
!! @author Charles F. Vardeman II |
6 |
!! @author Matthew Meineke |
7 |
!! @version $Id: do_Forces.F90,v 1.8 2003-03-10 19:26:45 chuckv Exp $, $Date: 2003-03-10 19:26:45 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $ |
8 |
|
9 |
|
10 |
|
11 |
module do_Forces |
12 |
use simulation |
13 |
use definitions |
14 |
use forceGlobals |
15 |
use atype_typedefs |
16 |
use neighborLists |
17 |
|
18 |
|
19 |
use lj_FF |
20 |
use sticky_FF |
21 |
use dp_FF |
22 |
use gb_FF |
23 |
|
24 |
#ifdef IS_MPI |
25 |
use mpiSimulation |
26 |
#endif |
27 |
implicit none |
28 |
PRIVATE |
29 |
|
30 |
public :: do_force_loop |
31 |
|
32 |
contains |
33 |
|
34 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
35 |
!-------------------------------------------------------------> |
36 |
subroutine do_force_loop(q,A,mu,u_l,f,t,tau,potE,do_pot,FFerror) |
37 |
!! Position array provided by C, dimensioned by getNlocal |
38 |
real ( kind = dp ), dimension(3,getNlocal()) :: q |
39 |
!! Rotation Matrix for each long range particle in simulation. |
40 |
real( kind = dp), dimension(9,getNlocal()) :: A |
41 |
|
42 |
!! Magnitude dipole moment |
43 |
real( kind = dp ), dimension(3,getNlocal()) :: mu |
44 |
!! Unit vectors for dipoles (lab frame) |
45 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
46 |
!! Force array provided by C, dimensioned by getNlocal |
47 |
real ( kind = dp ), dimension(3,getNlocal()) :: f |
48 |
!! Torsion array provided by C, dimensioned by getNlocal |
49 |
real( kind = dp ), dimension(3,getNlocal()) :: t |
50 |
|
51 |
!! Stress Tensor |
52 |
real( kind = dp), dimension(9) :: tau |
53 |
|
54 |
real ( kind = dp ) :: potE |
55 |
logical ( kind = 2) :: do_pot |
56 |
integer :: FFerror |
57 |
|
58 |
|
59 |
type(atype), pointer :: Atype_i |
60 |
type(atype), pointer :: Atype_j |
61 |
|
62 |
|
63 |
|
64 |
|
65 |
|
66 |
|
67 |
#ifdef IS_MPI |
68 |
real( kind = DP ) :: pot_local |
69 |
|
70 |
!! Local arrays needed for MPI |
71 |
|
72 |
#endif |
73 |
|
74 |
|
75 |
|
76 |
real( kind = DP ) :: pe |
77 |
logical :: update_nlist |
78 |
|
79 |
|
80 |
integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2 |
81 |
integer :: nlist |
82 |
integer :: j_start |
83 |
|
84 |
real( kind = DP ) :: r_ij, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp |
85 |
|
86 |
real( kind = DP ) :: rx_ij, ry_ij, rz_ij, rijsq |
87 |
real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut |
88 |
|
89 |
|
90 |
|
91 |
! a rig that need to be fixed. |
92 |
#ifdef IS_MPI |
93 |
real( kind = dp ) :: pe_local |
94 |
integer :: nlocal |
95 |
#endif |
96 |
integer :: nrow |
97 |
integer :: ncol |
98 |
integer :: natoms |
99 |
integer :: neighborListSize |
100 |
integer :: listerror |
101 |
!! should we calculate the stress tensor |
102 |
logical :: do_stress = .false. |
103 |
|
104 |
|
105 |
FFerror = 0 |
106 |
|
107 |
! Make sure we are properly initialized. |
108 |
if (.not. isFFInit) then |
109 |
write(default_error,*) "ERROR: lj_FF has not been properly initialized" |
110 |
FFerror = -1 |
111 |
return |
112 |
endif |
113 |
#ifdef IS_MPI |
114 |
if (.not. isMPISimSet()) then |
115 |
write(default_error,*) "ERROR: mpiSimulation has not been properly initialized" |
116 |
FFerror = -1 |
117 |
return |
118 |
endif |
119 |
#endif |
120 |
|
121 |
!! initialize local variables |
122 |
natoms = getNlocal() |
123 |
call getRcut(rcut,rcut2=rcutsq) |
124 |
call getRlist(rlist,rlistsq) |
125 |
|
126 |
!! Find ensemble |
127 |
if (isEnsemble("NPT")) do_stress = .true. |
128 |
!! set to wrap |
129 |
if (isPBC()) wrap = .true. |
130 |
|
131 |
|
132 |
|
133 |
|
134 |
!! See if we need to update neighbor lists |
135 |
call check(q,update_nlist) |
136 |
|
137 |
!--------------WARNING........................... |
138 |
! Zero variables, NOTE:::: Forces are zeroed in C |
139 |
! Zeroing them here could delete previously computed |
140 |
! Forces. |
141 |
!------------------------------------------------ |
142 |
call zero_module_variables() |
143 |
|
144 |
|
145 |
! communicate MPI positions |
146 |
#ifdef IS_MPI |
147 |
call gather(q,qRow,plan_row3d) |
148 |
call gather(q,qCol,plan_col3d) |
149 |
|
150 |
call gather(mu,muRow,plan_row3d) |
151 |
call gather(mu,muCol,plan_col3d) |
152 |
|
153 |
call gather(u_l,u_lRow,plan_row3d) |
154 |
call gather(u_l,u_lCol,plan_col3d) |
155 |
|
156 |
call gather(A,ARow,plan_row_rotation) |
157 |
call gather(A,ACol,plan_col_rotation) |
158 |
#endif |
159 |
|
160 |
|
161 |
#ifdef IS_MPI |
162 |
|
163 |
if (update_nlist) then |
164 |
|
165 |
! save current configuration, contruct neighbor list, |
166 |
! and calculate forces |
167 |
call save_neighborList(q) |
168 |
|
169 |
neighborListSize = getNeighborListSize() |
170 |
nlist = 0 |
171 |
|
172 |
nrow = getNrow(plan_row) |
173 |
ncol = getNcol(plan_col) |
174 |
nlocal = getNlocal() |
175 |
|
176 |
do i = 1, nrow |
177 |
point(i) = nlist + 1 |
178 |
Atype_i => identPtrListRow(i)%this |
179 |
|
180 |
inner: do j = 1, ncol |
181 |
Atype_j => identPtrListColumn(j)%this |
182 |
|
183 |
call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),& |
184 |
rxij,ryij,rzij,rijsq,r) |
185 |
|
186 |
! skip the loop if the atoms are identical |
187 |
if (mpi_cycle_jLoop(i,j)) cycle inner: |
188 |
|
189 |
if (rijsq < rlistsq) then |
190 |
|
191 |
nlist = nlist + 1 |
192 |
|
193 |
if (nlist > neighborListSize) then |
194 |
call expandList(listerror) |
195 |
if (listerror /= 0) then |
196 |
FFerror = -1 |
197 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
198 |
return |
199 |
end if |
200 |
endif |
201 |
|
202 |
list(nlist) = j |
203 |
|
204 |
|
205 |
if (rijsq < rcutsq) then |
206 |
call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij) |
207 |
endif |
208 |
endif |
209 |
enddo inner |
210 |
enddo |
211 |
|
212 |
point(nrow + 1) = nlist + 1 |
213 |
|
214 |
else !! (update) |
215 |
|
216 |
! use the list to find the neighbors |
217 |
do i = 1, nrow |
218 |
JBEG = POINT(i) |
219 |
JEND = POINT(i+1) - 1 |
220 |
! check thiat molecule i has neighbors |
221 |
if (jbeg .le. jend) then |
222 |
|
223 |
Atype_i => identPtrListRow(i)%this |
224 |
do jnab = jbeg, jend |
225 |
j = list(jnab) |
226 |
Atype_j = identPtrListColumn(j)%this |
227 |
call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),& |
228 |
rxij,ryij,rzij,rijsq,r) |
229 |
|
230 |
call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij) |
231 |
enddo |
232 |
endif |
233 |
enddo |
234 |
endif |
235 |
|
236 |
#else |
237 |
|
238 |
if (update_nlist) then |
239 |
|
240 |
! save current configuration, contruct neighbor list, |
241 |
! and calculate forces |
242 |
call save_neighborList(q) |
243 |
|
244 |
neighborListSize = getNeighborListSize() |
245 |
nlist = 0 |
246 |
|
247 |
|
248 |
do i = 1, natoms-1 |
249 |
point(i) = nlist + 1 |
250 |
Atype_i => identPtrList(i)%this |
251 |
|
252 |
inner: do j = i+1, natoms |
253 |
Atype_j => identPtrList(j)%this |
254 |
call get_interatomic_vector(i,j,q(:,i),q(:,j),& |
255 |
rxij,ryij,rzij,rijsq,r) |
256 |
|
257 |
if (rijsq < rlistsq) then |
258 |
|
259 |
nlist = nlist + 1 |
260 |
|
261 |
if (nlist > neighborListSize) then |
262 |
call expandList(listerror) |
263 |
if (listerror /= 0) then |
264 |
FFerror = -1 |
265 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
266 |
return |
267 |
end if |
268 |
endif |
269 |
|
270 |
list(nlist) = j |
271 |
|
272 |
|
273 |
if (rijsq < rcutsq) then |
274 |
call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij) |
275 |
endif |
276 |
endif |
277 |
enddo inner |
278 |
enddo |
279 |
|
280 |
point(natoms) = nlist + 1 |
281 |
|
282 |
else !! (update) |
283 |
|
284 |
! use the list to find the neighbors |
285 |
do i = 1, nrow |
286 |
JBEG = POINT(i) |
287 |
JEND = POINT(i+1) - 1 |
288 |
! check thiat molecule i has neighbors |
289 |
if (jbeg .le. jend) then |
290 |
|
291 |
Atype_i => identPtrList(i)%this |
292 |
do jnab = jbeg, jend |
293 |
j = list(jnab) |
294 |
Atype_j = identPtrList(j)%this |
295 |
call get_interatomic_vector(i,j,q(:,i),q(:,j),& |
296 |
rxij,ryij,rzij,rijsq,r) |
297 |
call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij) |
298 |
enddo |
299 |
endif |
300 |
enddo |
301 |
endif |
302 |
|
303 |
#endif |
304 |
|
305 |
|
306 |
#ifdef IS_MPI |
307 |
!! distribute all reaction field stuff (or anything for post-pair): |
308 |
call scatter(rflRow,rflTemp1,plan_row3d) |
309 |
call scatter(rflCol,rflTemp2,plan_col3d) |
310 |
do i = 1,nlocal |
311 |
rflTemp(1:3,i) = rflTemp1(1:3,i) + rflTemp2(1:3,i) |
312 |
end do |
313 |
#endif |
314 |
|
315 |
! This is the post-pair loop: |
316 |
#ifdef IS_MPI |
317 |
|
318 |
if (system_has_postpair_atoms) then |
319 |
do i = 1, nlocal |
320 |
Atype_i => identPtrListRow(i)%this |
321 |
call do_postpair(i, Atype_i) |
322 |
enddo |
323 |
endif |
324 |
|
325 |
#else |
326 |
|
327 |
if (system_has_postpair_atoms) then |
328 |
do i = 1, natoms |
329 |
Atype_i => identPtr(i)%this |
330 |
call do_postpair(i, Atype_i) |
331 |
enddo |
332 |
endif |
333 |
|
334 |
#endif |
335 |
|
336 |
|
337 |
|
338 |
|
339 |
#ifdef IS_MPI |
340 |
!!distribute forces |
341 |
|
342 |
call scatter(fRow,fTemp1,plan_row3d) |
343 |
call scatter(fCol,fTemp2,plan_col3d) |
344 |
|
345 |
|
346 |
do i = 1,nlocal |
347 |
fTemp(1:3,i) = fTemp1(1:3,i) + fTemp2(1:3,i) |
348 |
end do |
349 |
|
350 |
if (do_torque) then |
351 |
call scatter(tRow,tTemp1,plan_row3d) |
352 |
call scatter(tCol,tTemp2,plan_col3d) |
353 |
|
354 |
do i = 1,nlocal |
355 |
tTemp(1:3,i) = tTemp1(1:3,i) + tTemp2(1:3,i) |
356 |
end do |
357 |
endif |
358 |
|
359 |
if (do_pot) then |
360 |
! scatter/gather pot_row into the members of my column |
361 |
call scatter(eRow,eTemp,plan_row) |
362 |
|
363 |
! scatter/gather pot_local into all other procs |
364 |
! add resultant to get total pot |
365 |
do i = 1, nlocal |
366 |
pe_local = pe_local + eTemp(i) |
367 |
enddo |
368 |
|
369 |
eTemp = 0.0E0_DP |
370 |
call scatter(eCol,eTemp,plan_col) |
371 |
do i = 1, nlocal |
372 |
pe_local = pe_local + eTemp(i) |
373 |
enddo |
374 |
|
375 |
pe = pe_local |
376 |
endif |
377 |
#else |
378 |
! Copy local array into return array for c |
379 |
f = f+fTemp |
380 |
t = t+tTemp |
381 |
#endif |
382 |
|
383 |
potE = pe |
384 |
|
385 |
|
386 |
if (do_stress) then |
387 |
#ifdef IS_MPI |
388 |
mpi_allreduce = (tau,tauTemp,9,mpi_double_precision,mpi_sum, & |
389 |
mpi_comm_world,mpi_err) |
390 |
#else |
391 |
tau = tauTemp |
392 |
#endif |
393 |
endif |
394 |
|
395 |
end subroutine do_force_loop |
396 |
|
397 |
|
398 |
|
399 |
|
400 |
|
401 |
|
402 |
|
403 |
|
404 |
|
405 |
|
406 |
!! Calculate any pre-force loop components and update nlist if necessary. |
407 |
subroutine do_preForce(updateNlist) |
408 |
logical, intent(inout) :: updateNlist |
409 |
|
410 |
|
411 |
|
412 |
end subroutine do_preForce |
413 |
|
414 |
|
415 |
|
416 |
|
417 |
|
418 |
|
419 |
|
420 |
|
421 |
|
422 |
|
423 |
|
424 |
|
425 |
|
426 |
!! Calculate any post force loop components, i.e. reaction field, etc. |
427 |
subroutine do_postForce() |
428 |
|
429 |
|
430 |
|
431 |
end subroutine do_postForce |
432 |
|
433 |
|
434 |
|
435 |
|
436 |
|
437 |
|
438 |
|
439 |
|
440 |
|
441 |
|
442 |
|
443 |
|
444 |
|
445 |
|
446 |
|
447 |
|
448 |
|
449 |
subroutine do_pair(atype_i,atype_j,i,j,r_ij,rx_ij,ry_ij,rz_ij) |
450 |
type (atype ), pointer, intent(inout) :: atype_i |
451 |
type (atype ), pointer, intent(inout) :: atype_j |
452 |
integer :: i |
453 |
integer :: j |
454 |
real ( kind = dp ), intent(inout) :: rx_ij |
455 |
real ( kind = dp ), intent(inout) :: ry_ij |
456 |
real ( kind = dp ), intent(inout) :: rz_ij |
457 |
|
458 |
|
459 |
real( kind = dp ) :: fx = 0.0_dp |
460 |
real( kind = dp ) :: fy = 0.0_dp |
461 |
real( kind = dp ) :: fz = 0.0_dp |
462 |
|
463 |
real( kind = dp ) :: drdx = 0.0_dp |
464 |
real( kind = dp ) :: drdy = 0.0_dp |
465 |
real( kind = dp ) :: drdz = 0.0_dp |
466 |
|
467 |
|
468 |
#ifdef IS_MPI |
469 |
|
470 |
if (Atype_i%is_LJ .and. Atype_j%is_LJ) then |
471 |
call getLJForce(r,pot,dudr,ljAtype_i,ljAtype_j,fx,fy,fz) |
472 |
endif |
473 |
|
474 |
if (Atype_i%is_dp .and. Atype_j%is_dp) then |
475 |
|
476 |
call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, & |
477 |
ulRow(:,i), ulCol(:,j), rt, rrf, pot) |
478 |
|
479 |
if (do_reaction_field) then |
480 |
call accumulate_rf(i, j, r_ij, rflRow(:,i), rflCol(:j), & |
481 |
ulRow(:i), ulCol(:,j), rt, rrf) |
482 |
endif |
483 |
|
484 |
endif |
485 |
|
486 |
if (Atype_i%is_sticky .and. Atype_j%is_sticky) then |
487 |
call getstickyforce(r, pot, dudr, Atype_i, Atype_j) |
488 |
endif |
489 |
|
490 |
#else |
491 |
|
492 |
if (Atype_i%is_LJ .and. Atype_j%is_LJ) then |
493 |
call getLJForce(r,pot,dudr,ljAtype_i,ljAtype_j,fx,fy,fz) |
494 |
endif |
495 |
|
496 |
if (Atype_i%is_dp .and. Atype_j%is_dp) then |
497 |
call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, & |
498 |
ul(:,i), ul(:,j), rt, rrf, pot) |
499 |
|
500 |
if (do_reaction_field) then |
501 |
call accumulate_rf(i, j, r_ij, rfl(:,i), rfl(:j), & |
502 |
ul(:,i), ul(:,j), rt, rrf) |
503 |
endif |
504 |
|
505 |
endif |
506 |
|
507 |
if (Atype_i%is_sticky .and. Atype_j%is_sticky) then |
508 |
call getstickyforce(r,pot,dudr, Atype_i, Atype_j) |
509 |
endif |
510 |
|
511 |
#endif |
512 |
|
513 |
|
514 |
#ifdef IS_MPI |
515 |
eRow(i) = eRow(i) + pot*0.5 |
516 |
eCol(i) = eCol(i) + pot*0.5 |
517 |
#else |
518 |
pe = pe + pot |
519 |
#endif |
520 |
|
521 |
drdx = -rxij / r |
522 |
drdy = -ryij / r |
523 |
drdz = -rzij / r |
524 |
|
525 |
fx = dudr * drdx |
526 |
fy = dudr * drdy |
527 |
fz = dudr * drdz |
528 |
|
529 |
#ifdef IS_MPI |
530 |
fCol(1,j) = fCol(1,j) - fx |
531 |
fCol(2,j) = fCol(2,j) - fy |
532 |
fCol(3,j) = fCol(3,j) - fz |
533 |
|
534 |
fRow(1,j) = fRow(1,j) + fx |
535 |
fRow(2,j) = fRow(2,j) + fy |
536 |
fRow(3,j) = fRow(3,j) + fz |
537 |
#else |
538 |
fTemp(1,j) = fTemp(1,j) - fx |
539 |
fTemp(2,j) = fTemp(2,j) - fy |
540 |
fTemp(3,j) = fTemp(3,j) - fz |
541 |
fTemp(1,i) = fTemp(1,i) + fx |
542 |
fTemp(2,i) = fTemp(2,i) + fy |
543 |
fTemp(3,i) = fTemp(3,i) + fz |
544 |
#endif |
545 |
|
546 |
if (do_stress) then |
547 |
tauTemp(1) = tauTemp(1) + fx * rxij |
548 |
tauTemp(2) = tauTemp(2) + fx * ryij |
549 |
tauTemp(3) = tauTemp(3) + fx * rzij |
550 |
tauTemp(4) = tauTemp(4) + fy * rxij |
551 |
tauTemp(5) = tauTemp(5) + fy * ryij |
552 |
tauTemp(6) = tauTemp(6) + fy * rzij |
553 |
tauTemp(7) = tauTemp(7) + fz * rxij |
554 |
tauTemp(8) = tauTemp(8) + fz * ryij |
555 |
tauTemp(9) = tauTemp(9) + fz * rzij |
556 |
endif |
557 |
|
558 |
|
559 |
|
560 |
end subroutine do_pair |
561 |
|
562 |
|
563 |
|
564 |
|
565 |
|
566 |
|
567 |
|
568 |
|
569 |
|
570 |
|
571 |
|
572 |
|
573 |
|
574 |
|
575 |
|
576 |
|
577 |
subroutine get_interatomic_vector(q_i,q_j,rx_ij,ry_ij,rz_ij,r_sq,r_ij) |
578 |
!---------------- Arguments------------------------------- |
579 |
!! index i |
580 |
|
581 |
!! Position array |
582 |
real (kind = dp), dimension(3) :: q_i |
583 |
real (kind = dp), dimension(3) :: q_j |
584 |
!! x component of vector between i and j |
585 |
real ( kind = dp ), intent(out) :: rx_ij |
586 |
!! y component of vector between i and j |
587 |
real ( kind = dp ), intent(out) :: ry_ij |
588 |
!! z component of vector between i and j |
589 |
real ( kind = dp ), intent(out) :: rz_ij |
590 |
!! magnitude of r squared |
591 |
real ( kind = dp ), intent(out) :: r_sq |
592 |
!! magnitude of vector r between atoms i and j. |
593 |
real ( kind = dp ), intent(out) :: r_ij |
594 |
!! wrap into periodic box. |
595 |
logical, intent(in) :: wrap |
596 |
|
597 |
!--------------- Local Variables--------------------------- |
598 |
!! Distance between i and j |
599 |
real( kind = dp ) :: d(3) |
600 |
!---------------- END DECLARATIONS------------------------- |
601 |
|
602 |
|
603 |
! Find distance between i and j |
604 |
d(1:3) = q_i(1:3) - q_j(1:3) |
605 |
|
606 |
! Wrap back into periodic box if necessary |
607 |
if ( wrap ) then |
608 |
d(1:3) = d(1:3) - thisSim%box(1:3) * sign(1.0_dp,thisSim%box(1:3)) * & |
609 |
int(abs(d(1:3)/thisSim%box(1:3) + 0.5_dp) |
610 |
end if |
611 |
|
612 |
! Find Magnitude of the vector |
613 |
r_sq = dot_product(d,d) |
614 |
r_ij = sqrt(r_sq) |
615 |
|
616 |
! Set each component for force calculation |
617 |
rx_ij = d(1) |
618 |
ry_ij = d(2) |
619 |
rz_ij = d(3) |
620 |
|
621 |
|
622 |
end subroutine get_interatomic_vector |
623 |
|
624 |
subroutine zero_module_variables() |
625 |
|
626 |
#ifndef IS_MPI |
627 |
|
628 |
pe = 0.0E0_DP |
629 |
tauTemp = 0.0_dp |
630 |
fTemp = 0.0_dp |
631 |
tTemp = 0.0_dp |
632 |
#else |
633 |
qRow = 0.0_dp |
634 |
qCol = 0.0_dp |
635 |
|
636 |
muRow = 0.0_dp |
637 |
muCol = 0.0_dp |
638 |
|
639 |
u_lRow = 0.0_dp |
640 |
u_lCol = 0.0_dp |
641 |
|
642 |
ARow = 0.0_dp |
643 |
ACol = 0.0_dp |
644 |
|
645 |
fRow = 0.0_dp |
646 |
fCol = 0.0_dp |
647 |
|
648 |
|
649 |
tRow = 0.0_dp |
650 |
tCol = 0.0_dp |
651 |
|
652 |
|
653 |
|
654 |
eRow = 0.0_dp |
655 |
eCol = 0.0_dp |
656 |
eTemp = 0.0_dp |
657 |
#endif |
658 |
|
659 |
end subroutine zero_module_variables |
660 |
|
661 |
|
662 |
!! Function to properly build neighbor lists in MPI using newtons 3rd law. |
663 |
!! We don't want 2 processors doing the same i j pair twice. |
664 |
!! Also checks to see if i and j are the same particle. |
665 |
function checkExcludes(atom1,atom2) result(do_cycle) |
666 |
!--------------- Arguments-------------------------- |
667 |
! Index i |
668 |
integer,intent(in) :: atom1 |
669 |
! Index j |
670 |
integer,intent(in), optional :: atom2 |
671 |
! Result do_cycle |
672 |
logical :: do_cycle |
673 |
!--------------- Local variables-------------------- |
674 |
integer :: tag_i |
675 |
integer :: tag_j |
676 |
integer :: i |
677 |
!--------------- END DECLARATIONS------------------ |
678 |
do_cycle = .false. |
679 |
|
680 |
#ifdef IS_MPI |
681 |
tag_i = tagRow(atom1) |
682 |
#else |
683 |
tag_i = tag(atom1) |
684 |
#endif |
685 |
|
686 |
!! Check global excludes first |
687 |
if (.not. present(atom2)) then |
688 |
do i = 1,nGlobalExcludes |
689 |
if (excludeGlobal(i) == tag_i) then |
690 |
do_cycle = .true. |
691 |
return |
692 |
end if |
693 |
end do |
694 |
return !! return after checking globals |
695 |
end if |
696 |
|
697 |
!! we return if j not present here. |
698 |
tag_j = tagColumn(j) |
699 |
|
700 |
|
701 |
|
702 |
if (tag_i == tag_j) then |
703 |
do_cycle = .true. |
704 |
return |
705 |
end if |
706 |
|
707 |
if (tag_i < tag_j) then |
708 |
if (mod(tag_i + tag_j,2) == 0) do_cycle = .true. |
709 |
return |
710 |
else |
711 |
if (mod(tag_i + tag_j,2) == 1) do_cycle = .true. |
712 |
endif |
713 |
|
714 |
|
715 |
|
716 |
do i = 1, nLocalExcludes |
717 |
if (tag_i = excludes(1,i) .and. excludes(2,i) < 0) then |
718 |
do_cycle = .true. |
719 |
return |
720 |
end if |
721 |
end do |
722 |
|
723 |
|
724 |
end function checkExcludes |
725 |
|
726 |
|
727 |
end module do_Forces |