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