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