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