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.24 2003-07-23 22:13:59 chuckv Exp $, $Date: 2003-07-23 22:13:59 $, $Name: not supported by cvs2svn $, $Revision: 1.24 $ |
8 |
|
9 |
module do_Forces |
10 |
use force_globals |
11 |
use simulation |
12 |
use definitions |
13 |
use atype_module |
14 |
use neighborLists |
15 |
use lj |
16 |
use sticky_pair |
17 |
use dipole_dipole |
18 |
use reaction_field |
19 |
use gb_pair |
20 |
use vector_class |
21 |
#ifdef IS_MPI |
22 |
use mpiSimulation |
23 |
#endif |
24 |
|
25 |
implicit none |
26 |
PRIVATE |
27 |
|
28 |
#define __FORTRAN90 |
29 |
#include "fForceField.h" |
30 |
|
31 |
logical, save :: do_forces_initialized = .false., haveRlist = .false. |
32 |
logical, save :: havePolicies = .false. |
33 |
logical, save :: FF_uses_LJ |
34 |
logical, save :: FF_uses_sticky |
35 |
logical, save :: FF_uses_dipoles |
36 |
logical, save :: FF_uses_RF |
37 |
logical, save :: FF_uses_GB |
38 |
logical, save :: FF_uses_EAM |
39 |
|
40 |
real(kind=dp), save :: rlist, rlistsq |
41 |
|
42 |
public :: init_FF |
43 |
public :: do_force_loop |
44 |
public :: setRlistDF |
45 |
|
46 |
contains |
47 |
|
48 |
subroutine setRlistDF( this_rlist ) |
49 |
|
50 |
real(kind=dp) :: this_rlist |
51 |
|
52 |
rlist = this_rlist |
53 |
rlistsq = rlist * rlist |
54 |
|
55 |
haveRlist = .true. |
56 |
if( havePolicies ) do_forces_initialized = .true. |
57 |
|
58 |
end subroutine setRlistDF |
59 |
|
60 |
subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat) |
61 |
|
62 |
integer, intent(in) :: LJMIXPOLICY |
63 |
logical, intent(in) :: use_RF_c |
64 |
|
65 |
integer, intent(out) :: thisStat |
66 |
integer :: my_status, nMatches |
67 |
integer, pointer :: MatchList(:) => null() |
68 |
real(kind=dp) :: rcut, rrf, rt, dielect |
69 |
|
70 |
!! assume things are copacetic, unless they aren't |
71 |
thisStat = 0 |
72 |
|
73 |
!! Fortran's version of a cast: |
74 |
FF_uses_RF = use_RF_c |
75 |
|
76 |
!! init_FF is called *after* all of the atom types have been |
77 |
!! defined in atype_module using the new_atype subroutine. |
78 |
!! |
79 |
!! this will scan through the known atypes and figure out what |
80 |
!! interactions are used by the force field. |
81 |
|
82 |
FF_uses_LJ = .false. |
83 |
FF_uses_sticky = .false. |
84 |
FF_uses_dipoles = .false. |
85 |
FF_uses_GB = .false. |
86 |
FF_uses_EAM = .false. |
87 |
|
88 |
call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList) |
89 |
if (nMatches .gt. 0) FF_uses_LJ = .true. |
90 |
|
91 |
call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList) |
92 |
if (nMatches .gt. 0) FF_uses_dipoles = .true. |
93 |
|
94 |
call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, & |
95 |
MatchList) |
96 |
if (nMatches .gt. 0) FF_uses_Sticky = .true. |
97 |
|
98 |
call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList) |
99 |
if (nMatches .gt. 0) FF_uses_GB = .true. |
100 |
|
101 |
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
102 |
if (nMatches .gt. 0) FF_uses_EAM = .true. |
103 |
|
104 |
!! check to make sure the FF_uses_RF setting makes sense |
105 |
|
106 |
if (FF_uses_dipoles) then |
107 |
if (FF_uses_RF) then |
108 |
dielect = getDielect() |
109 |
call initialize_rf(dielect) |
110 |
endif |
111 |
else |
112 |
if (FF_uses_RF) then |
113 |
write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' |
114 |
thisStat = -1 |
115 |
return |
116 |
endif |
117 |
endif |
118 |
|
119 |
if (FF_uses_LJ) then |
120 |
|
121 |
select case (LJMIXPOLICY) |
122 |
case (LB_MIXING_RULE) |
123 |
call init_lj_FF(LB_MIXING_RULE, my_status) |
124 |
case (EXPLICIT_MIXING_RULE) |
125 |
call init_lj_FF(EXPLICIT_MIXING_RULE, my_status) |
126 |
case default |
127 |
write(default_error,*) 'unknown LJ Mixing Policy!' |
128 |
thisStat = -1 |
129 |
return |
130 |
end select |
131 |
if (my_status /= 0) then |
132 |
thisStat = -1 |
133 |
return |
134 |
end if |
135 |
endif |
136 |
|
137 |
if (FF_uses_sticky) then |
138 |
call check_sticky_FF(my_status) |
139 |
if (my_status /= 0) then |
140 |
thisStat = -1 |
141 |
return |
142 |
end if |
143 |
endif |
144 |
|
145 |
if (FF_uses_GB) then |
146 |
call check_gb_pair_FF(my_status) |
147 |
if (my_status .ne. 0) then |
148 |
thisStat = -1 |
149 |
return |
150 |
endif |
151 |
endif |
152 |
|
153 |
if (FF_uses_GB .and. FF_uses_LJ) then |
154 |
endif |
155 |
if (.not. do_forces_initialized) then |
156 |
!! Create neighbor lists |
157 |
call expandNeighborList(getNlocal(), my_status) |
158 |
if (my_Status /= 0) then |
159 |
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
160 |
thisStat = -1 |
161 |
return |
162 |
endif |
163 |
endif |
164 |
|
165 |
havePolicies = .true. |
166 |
if( haveRlist ) do_forces_initialized = .true. |
167 |
|
168 |
end subroutine init_FF |
169 |
|
170 |
|
171 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
172 |
!-------------------------------------------------------------> |
173 |
subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, & |
174 |
error) |
175 |
!! Position array provided by C, dimensioned by getNlocal |
176 |
real ( kind = dp ), dimension(3,getNlocal()) :: q |
177 |
!! Rotation Matrix for each long range particle in simulation. |
178 |
real( kind = dp), dimension(9,getNlocal()) :: A |
179 |
!! Unit vectors for dipoles (lab frame) |
180 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
181 |
!! Force array provided by C, dimensioned by getNlocal |
182 |
real ( kind = dp ), dimension(3,getNlocal()) :: f |
183 |
!! Torsion array provided by C, dimensioned by getNlocal |
184 |
real( kind = dp ), dimension(3,getNlocal()) :: t |
185 |
!! Stress Tensor |
186 |
real( kind = dp), dimension(9) :: tau |
187 |
real ( kind = dp ) :: pot |
188 |
logical ( kind = 2) :: do_pot_c, do_stress_c |
189 |
logical :: do_pot |
190 |
logical :: do_stress |
191 |
#ifdef IS_MPI |
192 |
real( kind = DP ) :: pot_local |
193 |
integer :: nrow |
194 |
integer :: ncol |
195 |
#endif |
196 |
integer :: nlocal |
197 |
integer :: natoms |
198 |
logical :: update_nlist |
199 |
integer :: i, j, jbeg, jend, jnab |
200 |
integer :: nlist |
201 |
real( kind = DP ) :: rijsq |
202 |
real(kind=dp),dimension(3) :: d |
203 |
real(kind=dp) :: rfpot, mu_i, virial |
204 |
integer :: me_i |
205 |
logical :: is_dp_i |
206 |
integer :: neighborListSize |
207 |
integer :: listerror, error |
208 |
integer :: localError |
209 |
|
210 |
real(kind=dp) :: listSkin = 1.0 |
211 |
|
212 |
|
213 |
!! initialize local variables |
214 |
|
215 |
#ifdef IS_MPI |
216 |
pot_local = 0.0_dp |
217 |
nlocal = getNlocal() |
218 |
nrow = getNrow(plan_row) |
219 |
ncol = getNcol(plan_col) |
220 |
#else |
221 |
nlocal = getNlocal() |
222 |
natoms = nlocal |
223 |
#endif |
224 |
|
225 |
call check_initialization(localError) |
226 |
if ( localError .ne. 0 ) then |
227 |
error = -1 |
228 |
return |
229 |
end if |
230 |
call zero_work_arrays() |
231 |
|
232 |
do_pot = do_pot_c |
233 |
do_stress = do_stress_c |
234 |
|
235 |
! Gather all information needed by all force loops: |
236 |
|
237 |
#ifdef IS_MPI |
238 |
|
239 |
call gather(q,q_Row,plan_row3d) |
240 |
call gather(q,q_Col,plan_col3d) |
241 |
|
242 |
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
243 |
call gather(u_l,u_l_Row,plan_row3d) |
244 |
call gather(u_l,u_l_Col,plan_col3d) |
245 |
|
246 |
call gather(A,A_Row,plan_row_rotation) |
247 |
call gather(A,A_Col,plan_col_rotation) |
248 |
endif |
249 |
|
250 |
#endif |
251 |
|
252 |
if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then |
253 |
!! See if we need to update neighbor lists |
254 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
255 |
!! if_mpi_gather_stuff_for_prepair |
256 |
!! do_prepair_loop_if_needed |
257 |
!! if_mpi_scatter_stuff_from_prepair |
258 |
!! if_mpi_gather_stuff_from_prepair_to_main_loop |
259 |
|
260 |
!--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>> |
261 |
#ifdef IS_MPI |
262 |
|
263 |
if (update_nlist) then |
264 |
|
265 |
!! save current configuration, construct neighbor list, |
266 |
!! and calculate forces |
267 |
call saveNeighborList(nlocal, q) |
268 |
|
269 |
neighborListSize = size(list) |
270 |
nlist = 0 |
271 |
|
272 |
do i = 1, nrow |
273 |
point(i) = nlist + 1 |
274 |
|
275 |
prepair_inner: do j = 1, ncol |
276 |
|
277 |
if (skipThisPair(i,j)) cycle prepair_inner |
278 |
|
279 |
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
280 |
|
281 |
if (rijsq < rlistsq) then |
282 |
|
283 |
nlist = nlist + 1 |
284 |
|
285 |
if (nlist > neighborListSize) then |
286 |
call expandNeighborList(nlocal, listerror) |
287 |
if (listerror /= 0) then |
288 |
error = -1 |
289 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
290 |
return |
291 |
end if |
292 |
neighborListSize = size(list) |
293 |
endif |
294 |
|
295 |
list(nlist) = j |
296 |
call do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot_local) |
297 |
endif |
298 |
enddo prepair_inner |
299 |
enddo |
300 |
|
301 |
point(nrow + 1) = nlist + 1 |
302 |
|
303 |
else !! (of update_check) |
304 |
|
305 |
! use the list to find the neighbors |
306 |
do i = 1, nrow |
307 |
JBEG = POINT(i) |
308 |
JEND = POINT(i+1) - 1 |
309 |
! check thiat molecule i has neighbors |
310 |
if (jbeg .le. jend) then |
311 |
|
312 |
do jnab = jbeg, jend |
313 |
j = list(jnab) |
314 |
|
315 |
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
316 |
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
317 |
u_l, A, f, t, pot_local) |
318 |
|
319 |
enddo |
320 |
endif |
321 |
enddo |
322 |
endif |
323 |
|
324 |
#else |
325 |
|
326 |
if (update_nlist) then |
327 |
|
328 |
! save current configuration, contruct neighbor list, |
329 |
! and calculate forces |
330 |
call saveNeighborList(natoms, q) |
331 |
|
332 |
neighborListSize = size(list) |
333 |
|
334 |
nlist = 0 |
335 |
|
336 |
do i = 1, natoms-1 |
337 |
point(i) = nlist + 1 |
338 |
|
339 |
prepair_inner: do j = i+1, natoms |
340 |
|
341 |
if (skipThisPair(i,j)) cycle prepair_inner |
342 |
|
343 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
344 |
|
345 |
|
346 |
if (rijsq < rlistsq) then |
347 |
|
348 |
nlist = nlist + 1 |
349 |
|
350 |
if (nlist > neighborListSize) then |
351 |
call expandNeighborList(natoms, listerror) |
352 |
if (listerror /= 0) then |
353 |
error = -1 |
354 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
355 |
return |
356 |
end if |
357 |
neighborListSize = size(list) |
358 |
endif |
359 |
|
360 |
list(nlist) = j |
361 |
|
362 |
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
363 |
u_l, A, f, t, pot) |
364 |
|
365 |
endif |
366 |
enddo prepair_inner |
367 |
enddo |
368 |
|
369 |
point(natoms) = nlist + 1 |
370 |
|
371 |
else !! (update) |
372 |
|
373 |
! use the list to find the neighbors |
374 |
do i = 1, natoms-1 |
375 |
JBEG = POINT(i) |
376 |
JEND = POINT(i+1) - 1 |
377 |
! check thiat molecule i has neighbors |
378 |
if (jbeg .le. jend) then |
379 |
|
380 |
do jnab = jbeg, jend |
381 |
j = list(jnab) |
382 |
|
383 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
384 |
call do_prepair(i, j, rijsq, d, do_pot, do_stress, & |
385 |
u_l, A, f, t, pot) |
386 |
|
387 |
enddo |
388 |
endif |
389 |
enddo |
390 |
endif |
391 |
#endif |
392 |
!! Do rest of preforce calculations |
393 |
call do_preforce(nlocal,pot) |
394 |
else |
395 |
!! See if we need to update neighbor lists for non pre-pair |
396 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
397 |
endif |
398 |
|
399 |
|
400 |
|
401 |
|
402 |
|
403 |
!---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>>>>>>>>>>>>>>>> |
404 |
|
405 |
|
406 |
|
407 |
|
408 |
|
409 |
#ifdef IS_MPI |
410 |
|
411 |
if (update_nlist) then |
412 |
|
413 |
!! save current configuration, construct neighbor list, |
414 |
!! and calculate forces |
415 |
call saveNeighborList(nlocal, q) |
416 |
|
417 |
neighborListSize = size(list) |
418 |
nlist = 0 |
419 |
|
420 |
do i = 1, nrow |
421 |
point(i) = nlist + 1 |
422 |
|
423 |
inner: do j = 1, ncol |
424 |
|
425 |
if (skipThisPair(i,j)) cycle inner |
426 |
|
427 |
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
428 |
|
429 |
if (rijsq < rlistsq) then |
430 |
|
431 |
nlist = nlist + 1 |
432 |
|
433 |
if (nlist > neighborListSize) then |
434 |
call expandNeighborList(nlocal, listerror) |
435 |
if (listerror /= 0) then |
436 |
error = -1 |
437 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
438 |
return |
439 |
end if |
440 |
neighborListSize = size(list) |
441 |
endif |
442 |
|
443 |
list(nlist) = j |
444 |
|
445 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
446 |
u_l, A, f, t, pot_local) |
447 |
|
448 |
endif |
449 |
enddo inner |
450 |
enddo |
451 |
|
452 |
point(nrow + 1) = nlist + 1 |
453 |
|
454 |
else !! (of update_check) |
455 |
|
456 |
! use the list to find the neighbors |
457 |
do i = 1, nrow |
458 |
JBEG = POINT(i) |
459 |
JEND = POINT(i+1) - 1 |
460 |
! check thiat molecule i has neighbors |
461 |
if (jbeg .le. jend) then |
462 |
|
463 |
do jnab = jbeg, jend |
464 |
j = list(jnab) |
465 |
|
466 |
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
467 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
468 |
u_l, A, f, t, pot_local) |
469 |
|
470 |
enddo |
471 |
endif |
472 |
enddo |
473 |
endif |
474 |
|
475 |
#else |
476 |
|
477 |
if (update_nlist) then |
478 |
|
479 |
! save current configuration, contruct neighbor list, |
480 |
! and calculate forces |
481 |
call saveNeighborList(natoms, q) |
482 |
|
483 |
neighborListSize = size(list) |
484 |
|
485 |
nlist = 0 |
486 |
|
487 |
do i = 1, natoms-1 |
488 |
point(i) = nlist + 1 |
489 |
|
490 |
inner: do j = i+1, natoms |
491 |
|
492 |
if (skipThisPair(i,j)) cycle inner |
493 |
|
494 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
495 |
|
496 |
|
497 |
if (rijsq < rlistsq) then |
498 |
|
499 |
nlist = nlist + 1 |
500 |
|
501 |
if (nlist > neighborListSize) then |
502 |
call expandNeighborList(natoms, listerror) |
503 |
if (listerror /= 0) then |
504 |
error = -1 |
505 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
506 |
return |
507 |
end if |
508 |
neighborListSize = size(list) |
509 |
endif |
510 |
|
511 |
list(nlist) = j |
512 |
|
513 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
514 |
u_l, A, f, t, pot) |
515 |
|
516 |
endif |
517 |
enddo inner |
518 |
enddo |
519 |
|
520 |
point(natoms) = nlist + 1 |
521 |
|
522 |
else !! (update) |
523 |
|
524 |
! use the list to find the neighbors |
525 |
do i = 1, natoms-1 |
526 |
JBEG = POINT(i) |
527 |
JEND = POINT(i+1) - 1 |
528 |
! check thiat molecule i has neighbors |
529 |
if (jbeg .le. jend) then |
530 |
|
531 |
do jnab = jbeg, jend |
532 |
j = list(jnab) |
533 |
|
534 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
535 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
536 |
u_l, A, f, t, pot) |
537 |
|
538 |
enddo |
539 |
endif |
540 |
enddo |
541 |
endif |
542 |
|
543 |
#endif |
544 |
|
545 |
! phew, done with main loop. |
546 |
|
547 |
#ifdef IS_MPI |
548 |
!!distribute forces |
549 |
|
550 |
f_temp = 0.0_dp |
551 |
call scatter(f_Row,f_temp,plan_row3d) |
552 |
do i = 1,nlocal |
553 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
554 |
end do |
555 |
|
556 |
f_temp = 0.0_dp |
557 |
call scatter(f_Col,f_temp,plan_col3d) |
558 |
do i = 1,nlocal |
559 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
560 |
end do |
561 |
|
562 |
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
563 |
t_temp = 0.0_dp |
564 |
call scatter(t_Row,t_temp,plan_row3d) |
565 |
do i = 1,nlocal |
566 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
567 |
end do |
568 |
t_temp = 0.0_dp |
569 |
call scatter(t_Col,t_temp,plan_col3d) |
570 |
|
571 |
do i = 1,nlocal |
572 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
573 |
end do |
574 |
endif |
575 |
|
576 |
if (do_pot) then |
577 |
! scatter/gather pot_row into the members of my column |
578 |
call scatter(pot_Row, pot_Temp, plan_row) |
579 |
|
580 |
! scatter/gather pot_local into all other procs |
581 |
! add resultant to get total pot |
582 |
do i = 1, nlocal |
583 |
pot_local = pot_local + pot_Temp(i) |
584 |
enddo |
585 |
|
586 |
pot_Temp = 0.0_DP |
587 |
|
588 |
call scatter(pot_Col, pot_Temp, plan_col) |
589 |
do i = 1, nlocal |
590 |
pot_local = pot_local + pot_Temp(i) |
591 |
enddo |
592 |
|
593 |
endif |
594 |
#endif |
595 |
|
596 |
if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then |
597 |
|
598 |
if (FF_uses_RF .and. SimUsesRF()) then |
599 |
|
600 |
#ifdef IS_MPI |
601 |
call scatter(rf_Row,rf,plan_row3d) |
602 |
call scatter(rf_Col,rf_Temp,plan_col3d) |
603 |
do i = 1,nlocal |
604 |
rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i) |
605 |
end do |
606 |
#endif |
607 |
|
608 |
do i = 1, getNlocal() |
609 |
|
610 |
rfpot = 0.0_DP |
611 |
#ifdef IS_MPI |
612 |
me_i = atid_row(i) |
613 |
#else |
614 |
me_i = atid(i) |
615 |
#endif |
616 |
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
617 |
if ( is_DP_i ) then |
618 |
call getElementProperty(atypes, me_i, "dipole_moment", mu_i) |
619 |
!! The reaction field needs to include a self contribution |
620 |
!! to the field: |
621 |
call accumulate_self_rf(i, mu_i, u_l) |
622 |
!! Get the reaction field contribution to the |
623 |
!! potential and torques: |
624 |
call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot) |
625 |
#ifdef IS_MPI |
626 |
pot_local = pot_local + rfpot |
627 |
#else |
628 |
pot = pot + rfpot |
629 |
|
630 |
#endif |
631 |
endif |
632 |
enddo |
633 |
endif |
634 |
endif |
635 |
|
636 |
|
637 |
#ifdef IS_MPI |
638 |
|
639 |
if (do_pot) then |
640 |
pot = pot + pot_local |
641 |
!! we assume the c code will do the allreduce to get the total potential |
642 |
!! we could do it right here if we needed to... |
643 |
endif |
644 |
|
645 |
if (do_stress) then |
646 |
call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, & |
647 |
mpi_comm_world,mpi_err) |
648 |
call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, & |
649 |
mpi_comm_world,mpi_err) |
650 |
endif |
651 |
|
652 |
#else |
653 |
|
654 |
if (do_stress) then |
655 |
tau = tau_Temp |
656 |
virial = virial_Temp |
657 |
endif |
658 |
|
659 |
#endif |
660 |
|
661 |
end subroutine do_force_loop |
662 |
|
663 |
subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
664 |
|
665 |
real( kind = dp ) :: pot |
666 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
667 |
real (kind=dp), dimension(9,getNlocal()) :: A |
668 |
real (kind=dp), dimension(3,getNlocal()) :: f |
669 |
real (kind=dp), dimension(3,getNlocal()) :: t |
670 |
|
671 |
logical, intent(inout) :: do_pot, do_stress |
672 |
integer, intent(in) :: i, j |
673 |
real ( kind = dp ), intent(inout) :: rijsq |
674 |
real ( kind = dp ) :: r |
675 |
real ( kind = dp ), intent(inout) :: d(3) |
676 |
logical :: is_LJ_i, is_LJ_j |
677 |
logical :: is_DP_i, is_DP_j |
678 |
logical :: is_GB_i, is_GB_j |
679 |
logical :: is_EAM_i,is_EAM_j |
680 |
logical :: is_Sticky_i, is_Sticky_j |
681 |
integer :: me_i, me_j |
682 |
|
683 |
r = sqrt(rijsq) |
684 |
|
685 |
#ifdef IS_MPI |
686 |
if (tagRow(i) .eq. tagColumn(j)) then |
687 |
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
688 |
endif |
689 |
|
690 |
me_i = atid_row(i) |
691 |
me_j = atid_col(j) |
692 |
|
693 |
#else |
694 |
|
695 |
me_i = atid(i) |
696 |
me_j = atid(j) |
697 |
|
698 |
#endif |
699 |
|
700 |
if (FF_uses_LJ .and. SimUsesLJ()) then |
701 |
call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i) |
702 |
call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j) |
703 |
|
704 |
if ( is_LJ_i .and. is_LJ_j ) & |
705 |
call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
706 |
endif |
707 |
|
708 |
if (FF_uses_dipoles .and. SimUsesDipoles()) then |
709 |
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
710 |
call getElementProperty(atypes, me_j, "is_DP", is_DP_j) |
711 |
|
712 |
if ( is_DP_i .and. is_DP_j ) then |
713 |
|
714 |
call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, & |
715 |
do_pot, do_stress) |
716 |
if (FF_uses_RF .and. SimUsesRF()) then |
717 |
call accumulate_rf(i, j, r, u_l) |
718 |
call rf_correct_forces(i, j, d, r, u_l, f, do_stress) |
719 |
endif |
720 |
|
721 |
endif |
722 |
endif |
723 |
|
724 |
if (FF_uses_Sticky .and. SimUsesSticky()) then |
725 |
|
726 |
call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i) |
727 |
call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j) |
728 |
|
729 |
if ( is_Sticky_i .and. is_Sticky_j ) then |
730 |
call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, & |
731 |
do_pot, do_stress) |
732 |
endif |
733 |
endif |
734 |
|
735 |
|
736 |
if (FF_uses_GB .and. SimUsesGB()) then |
737 |
|
738 |
call getElementProperty(atypes, me_i, "is_GB", is_GB_i) |
739 |
call getElementProperty(atypes, me_j, "is_GB", is_GB_j) |
740 |
|
741 |
if ( is_GB_i .and. is_GB_j ) then |
742 |
call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, & |
743 |
do_pot, do_stress) |
744 |
endif |
745 |
endif |
746 |
|
747 |
|
748 |
|
749 |
if (FF_uses_EAM .and. SimUsesEAM()) then |
750 |
call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i) |
751 |
call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j) |
752 |
|
753 |
if ( is_EAM_i .and. is_EAM_j ) & |
754 |
call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
755 |
endif |
756 |
|
757 |
|
758 |
|
759 |
|
760 |
end subroutine do_pair |
761 |
|
762 |
|
763 |
|
764 |
subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
765 |
real( kind = dp ) :: pot |
766 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
767 |
real (kind=dp), dimension(9,getNlocal()) :: A |
768 |
real (kind=dp), dimension(3,getNlocal()) :: f |
769 |
real (kind=dp), dimension(3,getNlocal()) :: t |
770 |
|
771 |
logical, intent(inout) :: do_pot, do_stress |
772 |
integer, intent(in) :: i, j |
773 |
real ( kind = dp ), intent(inout) :: rijsq |
774 |
real ( kind = dp ) :: r |
775 |
real ( kind = dp ), intent(inout) :: d(3) |
776 |
|
777 |
logical :: is_EAM_i, is_EAM_j |
778 |
|
779 |
integer :: me_i, me_j |
780 |
|
781 |
r = sqrt(rijsq) |
782 |
|
783 |
#ifdef IS_MPI |
784 |
if (tagRow(i) .eq. tagColumn(j)) then |
785 |
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
786 |
endif |
787 |
|
788 |
me_i = atid_row(i) |
789 |
me_j = atid_col(j) |
790 |
|
791 |
#else |
792 |
|
793 |
me_i = atid(i) |
794 |
me_j = atid(j) |
795 |
|
796 |
#endif |
797 |
|
798 |
if (FF_uses_EAM .and. SimUsesEAM()) then |
799 |
call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i) |
800 |
call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j) |
801 |
|
802 |
if ( is_EAM_i .and. is_EAM_j ) & |
803 |
call calc_EAM_prepair_rho(i, j, d, r, rijsq ) |
804 |
endif |
805 |
end subroutine do_prepair |
806 |
|
807 |
|
808 |
|
809 |
|
810 |
subroutine do_preforce(nlocal,pot) |
811 |
integer :: nlocal |
812 |
real( kind = dp ) :: pot |
813 |
|
814 |
if (FF_uses_EAM .and. SimUsesEAM()) then |
815 |
call calc_EAM_preforce_Frho(nlocal,pot) |
816 |
endif |
817 |
|
818 |
|
819 |
end subroutine do_preforce |
820 |
|
821 |
|
822 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
823 |
|
824 |
real (kind = dp), dimension(3) :: q_i |
825 |
real (kind = dp), dimension(3) :: q_j |
826 |
real ( kind = dp ), intent(out) :: r_sq |
827 |
real( kind = dp ) :: d(3), scaled(3) |
828 |
integer i |
829 |
|
830 |
d(1:3) = q_j(1:3) - q_i(1:3) |
831 |
|
832 |
! Wrap back into periodic box if necessary |
833 |
if ( SimUsesPBC() ) then |
834 |
|
835 |
if( .not.boxIsOrthorhombic ) then |
836 |
! calc the scaled coordinates. |
837 |
|
838 |
scaled = matmul(HmatInv, d) |
839 |
|
840 |
! wrap the scaled coordinates |
841 |
|
842 |
scaled = scaled - anint(scaled) |
843 |
|
844 |
|
845 |
! calc the wrapped real coordinates from the wrapped scaled |
846 |
! coordinates |
847 |
|
848 |
d = matmul(Hmat,scaled) |
849 |
|
850 |
else |
851 |
! calc the scaled coordinates. |
852 |
|
853 |
do i = 1, 3 |
854 |
scaled(i) = d(i) * HmatInv(i,i) |
855 |
|
856 |
! wrap the scaled coordinates |
857 |
|
858 |
scaled(i) = scaled(i) - anint(scaled(i)) |
859 |
|
860 |
! calc the wrapped real coordinates from the wrapped scaled |
861 |
! coordinates |
862 |
|
863 |
d(i) = scaled(i)*Hmat(i,i) |
864 |
enddo |
865 |
endif |
866 |
|
867 |
endif |
868 |
|
869 |
r_sq = dot_product(d,d) |
870 |
|
871 |
end subroutine get_interatomic_vector |
872 |
|
873 |
subroutine check_initialization(error) |
874 |
integer, intent(out) :: error |
875 |
|
876 |
error = 0 |
877 |
! Make sure we are properly initialized. |
878 |
if (.not. do_forces_initialized) then |
879 |
error = -1 |
880 |
return |
881 |
endif |
882 |
|
883 |
#ifdef IS_MPI |
884 |
if (.not. isMPISimSet()) then |
885 |
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
886 |
error = -1 |
887 |
return |
888 |
endif |
889 |
#endif |
890 |
|
891 |
return |
892 |
end subroutine check_initialization |
893 |
|
894 |
|
895 |
subroutine zero_work_arrays() |
896 |
|
897 |
#ifdef IS_MPI |
898 |
|
899 |
q_Row = 0.0_dp |
900 |
q_Col = 0.0_dp |
901 |
|
902 |
u_l_Row = 0.0_dp |
903 |
u_l_Col = 0.0_dp |
904 |
|
905 |
A_Row = 0.0_dp |
906 |
A_Col = 0.0_dp |
907 |
|
908 |
f_Row = 0.0_dp |
909 |
f_Col = 0.0_dp |
910 |
f_Temp = 0.0_dp |
911 |
|
912 |
t_Row = 0.0_dp |
913 |
t_Col = 0.0_dp |
914 |
t_Temp = 0.0_dp |
915 |
|
916 |
pot_Row = 0.0_dp |
917 |
pot_Col = 0.0_dp |
918 |
pot_Temp = 0.0_dp |
919 |
|
920 |
rf_Row = 0.0_dp |
921 |
rf_Col = 0.0_dp |
922 |
rf_Temp = 0.0_dp |
923 |
|
924 |
#endif |
925 |
|
926 |
rf = 0.0_dp |
927 |
tau_Temp = 0.0_dp |
928 |
virial_Temp = 0.0_dp |
929 |
end subroutine zero_work_arrays |
930 |
|
931 |
function skipThisPair(atom1, atom2) result(skip_it) |
932 |
integer, intent(in) :: atom1 |
933 |
integer, intent(in), optional :: atom2 |
934 |
logical :: skip_it |
935 |
integer :: unique_id_1, unique_id_2 |
936 |
integer :: me_i,me_j |
937 |
integer :: i |
938 |
|
939 |
skip_it = .false. |
940 |
|
941 |
!! there are a number of reasons to skip a pair or a particle |
942 |
!! mostly we do this to exclude atoms who are involved in short |
943 |
!! range interactions (bonds, bends, torsions), but we also need |
944 |
!! to exclude some overcounted interactions that result from |
945 |
!! the parallel decomposition |
946 |
|
947 |
#ifdef IS_MPI |
948 |
!! in MPI, we have to look up the unique IDs for each atom |
949 |
unique_id_1 = tagRow(atom1) |
950 |
#else |
951 |
!! in the normal loop, the atom numbers are unique |
952 |
unique_id_1 = atom1 |
953 |
#endif |
954 |
|
955 |
!! We were called with only one atom, so just check the global exclude |
956 |
!! list for this atom |
957 |
if (.not. present(atom2)) then |
958 |
do i = 1, nExcludes_global |
959 |
if (excludesGlobal(i) == unique_id_1) then |
960 |
skip_it = .true. |
961 |
return |
962 |
end if |
963 |
end do |
964 |
return |
965 |
end if |
966 |
|
967 |
#ifdef IS_MPI |
968 |
unique_id_2 = tagColumn(atom2) |
969 |
#else |
970 |
unique_id_2 = atom2 |
971 |
#endif |
972 |
|
973 |
#ifdef IS_MPI |
974 |
!! this situation should only arise in MPI simulations |
975 |
if (unique_id_1 == unique_id_2) then |
976 |
skip_it = .true. |
977 |
return |
978 |
end if |
979 |
|
980 |
!! this prevents us from doing the pair on multiple processors |
981 |
if (unique_id_1 < unique_id_2) then |
982 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
983 |
skip_it = .true. |
984 |
return |
985 |
endif |
986 |
else |
987 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
988 |
skip_it = .true. |
989 |
return |
990 |
endif |
991 |
endif |
992 |
#endif |
993 |
|
994 |
!! the rest of these situations can happen in all simulations: |
995 |
do i = 1, nExcludes_global |
996 |
if ((excludesGlobal(i) == unique_id_1) .or. & |
997 |
(excludesGlobal(i) == unique_id_2)) then |
998 |
skip_it = .true. |
999 |
return |
1000 |
endif |
1001 |
enddo |
1002 |
|
1003 |
do i = 1, nExcludes_local |
1004 |
if (excludesLocal(1,i) == unique_id_1) then |
1005 |
if (excludesLocal(2,i) == unique_id_2) then |
1006 |
skip_it = .true. |
1007 |
return |
1008 |
endif |
1009 |
else |
1010 |
if (excludesLocal(1,i) == unique_id_2) then |
1011 |
if (excludesLocal(2,i) == unique_id_1) then |
1012 |
skip_it = .true. |
1013 |
return |
1014 |
endif |
1015 |
endif |
1016 |
endif |
1017 |
end do |
1018 |
|
1019 |
return |
1020 |
end function skipThisPair |
1021 |
|
1022 |
function FF_UsesDirectionalAtoms() result(doesit) |
1023 |
logical :: doesit |
1024 |
doesit = FF_uses_dipoles .or. FF_uses_sticky .or. & |
1025 |
FF_uses_GB .or. FF_uses_RF |
1026 |
end function FF_UsesDirectionalAtoms |
1027 |
|
1028 |
function FF_RequiresPrepairCalc() result(doesit) |
1029 |
logical :: doesit |
1030 |
doesit = FF_uses_EAM |
1031 |
end function FF_RequiresPrepairCalc |
1032 |
|
1033 |
function FF_RequiresPostpairCalc() result(doesit) |
1034 |
logical :: doesit |
1035 |
doesit = FF_uses_RF |
1036 |
end function FF_RequiresPostpairCalc |
1037 |
|
1038 |
end module do_Forces |