1 |
!! |
2 |
!! Copyright (c) 2005, 2009 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
|
42 |
!! doForces.F90 |
43 |
!! module doForces |
44 |
!! Calculates Long Range forces. |
45 |
|
46 |
!! @author Charles F. Vardeman II |
47 |
!! @author Matthew Meineke |
48 |
!! @version $Id$, $Date$, $Name: not supported by cvs2svn $, $Revision$ |
49 |
|
50 |
|
51 |
module doForces |
52 |
use force_globals |
53 |
use simulation |
54 |
use definitions |
55 |
use neighborLists |
56 |
use vector_class |
57 |
use status |
58 |
use ISO_C_BINDING |
59 |
|
60 |
#ifdef IS_MPI |
61 |
use mpiSimulation |
62 |
#endif |
63 |
|
64 |
implicit none |
65 |
PRIVATE |
66 |
|
67 |
real(kind=dp), external :: get_cutoff |
68 |
|
69 |
#define __FORTRAN90 |
70 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
71 |
|
72 |
INTEGER, PARAMETER:: PREPAIR_LOOP = 1 |
73 |
INTEGER, PARAMETER:: PAIR_LOOP = 2 |
74 |
|
75 |
logical, save :: haveNeighborList = .false. |
76 |
logical, save :: haveSIMvariables = .false. |
77 |
logical, save :: haveCutoffs = .false. |
78 |
logical, save :: haveSkinThickness = .false. |
79 |
|
80 |
logical, save :: SIM_uses_DirectionalAtoms |
81 |
logical, save :: SIM_uses_MetallicAtoms |
82 |
logical, save :: SIM_requires_skip_correction |
83 |
logical, save :: SIM_requires_self_correction |
84 |
logical, save :: SIM_uses_PBC |
85 |
logical, save :: SIM_uses_AtomicVirial |
86 |
|
87 |
real(kind=dp), save :: rCut, rSwitch, rList, rListSq, rCutSq |
88 |
real(kind=dp), save :: skinThickness |
89 |
logical, save :: shifted_pot |
90 |
logical, save :: shifted_force |
91 |
|
92 |
public :: init_FF |
93 |
public :: setCutoffs |
94 |
public :: setSkinThickness |
95 |
public :: do_force_loop |
96 |
|
97 |
#ifdef PROFILE |
98 |
public :: getforcetime |
99 |
real, save :: forceTime = 0 |
100 |
real :: forceTimeInitial, forceTimeFinal |
101 |
integer :: nLoops |
102 |
#endif |
103 |
|
104 |
contains |
105 |
|
106 |
subroutine setCutoffs(rc, rsw, sp, sf) |
107 |
|
108 |
real(kind=dp),intent(in) :: rc, rsw |
109 |
integer, intent(in) :: sp, sf |
110 |
character(len = statusMsgSize) :: errMsg |
111 |
integer :: localError |
112 |
|
113 |
rCut = rc |
114 |
rSwitch = rsw |
115 |
rCutsq = rCut * rCut |
116 |
|
117 |
if (haveSkinThickness) then |
118 |
rList = rCut + skinThickness |
119 |
rListSq = rListSq |
120 |
endif |
121 |
|
122 |
if (sp .ne. 0) then |
123 |
shifted_pot = .true. |
124 |
else |
125 |
shifted_pot = .false. |
126 |
endif |
127 |
if (sf .ne. 0) then |
128 |
shifted_force = .true. |
129 |
else |
130 |
shifted_force = .false. |
131 |
endif |
132 |
|
133 |
if (abs(rCut - rSwitch) .lt. 0.0001) then |
134 |
if (shifted_force) then |
135 |
write(errMsg, *) & |
136 |
'cutoffRadius and switchingRadius are set to the', newline & |
137 |
// tab, 'same value. OpenMD will use shifted force', newline & |
138 |
// tab, 'potentials instead of switching functions.' |
139 |
|
140 |
call handleInfo("setCutoffs", errMsg) |
141 |
else |
142 |
write(errMsg, *) & |
143 |
'cutoffRadius and switchingRadius are set to the', newline & |
144 |
// tab, 'same value. OpenMD will use shifted', newline & |
145 |
// tab, 'potentials instead of switching functions.' |
146 |
|
147 |
call handleInfo("setCutoffs", errMsg) |
148 |
|
149 |
shifted_pot = .true. |
150 |
endif |
151 |
|
152 |
endif |
153 |
|
154 |
localError = 0 |
155 |
call set_switch(rSwitch, rCut) |
156 |
call setHmatDangerousRcutValue(rCut) |
157 |
|
158 |
haveCutoffs = .true. |
159 |
end subroutine setCutoffs |
160 |
|
161 |
subroutine setSkinThickness( thisSkin ) |
162 |
|
163 |
real(kind=dp), intent(in) :: thisSkin |
164 |
|
165 |
skinThickness = thisSkin |
166 |
haveSkinThickness = .true. |
167 |
|
168 |
if (haveCutoffs) then |
169 |
rList = rCut + skinThickness |
170 |
rListSq = rList * rList |
171 |
endif |
172 |
|
173 |
end subroutine setSkinThickness |
174 |
|
175 |
subroutine setSimVariables() |
176 |
SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms() |
177 |
SIM_uses_MetallicAtoms = SimUsesMetallicAtoms() |
178 |
SIM_requires_skip_correction = SimRequiresSkipCorrection() |
179 |
SIM_requires_self_correction = SimRequiresSelfCorrection() |
180 |
SIM_uses_PBC = SimUsesPBC() |
181 |
SIM_uses_AtomicVirial = SimUsesAtomicVirial() |
182 |
|
183 |
haveSIMvariables = .true. |
184 |
return |
185 |
end subroutine setSimVariables |
186 |
|
187 |
subroutine doReadyCheck(error) |
188 |
integer, intent(out) :: error |
189 |
integer :: myStatus |
190 |
|
191 |
error = 0 |
192 |
|
193 |
if (.not. haveSIMvariables) then |
194 |
call setSimVariables() |
195 |
endif |
196 |
|
197 |
if (.not. haveCutoffs) then |
198 |
write(default_error, *) 'cutoffs have not been set in doForces!' |
199 |
error = -1 |
200 |
return |
201 |
endif |
202 |
|
203 |
if (.not. haveSkinThickness) then |
204 |
write(default_error, *) 'skin thickness has not been set in doForces!' |
205 |
error = -1 |
206 |
return |
207 |
endif |
208 |
|
209 |
if (.not. haveNeighborList) then |
210 |
write(default_error, *) 'neighbor list has not been initialized in doForces!' |
211 |
error = -1 |
212 |
return |
213 |
end if |
214 |
|
215 |
|
216 |
#ifdef IS_MPI |
217 |
if (.not. isMPISimSet()) then |
218 |
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
219 |
error = -1 |
220 |
return |
221 |
endif |
222 |
#endif |
223 |
return |
224 |
end subroutine doReadyCheck |
225 |
|
226 |
subroutine init_FF(thisStat) |
227 |
|
228 |
integer :: my_status |
229 |
integer, intent(out) :: thisStat |
230 |
|
231 |
!! assume things are copacetic, unless they aren't |
232 |
thisStat = 0 |
233 |
|
234 |
if (.not. haveNeighborList) then |
235 |
!! Create neighbor lists |
236 |
call expandNeighborList(nLocal, my_status) |
237 |
if (my_Status /= 0) then |
238 |
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
239 |
thisStat = -1 |
240 |
return |
241 |
endif |
242 |
haveNeighborList = .true. |
243 |
endif |
244 |
|
245 |
end subroutine init_FF |
246 |
|
247 |
|
248 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
249 |
!-------------------------------------------------------------> |
250 |
subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, & |
251 |
error) |
252 |
!! Position array provided by C, dimensioned by getNlocal |
253 |
real ( kind = dp ), dimension(3, nLocal) :: q |
254 |
!! molecular center-of-mass position array |
255 |
real ( kind = dp ), dimension(3, nGroups) :: q_group |
256 |
!! Rotation Matrix for each long range particle in simulation. |
257 |
real( kind = dp), dimension(9, nLocal) :: A |
258 |
!! Unit vectors for dipoles (lab frame) |
259 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
260 |
!! Force array provided by C, dimensioned by getNlocal |
261 |
real ( kind = dp ), dimension(3,nLocal) :: f |
262 |
!! Torsion array provided by C, dimensioned by getNlocal |
263 |
real( kind = dp ), dimension(3,nLocal) :: t |
264 |
|
265 |
!! Stress Tensor |
266 |
real( kind = dp), dimension(9) :: tau |
267 |
real ( kind = dp ),dimension(LR_POT_TYPES) :: pot |
268 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
269 |
real( kind = dp ), dimension(nLocal) :: skipped_charge |
270 |
|
271 |
logical :: in_switching_region |
272 |
#ifdef IS_MPI |
273 |
real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local |
274 |
integer :: nAtomsInRow |
275 |
integer :: nAtomsInCol |
276 |
integer :: nprocs |
277 |
integer :: nGroupsInRow |
278 |
integer :: nGroupsInCol |
279 |
#endif |
280 |
integer :: natoms |
281 |
logical :: update_nlist |
282 |
integer :: i, j, jstart, jend, jnab |
283 |
integer :: istart, iend |
284 |
integer :: ia, jb, atom1, atom2 |
285 |
integer :: nlist |
286 |
real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij |
287 |
real( kind = DP ) :: sw, dswdr, swderiv, mf |
288 |
real( kind = DP ) :: rVal |
289 |
real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag |
290 |
real(kind=dp) :: rfpot, mu_i |
291 |
real(kind=dp):: rCut |
292 |
integer :: n_in_i, n_in_j, iG, j1 |
293 |
logical :: is_dp_i |
294 |
integer :: neighborListSize |
295 |
integer :: listerror, error |
296 |
integer :: localError |
297 |
integer :: propPack_i, propPack_j |
298 |
integer :: loopStart, loopEnd, loop |
299 |
integer :: i1, topoDist |
300 |
|
301 |
real(kind=dp) :: skch |
302 |
|
303 |
!! initialize local variables |
304 |
|
305 |
#ifdef IS_MPI |
306 |
pot_local = 0.0_dp |
307 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
308 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
309 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
310 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
311 |
#else |
312 |
natoms = nlocal |
313 |
#endif |
314 |
|
315 |
call doReadyCheck(localError) |
316 |
if ( localError .ne. 0 ) then |
317 |
call handleError("do_force_loop", "Not Initialized") |
318 |
error = -1 |
319 |
return |
320 |
end if |
321 |
call zero_work_arrays() |
322 |
|
323 |
! Gather all information needed by all force loops: |
324 |
|
325 |
#ifdef IS_MPI |
326 |
|
327 |
call gather(q, q_Row, plan_atom_row_3d) |
328 |
call gather(q, q_Col, plan_atom_col_3d) |
329 |
|
330 |
call gather(q_group, q_group_Row, plan_group_row_3d) |
331 |
call gather(q_group, q_group_Col, plan_group_col_3d) |
332 |
|
333 |
if (SIM_uses_DirectionalAtoms) then |
334 |
call gather(eFrame, eFrame_Row, plan_atom_row_rotation) |
335 |
call gather(eFrame, eFrame_Col, plan_atom_col_rotation) |
336 |
|
337 |
call gather(A, A_Row, plan_atom_row_rotation) |
338 |
call gather(A, A_Col, plan_atom_col_rotation) |
339 |
endif |
340 |
|
341 |
#endif |
342 |
|
343 |
!! Begin force loop timing: |
344 |
#ifdef PROFILE |
345 |
call cpu_time(forceTimeInitial) |
346 |
nloops = nloops + 1 |
347 |
#endif |
348 |
|
349 |
loopEnd = PAIR_LOOP |
350 |
if (SIM_uses_MetallicAtoms) then |
351 |
loopStart = PREPAIR_LOOP |
352 |
else |
353 |
loopStart = PAIR_LOOP |
354 |
endif |
355 |
|
356 |
do loop = loopStart, loopEnd |
357 |
|
358 |
! See if we need to update neighbor lists |
359 |
! (but only on the first time through): |
360 |
if (loop .eq. loopStart) then |
361 |
#ifdef IS_MPI |
362 |
call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, & |
363 |
update_nlist) |
364 |
#else |
365 |
call checkNeighborList(nGroups, q_group, skinThickness, & |
366 |
update_nlist) |
367 |
#endif |
368 |
endif |
369 |
|
370 |
if (update_nlist) then |
371 |
!! save current configuration and construct neighbor list |
372 |
#ifdef IS_MPI |
373 |
call saveNeighborList(nGroupsInRow, q_group_row) |
374 |
#else |
375 |
call saveNeighborList(nGroups, q_group) |
376 |
#endif |
377 |
neighborListSize = size(list) |
378 |
nlist = 0 |
379 |
endif |
380 |
|
381 |
istart = 1 |
382 |
#ifdef IS_MPI |
383 |
iend = nGroupsInRow |
384 |
#else |
385 |
iend = nGroups - 1 |
386 |
#endif |
387 |
outer: do i = istart, iend |
388 |
|
389 |
if (update_nlist) point(i) = nlist + 1 |
390 |
|
391 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
392 |
|
393 |
if (update_nlist) then |
394 |
#ifdef IS_MPI |
395 |
jstart = 1 |
396 |
jend = nGroupsInCol |
397 |
#else |
398 |
jstart = i+1 |
399 |
jend = nGroups |
400 |
#endif |
401 |
else |
402 |
jstart = point(i) |
403 |
jend = point(i+1) - 1 |
404 |
! make sure group i has neighbors |
405 |
if (jstart .gt. jend) cycle outer |
406 |
endif |
407 |
|
408 |
do jnab = jstart, jend |
409 |
if (update_nlist) then |
410 |
j = jnab |
411 |
else |
412 |
j = list(jnab) |
413 |
endif |
414 |
|
415 |
#ifdef IS_MPI |
416 |
call get_interatomic_vector(q_group_Row(:,i), & |
417 |
q_group_Col(:,j), d_grp, rgrpsq) |
418 |
#else |
419 |
call get_interatomic_vector(q_group(:,i), & |
420 |
q_group(:,j), d_grp, rgrpsq) |
421 |
#endif |
422 |
|
423 |
if (rgrpsq < rListsq) then |
424 |
if (update_nlist) then |
425 |
nlist = nlist + 1 |
426 |
|
427 |
if (nlist > neighborListSize) then |
428 |
#ifdef IS_MPI |
429 |
call expandNeighborList(nGroupsInRow, listerror) |
430 |
#else |
431 |
call expandNeighborList(nGroups, listerror) |
432 |
#endif |
433 |
if (listerror /= 0) then |
434 |
error = -1 |
435 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
436 |
return |
437 |
end if |
438 |
neighborListSize = size(list) |
439 |
endif |
440 |
|
441 |
list(nlist) = j |
442 |
endif |
443 |
|
444 |
if (rgrpsq < rCutsq) then |
445 |
|
446 |
if (loop .eq. PAIR_LOOP) then |
447 |
vij = 0.0_dp |
448 |
fij(1) = 0.0_dp |
449 |
fij(2) = 0.0_dp |
450 |
fij(3) = 0.0_dp |
451 |
endif |
452 |
|
453 |
call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region) |
454 |
|
455 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
456 |
|
457 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
458 |
|
459 |
atom1 = groupListRow(ia) |
460 |
|
461 |
inner: do jb = groupStartCol(j), groupStartCol(j+1)-1 |
462 |
|
463 |
atom2 = groupListCol(jb) |
464 |
|
465 |
if (skipThisPair(atom1, atom2)) cycle inner |
466 |
|
467 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
468 |
d_atm(1) = d_grp(1) |
469 |
d_atm(2) = d_grp(2) |
470 |
d_atm(3) = d_grp(3) |
471 |
ratmsq = rgrpsq |
472 |
else |
473 |
#ifdef IS_MPI |
474 |
call get_interatomic_vector(q_Row(:,atom1), & |
475 |
q_Col(:,atom2), d_atm, ratmsq) |
476 |
#else |
477 |
call get_interatomic_vector(q(:,atom1), & |
478 |
q(:,atom2), d_atm, ratmsq) |
479 |
#endif |
480 |
endif |
481 |
|
482 |
topoDist = getTopoDistance(atom1, atom2) |
483 |
|
484 |
if (loop .eq. PREPAIR_LOOP) then |
485 |
#ifdef IS_MPI |
486 |
call f_do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
487 |
rgrpsq, d_grp, rCut, & |
488 |
eFrame, A, f, t, pot_local) |
489 |
#else |
490 |
call f_do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
491 |
rgrpsq, d_grp, rCut, & |
492 |
eFrame, A, f, t, pot) |
493 |
#endif |
494 |
else |
495 |
#ifdef IS_MPI |
496 |
call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
497 |
eFrame, A, f, t, pot_local, particle_pot, vpair, & |
498 |
fpair, d_grp, rgrp, rCut, topoDist) |
499 |
! particle_pot will be accumulated from row & column |
500 |
! arrays later |
501 |
#else |
502 |
call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
503 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
504 |
fpair, d_grp, rgrp, rCut, topoDist) |
505 |
#endif |
506 |
vij = vij + vpair |
507 |
fij(1) = fij(1) + fpair(1) |
508 |
fij(2) = fij(2) + fpair(2) |
509 |
fij(3) = fij(3) + fpair(3) |
510 |
call add_stress_tensor(d_atm, fpair, tau) |
511 |
endif |
512 |
enddo inner |
513 |
enddo |
514 |
|
515 |
if (loop .eq. PAIR_LOOP) then |
516 |
if (in_switching_region) then |
517 |
swderiv = vij*dswdr/rgrp |
518 |
fg = swderiv*d_grp |
519 |
|
520 |
fij(1) = fij(1) + fg(1) |
521 |
fij(2) = fij(2) + fg(2) |
522 |
fij(3) = fij(3) + fg(3) |
523 |
|
524 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
525 |
call add_stress_tensor(d_atm, fg, tau) |
526 |
endif |
527 |
|
528 |
do ia=groupStartRow(i), groupStartRow(i+1)-1 |
529 |
atom1=groupListRow(ia) |
530 |
mf = mfactRow(atom1) |
531 |
! fg is the force on atom ia due to cutoff group's |
532 |
! presence in switching region |
533 |
fg = swderiv*d_grp*mf |
534 |
#ifdef IS_MPI |
535 |
f_Row(1,atom1) = f_Row(1,atom1) + fg(1) |
536 |
f_Row(2,atom1) = f_Row(2,atom1) + fg(2) |
537 |
f_Row(3,atom1) = f_Row(3,atom1) + fg(3) |
538 |
#else |
539 |
f(1,atom1) = f(1,atom1) + fg(1) |
540 |
f(2,atom1) = f(2,atom1) + fg(2) |
541 |
f(3,atom1) = f(3,atom1) + fg(3) |
542 |
#endif |
543 |
if (n_in_i .gt. 1) then |
544 |
if (SIM_uses_AtomicVirial) then |
545 |
! find the distance between the atom |
546 |
! and the center of the cutoff group: |
547 |
#ifdef IS_MPI |
548 |
call get_interatomic_vector(q_Row(:,atom1), & |
549 |
q_group_Row(:,i), dag, rag) |
550 |
#else |
551 |
call get_interatomic_vector(q(:,atom1), & |
552 |
q_group(:,i), dag, rag) |
553 |
#endif |
554 |
call add_stress_tensor(dag,fg,tau) |
555 |
endif |
556 |
endif |
557 |
enddo |
558 |
|
559 |
do jb=groupStartCol(j), groupStartCol(j+1)-1 |
560 |
atom2=groupListCol(jb) |
561 |
mf = mfactCol(atom2) |
562 |
! fg is the force on atom jb due to cutoff group's |
563 |
! presence in switching region |
564 |
fg = -swderiv*d_grp*mf |
565 |
#ifdef IS_MPI |
566 |
f_Col(1,atom2) = f_Col(1,atom2) + fg(1) |
567 |
f_Col(2,atom2) = f_Col(2,atom2) + fg(2) |
568 |
f_Col(3,atom2) = f_Col(3,atom2) + fg(3) |
569 |
#else |
570 |
f(1,atom2) = f(1,atom2) + fg(1) |
571 |
f(2,atom2) = f(2,atom2) + fg(2) |
572 |
f(3,atom2) = f(3,atom2) + fg(3) |
573 |
#endif |
574 |
if (n_in_j .gt. 1) then |
575 |
if (SIM_uses_AtomicVirial) then |
576 |
! find the distance between the atom |
577 |
! and the center of the cutoff group: |
578 |
#ifdef IS_MPI |
579 |
call get_interatomic_vector(q_Col(:,atom2), & |
580 |
q_group_Col(:,j), dag, rag) |
581 |
#else |
582 |
call get_interatomic_vector(q(:,atom2), & |
583 |
q_group(:,j), dag, rag) |
584 |
#endif |
585 |
call add_stress_tensor(dag,fg,tau) |
586 |
endif |
587 |
endif |
588 |
enddo |
589 |
endif |
590 |
!if (.not.SIM_uses_AtomicVirial) then |
591 |
! call add_stress_tensor(d_grp, fij, tau) |
592 |
!endif |
593 |
endif |
594 |
endif |
595 |
endif |
596 |
enddo |
597 |
|
598 |
enddo outer |
599 |
|
600 |
if (update_nlist) then |
601 |
#ifdef IS_MPI |
602 |
point(nGroupsInRow + 1) = nlist + 1 |
603 |
#else |
604 |
point(nGroups) = nlist + 1 |
605 |
#endif |
606 |
if (loop .eq. PREPAIR_LOOP) then |
607 |
! we just did the neighbor list update on the first |
608 |
! pass, so we don't need to do it |
609 |
! again on the second pass |
610 |
update_nlist = .false. |
611 |
endif |
612 |
endif |
613 |
|
614 |
if (loop .eq. PREPAIR_LOOP) then |
615 |
#ifdef IS_MPI |
616 |
call f_do_preforce(nlocal, pot_local, particle_pot) |
617 |
#else |
618 |
call f_do_preforce(nlocal, pot, particle_pot) |
619 |
#endif |
620 |
endif |
621 |
|
622 |
enddo |
623 |
|
624 |
!! Do timing |
625 |
#ifdef PROFILE |
626 |
call cpu_time(forceTimeFinal) |
627 |
forceTime = forceTime + forceTimeFinal - forceTimeInitial |
628 |
#endif |
629 |
|
630 |
#ifdef IS_MPI |
631 |
!!distribute forces |
632 |
|
633 |
f_temp = 0.0_dp |
634 |
call scatter(f_Row,f_temp,plan_atom_row_3d) |
635 |
do i = 1,nlocal |
636 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
637 |
end do |
638 |
|
639 |
f_temp = 0.0_dp |
640 |
call scatter(f_Col,f_temp,plan_atom_col_3d) |
641 |
do i = 1,nlocal |
642 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
643 |
end do |
644 |
|
645 |
if (SIM_uses_DirectionalAtoms) then |
646 |
t_temp = 0.0_dp |
647 |
call scatter(t_Row,t_temp,plan_atom_row_3d) |
648 |
do i = 1,nlocal |
649 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
650 |
end do |
651 |
t_temp = 0.0_dp |
652 |
call scatter(t_Col,t_temp,plan_atom_col_3d) |
653 |
|
654 |
do i = 1,nlocal |
655 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
656 |
end do |
657 |
endif |
658 |
|
659 |
! scatter/gather pot_row into the members of my column |
660 |
do i = 1,LR_POT_TYPES |
661 |
call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row) |
662 |
end do |
663 |
! scatter/gather pot_local into all other procs |
664 |
! add resultant to get total pot |
665 |
do i = 1, nlocal |
666 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) & |
667 |
+ pot_Temp(1:LR_POT_TYPES,i) |
668 |
enddo |
669 |
|
670 |
do i = 1,LR_POT_TYPES |
671 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
672 |
enddo |
673 |
|
674 |
pot_Temp = 0.0_DP |
675 |
|
676 |
do i = 1,LR_POT_TYPES |
677 |
call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col) |
678 |
end do |
679 |
|
680 |
do i = 1, nlocal |
681 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)& |
682 |
+ pot_Temp(1:LR_POT_TYPES,i) |
683 |
enddo |
684 |
|
685 |
do i = 1,LR_POT_TYPES |
686 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
687 |
enddo |
688 |
|
689 |
ppot_Temp = 0.0_DP |
690 |
|
691 |
call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row) |
692 |
do i = 1, nlocal |
693 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
694 |
enddo |
695 |
|
696 |
ppot_Temp = 0.0_DP |
697 |
|
698 |
call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col) |
699 |
do i = 1, nlocal |
700 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
701 |
enddo |
702 |
|
703 |
#endif |
704 |
|
705 |
if (SIM_requires_skip_correction) then |
706 |
do i = 1, nlocal |
707 |
|
708 |
do i1 = 1, nSkipsForLocalAtom(i) |
709 |
j = skipsForLocalAtom(i, i1) |
710 |
|
711 |
! prevent overcounting the skips |
712 |
if (i.lt.j) then |
713 |
|
714 |
call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq) |
715 |
rVal = sqrt(ratmsq) |
716 |
call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region) |
717 |
#ifdef IS_MPI |
718 |
call do_skip_correction(c_idents_local(i), c_idents_local(j), & |
719 |
d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, & |
720 |
1.0_dp, pot_local(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j)) |
721 |
# else |
722 |
call do_skip_correction(c_idents_local(i), c_idents_local(j), & |
723 |
d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, & |
724 |
1.0_dp, pot(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j)) |
725 |
#endif |
726 |
endif |
727 |
enddo |
728 |
enddo |
729 |
endif |
730 |
|
731 |
if (SIM_requires_self_correction) then |
732 |
|
733 |
do i = 1, nlocal |
734 |
|
735 |
#ifdef IS_MPI |
736 |
call do_self_correction(c_idents_local(i), eFrame(:,i), & |
737 |
skipped_charge(i), pot_local(ELECTROSTATIC_POT), t(:,i)) |
738 |
#else |
739 |
call do_self_correction(c_idents_local(i), eFrame(:,i), & |
740 |
skipped_charge(i), pot(ELECTROSTATIC_POT), t(:,i)) |
741 |
#endif |
742 |
enddo |
743 |
endif |
744 |
|
745 |
#ifdef IS_MPI |
746 |
#ifdef SINGLE_PRECISION |
747 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, & |
748 |
mpi_comm_world,mpi_err) |
749 |
#else |
750 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, & |
751 |
mpi_sum, mpi_comm_world,mpi_err) |
752 |
#endif |
753 |
#endif |
754 |
|
755 |
end subroutine do_force_loop |
756 |
|
757 |
subroutine f_do_pair(i, j, rijsq, d, sw, & |
758 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
759 |
fpair, d_grp, r_grp, rCut, topoDist) |
760 |
|
761 |
real( kind = dp ) :: vpair, sw |
762 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot, pairpot |
763 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
764 |
real( kind = dp ), dimension(3) :: fpair |
765 |
real( kind = dp ), dimension(nLocal) :: mfact |
766 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
767 |
real( kind = dp ), dimension(9,nLocal) :: A |
768 |
real( kind = dp ), dimension(3,nLocal) :: f |
769 |
real( kind = dp ), dimension(3,nLocal) :: t |
770 |
|
771 |
integer, intent(in) :: i, j |
772 |
real ( kind = dp ), intent(inout) :: rijsq |
773 |
real ( kind = dp ), intent(inout) :: r_grp |
774 |
real ( kind = dp ), intent(inout) :: d(3) |
775 |
real ( kind = dp ), intent(inout) :: d_grp(3) |
776 |
real ( kind = dp ), intent(inout) :: rCut |
777 |
integer, intent(inout) :: topoDist |
778 |
real ( kind = dp ) :: r, pair_pot |
779 |
real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx |
780 |
|
781 |
real( kind = dp), dimension(3) :: f1, t1, t2 |
782 |
real( kind = dp), dimension(9) :: A1, A2, eF1, eF2 |
783 |
real( kind = dp) :: dfrhodrho_i, dfrhodrho_j |
784 |
real( kind = dp) :: rho_i, rho_j |
785 |
real( kind = dp) :: fshift_i, fshift_j |
786 |
real( kind = dp) :: p_vdw, p_elect, p_hb, p_met |
787 |
integer :: id1, id2, idx |
788 |
integer :: k |
789 |
integer :: c_ident_i, c_ident_j |
790 |
|
791 |
integer :: iHash |
792 |
|
793 |
r = sqrt(rijsq) |
794 |
|
795 |
vpair = 0.0_dp |
796 |
fpair(1:3) = 0.0_dp |
797 |
|
798 |
p_vdw = 0.0 |
799 |
p_elect = 0.0 |
800 |
p_hb = 0.0 |
801 |
p_met = 0.0 |
802 |
|
803 |
f1(1:3) = 0.0 |
804 |
t1(1:3) = 0.0 |
805 |
t2(1:3) = 0.0 |
806 |
|
807 |
#ifdef IS_MPI |
808 |
c_ident_i = c_idents_row(i) |
809 |
c_ident_j = c_idents_col(j) |
810 |
|
811 |
A1(:) = A_Row(:,i) |
812 |
A2(:) = A_Col(:,j) |
813 |
eF1(:) = eFrame_Row(:,i) |
814 |
eF2(:) = eFrame_Col(:,j) |
815 |
|
816 |
dfrhodrho_i = dfrhodrho_row(i) |
817 |
dfrhodrho_j = dfrhodrho_col(j) |
818 |
rho_i = rho_row(i) |
819 |
rho_j = rho_col(j) |
820 |
#else |
821 |
c_ident_i = c_idents_local(i) |
822 |
c_ident_j = c_idents_local(j) |
823 |
|
824 |
A1(:) = A(:,i) |
825 |
A2(:) = A(:,j) |
826 |
eF1(:) = eFrame(:,i) |
827 |
eF2(:) = eFrame(:,j) |
828 |
|
829 |
dfrhodrho_i = dfrhodrho(i) |
830 |
dfrhodrho_j = dfrhodrho(j) |
831 |
rho_i = rho(i) |
832 |
rho_j = rho(j) |
833 |
#endif |
834 |
|
835 |
call do_pair(c_ident_i, c_ident_j, d, r, rijsq, sw, vpair, & |
836 |
topoDist, A1, A2, eF1, eF2, & |
837 |
pairpot, f1, t1, t2, & |
838 |
rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j) |
839 |
|
840 |
#ifdef IS_MPI |
841 |
id1 = AtomRowToGlobal(i) |
842 |
id2 = AtomColToGlobal(j) |
843 |
|
844 |
pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*pairpot(VDW_POT) |
845 |
pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*pairpot(VDW_POT) |
846 |
pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*pairpot(ELECTROSTATIC_POT) |
847 |
pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*pairpot(ELECTROSTATIC_POT) |
848 |
pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*pairpot(HB_POT) |
849 |
pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*pairpot(HB_POT) |
850 |
pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*pairpot(METALLIC_POT) |
851 |
pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*pairpot(METALLIC_POT) |
852 |
|
853 |
do idx = 1, 3 |
854 |
f_Row(idx,i) = f_Row(idx,i) + f1(idx) |
855 |
f_Col(idx,j) = f_Col(idx,j) - f1(idx) |
856 |
|
857 |
t_Row(idx,i) = t_Row(idx,i) + t1(idx) |
858 |
t_Col(idx,j) = t_Col(idx,j) + t2(idx) |
859 |
enddo |
860 |
! particle_pot is the difference between the full potential |
861 |
! and the full potential without the presence of a particular |
862 |
! particle (atom1). |
863 |
! |
864 |
! This reduces the density at other particle locations, so |
865 |
! we need to recompute the density at atom2 assuming atom1 |
866 |
! didn't contribute. This then requires recomputing the |
867 |
! density functional for atom2 as well. |
868 |
! |
869 |
! Most of the particle_pot heavy lifting comes from the |
870 |
! pair interaction, and will be handled by vpair. Parallel version. |
871 |
|
872 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
873 |
ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j |
874 |
ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i |
875 |
end if |
876 |
|
877 |
#else |
878 |
id1 = i |
879 |
id2 = j |
880 |
|
881 |
pot(VDW_POT) = pot(VDW_POT) + pairpot(VDW_POT) |
882 |
pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + pairpot(ELECTROSTATIC_POT) |
883 |
pot(HB_POT) = pot(HB_POT) + pairpot(HB_POT) |
884 |
pot(METALLIC_POT) = pot(METALLIC_POT) + pairpot(METALLIC_POT) |
885 |
|
886 |
do idx = 1, 3 |
887 |
f(idx,i) = f(idx,i) + f1(idx) |
888 |
f(idx,j) = f(idx,j) - f1(idx) |
889 |
|
890 |
t(idx,i) = t(idx,i) + t1(idx) |
891 |
t(idx,j) = t(idx,j) + t2(idx) |
892 |
enddo |
893 |
! particle_pot is the difference between the full potential |
894 |
! and the full potential without the presence of a particular |
895 |
! particle (atom1). |
896 |
! |
897 |
! This reduces the density at other particle locations, so |
898 |
! we need to recompute the density at atom2 assuming atom1 |
899 |
! didn't contribute. This then requires recomputing the |
900 |
! density functional for atom2 as well. |
901 |
! |
902 |
! Most of the particle_pot heavy lifting comes from the |
903 |
! pair interaction, and will be handled by vpair. NonParallel version. |
904 |
|
905 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
906 |
particle_pot(i) = particle_pot(i) - frho(j) + fshift_j |
907 |
particle_pot(j) = particle_pot(j) - frho(i) + fshift_i |
908 |
end if |
909 |
|
910 |
|
911 |
#endif |
912 |
|
913 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
914 |
|
915 |
fpair(1) = fpair(1) + f1(1) |
916 |
fpair(2) = fpair(2) + f1(2) |
917 |
fpair(3) = fpair(3) + f1(3) |
918 |
|
919 |
endif |
920 |
end subroutine f_do_pair |
921 |
|
922 |
subroutine f_do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, & |
923 |
eFrame, A, f, t, pot) |
924 |
|
925 |
real( kind = dp ) :: sw |
926 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
927 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
928 |
real (kind=dp), dimension(9,nLocal) :: A |
929 |
real (kind=dp), dimension(3,nLocal) :: f |
930 |
real (kind=dp), dimension(3,nLocal) :: t |
931 |
|
932 |
integer, intent(in) :: i, j |
933 |
real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut |
934 |
real ( kind = dp ) :: r, rc |
935 |
real ( kind = dp ), intent(inout) :: d(3), dc(3) |
936 |
real ( kind = dp ) :: rho_i_at_j, rho_j_at_i |
937 |
integer :: c_ident_i, c_ident_j |
938 |
|
939 |
r = sqrt(rijsq) |
940 |
|
941 |
#ifdef IS_MPI |
942 |
c_ident_i = c_idents_row(i) |
943 |
c_ident_j = c_idents_col(j) |
944 |
#else |
945 |
c_ident_i = c_idents_local(i) |
946 |
c_ident_j = c_idents_local(j) |
947 |
#endif |
948 |
rho_i_at_j = 0.0_dp |
949 |
rho_j_at_i = 0.0_dp |
950 |
|
951 |
call do_prepair(c_ident_i, c_ident_j, r, & |
952 |
rho_i_at_j, rho_j_at_i) |
953 |
|
954 |
#ifdef IS_MPI |
955 |
rho_col(j) = rho_col(j) + rho_i_at_j |
956 |
rho_row(i) = rho_row(i) + rho_j_at_i |
957 |
#else |
958 |
rho(j) = rho(j) + rho_i_at_j |
959 |
rho(i) = rho(i) + rho_j_at_i |
960 |
#endif |
961 |
|
962 |
end subroutine f_do_prepair |
963 |
|
964 |
|
965 |
subroutine f_do_preforce(nlocal, pot, particle_pot) |
966 |
integer :: nlocal |
967 |
real( kind = dp ),dimension(LR_POT_TYPES) :: pot |
968 |
real( kind = dp ),dimension(nlocal) :: particle_pot |
969 |
integer :: sc_err = 0 |
970 |
integer :: atom, c_ident1 |
971 |
|
972 |
if (SIM_uses_MetallicAtoms) then |
973 |
|
974 |
#ifdef IS_MPI |
975 |
call scatter(rho_row,rho,plan_atom_row,sc_err) |
976 |
if (sc_err /= 0 ) then |
977 |
call handleError("do_preforce()", "Error scattering rho_row into rho") |
978 |
endif |
979 |
call scatter(rho_col,rho_tmp,plan_atom_col,sc_err) |
980 |
if (sc_err /= 0 ) then |
981 |
call handleError("do_preforce()", "Error scattering rho_col into rho") |
982 |
endif |
983 |
rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) |
984 |
#endif |
985 |
|
986 |
|
987 |
do atom = 1, nlocal |
988 |
c_ident1 = c_idents_local(atom) |
989 |
|
990 |
call do_preforce(c_ident1, rho(atom), frho(atom), dfrhodrho(atom)) |
991 |
pot(METALLIC_POT) = pot(METALLIC_POT) + frho(atom) |
992 |
particle_pot(atom) = particle_pot(atom) + frho(atom) |
993 |
end do |
994 |
|
995 |
#ifdef IS_MPI |
996 |
!! communicate f(rho) and derivatives back into row and column arrays |
997 |
call gather(frho,frho_row,plan_atom_row, sc_err) |
998 |
if (sc_err /= 0) then |
999 |
call handleError("do_preforce()","MPI gather frho_row failure") |
1000 |
endif |
1001 |
call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err) |
1002 |
if (sc_err /= 0) then |
1003 |
call handleError("do_preforce()","MPI gather dfrhodrho_row failure") |
1004 |
endif |
1005 |
call gather(frho,frho_col,plan_atom_col, sc_err) |
1006 |
if (sc_err /= 0) then |
1007 |
call handleError("do_preforce()","MPI gather frho_col failure") |
1008 |
endif |
1009 |
call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err) |
1010 |
if (sc_err /= 0) then |
1011 |
call handleError("do_preforce()","MPI gather dfrhodrho_col failure") |
1012 |
endif |
1013 |
#endif |
1014 |
|
1015 |
end if |
1016 |
end subroutine f_do_preforce |
1017 |
|
1018 |
|
1019 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
1020 |
|
1021 |
real(kind = dp), dimension(3) :: q_i |
1022 |
real(kind = dp), dimension(3) :: q_j |
1023 |
real(kind = dp), intent(out) :: r_sq |
1024 |
real(kind = dp) :: d(3), scaled(3) |
1025 |
real(kind = dp) :: t |
1026 |
integer i |
1027 |
|
1028 |
d(1) = q_j(1) - q_i(1) |
1029 |
d(2) = q_j(2) - q_i(2) |
1030 |
d(3) = q_j(3) - q_i(3) |
1031 |
|
1032 |
! Wrap back into periodic box if necessary |
1033 |
if ( SIM_uses_PBC ) then |
1034 |
|
1035 |
if( .not.boxIsOrthorhombic ) then |
1036 |
! calc the scaled coordinates. |
1037 |
! scaled = matmul(HmatInv, d) |
1038 |
! unwrap the matmul and do things explicitly: |
1039 |
|
1040 |
scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3) |
1041 |
scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3) |
1042 |
scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3) |
1043 |
|
1044 |
! wrap the scaled coordinates (but don't use anint for speed) |
1045 |
|
1046 |
t = scaled(1) |
1047 |
if (t .ge. 0.0) then |
1048 |
scaled(1) = t - floor(t + 0.5) |
1049 |
else |
1050 |
scaled(1) = t + ceiling(t - 0.5) |
1051 |
endif |
1052 |
|
1053 |
t = scaled(2) |
1054 |
if (t .ge. 0.0) then |
1055 |
scaled(2) = t - floor(t + 0.5) |
1056 |
else |
1057 |
scaled(2) = t + ceiling(t - 0.5) |
1058 |
endif |
1059 |
|
1060 |
t = scaled(3) |
1061 |
if (t .ge. 0.0) then |
1062 |
scaled(3) = t - floor(t + 0.5) |
1063 |
else |
1064 |
scaled(3) = t + ceiling(t - 0.5) |
1065 |
endif |
1066 |
|
1067 |
! calc the wrapped real coordinates from the wrapped scaled |
1068 |
! coordinates |
1069 |
! d = matmul(Hmat,scaled) |
1070 |
! unwrap the matmul and do things explicitly: |
1071 |
|
1072 |
d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3) |
1073 |
d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3) |
1074 |
d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3) |
1075 |
|
1076 |
else |
1077 |
! calc the scaled coordinates. |
1078 |
|
1079 |
scaled(1) = d(1) * HmatInv(1,1) |
1080 |
scaled(2) = d(2) * HmatInv(2,2) |
1081 |
scaled(3) = d(3) * HmatInv(3,3) |
1082 |
|
1083 |
! wrap the scaled coordinates (but don't use anint for speed) |
1084 |
|
1085 |
t = scaled(1) |
1086 |
if (t .ge. 0.0) then |
1087 |
scaled(1) = t - floor(t + 0.5) |
1088 |
else |
1089 |
scaled(1) = t + ceiling(t - 0.5) |
1090 |
endif |
1091 |
|
1092 |
t = scaled(2) |
1093 |
if (t .ge. 0.0) then |
1094 |
scaled(2) = t - floor(t + 0.5) |
1095 |
else |
1096 |
scaled(2) = t + ceiling(t - 0.5) |
1097 |
endif |
1098 |
|
1099 |
t = scaled(3) |
1100 |
if (t .ge. 0.0) then |
1101 |
scaled(3) = t - floor(t + 0.5) |
1102 |
else |
1103 |
scaled(3) = t + ceiling(t - 0.5) |
1104 |
endif |
1105 |
|
1106 |
! calc the wrapped real coordinates from the wrapped scaled |
1107 |
! coordinates |
1108 |
|
1109 |
d(1) = scaled(1)*Hmat(1,1) |
1110 |
d(2) = scaled(2)*Hmat(2,2) |
1111 |
d(3) = scaled(3)*Hmat(3,3) |
1112 |
|
1113 |
endif |
1114 |
|
1115 |
endif |
1116 |
|
1117 |
r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3) |
1118 |
|
1119 |
end subroutine get_interatomic_vector |
1120 |
|
1121 |
subroutine zero_work_arrays() |
1122 |
|
1123 |
#ifdef IS_MPI |
1124 |
|
1125 |
q_Row = 0.0_dp |
1126 |
q_Col = 0.0_dp |
1127 |
|
1128 |
q_group_Row = 0.0_dp |
1129 |
q_group_Col = 0.0_dp |
1130 |
|
1131 |
eFrame_Row = 0.0_dp |
1132 |
eFrame_Col = 0.0_dp |
1133 |
|
1134 |
A_Row = 0.0_dp |
1135 |
A_Col = 0.0_dp |
1136 |
|
1137 |
f_Row = 0.0_dp |
1138 |
f_Col = 0.0_dp |
1139 |
f_Temp = 0.0_dp |
1140 |
|
1141 |
t_Row = 0.0_dp |
1142 |
t_Col = 0.0_dp |
1143 |
t_Temp = 0.0_dp |
1144 |
|
1145 |
pot_Row = 0.0_dp |
1146 |
pot_Col = 0.0_dp |
1147 |
pot_Temp = 0.0_dp |
1148 |
ppot_Temp = 0.0_dp |
1149 |
|
1150 |
frho_row = 0.0_dp |
1151 |
frho_col = 0.0_dp |
1152 |
rho_row = 0.0_dp |
1153 |
rho_col = 0.0_dp |
1154 |
rho_tmp = 0.0_dp |
1155 |
dfrhodrho_row = 0.0_dp |
1156 |
dfrhodrho_col = 0.0_dp |
1157 |
|
1158 |
#endif |
1159 |
rho = 0.0_dp |
1160 |
frho = 0.0_dp |
1161 |
dfrhodrho = 0.0_dp |
1162 |
|
1163 |
end subroutine zero_work_arrays |
1164 |
|
1165 |
function skipThisPair(atom1, atom2) result(skip_it) |
1166 |
integer, intent(in) :: atom1 |
1167 |
integer, intent(in), optional :: atom2 |
1168 |
logical :: skip_it |
1169 |
integer :: unique_id_1, unique_id_2 |
1170 |
integer :: i |
1171 |
|
1172 |
skip_it = .false. |
1173 |
|
1174 |
!! there are a number of reasons to skip a pair or a particle |
1175 |
!! mostly we do this to exclude atoms who are involved in short |
1176 |
!! range interactions (bonds, bends, torsions), but we also need |
1177 |
!! to exclude some overcounted interactions that result from |
1178 |
!! the parallel decomposition |
1179 |
|
1180 |
#ifdef IS_MPI |
1181 |
!! in MPI, we have to look up the unique IDs for each atom |
1182 |
unique_id_1 = AtomRowToGlobal(atom1) |
1183 |
unique_id_2 = AtomColToGlobal(atom2) |
1184 |
!! this situation should only arise in MPI simulations |
1185 |
if (unique_id_1 == unique_id_2) then |
1186 |
skip_it = .true. |
1187 |
return |
1188 |
end if |
1189 |
|
1190 |
!! this prevents us from doing the pair on multiple processors |
1191 |
if (unique_id_1 < unique_id_2) then |
1192 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
1193 |
skip_it = .true. |
1194 |
return |
1195 |
endif |
1196 |
else |
1197 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
1198 |
skip_it = .true. |
1199 |
return |
1200 |
endif |
1201 |
endif |
1202 |
#else |
1203 |
!! in the normal loop, the atom numbers are unique |
1204 |
unique_id_1 = atom1 |
1205 |
unique_id_2 = atom2 |
1206 |
#endif |
1207 |
|
1208 |
#ifdef IS_MPI |
1209 |
do i = 1, nSkipsForRowAtom(atom1) |
1210 |
if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then |
1211 |
skip_it = .true. |
1212 |
return |
1213 |
endif |
1214 |
end do |
1215 |
#else |
1216 |
do i = 1, nSkipsForLocalAtom(atom1) |
1217 |
if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then |
1218 |
skip_it = .true. |
1219 |
return |
1220 |
endif |
1221 |
end do |
1222 |
#endif |
1223 |
|
1224 |
return |
1225 |
end function skipThisPair |
1226 |
|
1227 |
function getTopoDistance(atom1, atom2) result(topoDist) |
1228 |
integer, intent(in) :: atom1 |
1229 |
integer, intent(in) :: atom2 |
1230 |
integer :: topoDist |
1231 |
integer :: unique_id_2 |
1232 |
integer :: i |
1233 |
|
1234 |
#ifdef IS_MPI |
1235 |
unique_id_2 = AtomColToGlobal(atom2) |
1236 |
#else |
1237 |
unique_id_2 = atom2 |
1238 |
#endif |
1239 |
|
1240 |
! zero is default for unconnected (i.e. normal) pair interactions |
1241 |
|
1242 |
topoDist = 0 |
1243 |
|
1244 |
do i = 1, nTopoPairsForAtom(atom1) |
1245 |
if (toposForAtom(atom1, i) .eq. unique_id_2) then |
1246 |
topoDist = topoDistance(atom1, i) |
1247 |
return |
1248 |
endif |
1249 |
end do |
1250 |
|
1251 |
return |
1252 |
end function getTopoDistance |
1253 |
|
1254 |
#ifdef PROFILE |
1255 |
function getforcetime() result(totalforcetime) |
1256 |
real(kind=dp) :: totalforcetime |
1257 |
totalforcetime = forcetime |
1258 |
end function getforcetime |
1259 |
#endif |
1260 |
|
1261 |
!! This cleans componets of force arrays belonging only to fortran |
1262 |
|
1263 |
subroutine add_stress_tensor(dpair, fpair, tau) |
1264 |
|
1265 |
real( kind = dp ), dimension(3), intent(in) :: dpair, fpair |
1266 |
real( kind = dp ), dimension(9), intent(inout) :: tau |
1267 |
|
1268 |
! because the d vector is the rj - ri vector, and |
1269 |
! because fx, fy, fz are the force on atom i, we need a |
1270 |
! negative sign here: |
1271 |
|
1272 |
tau(1) = tau(1) - dpair(1) * fpair(1) |
1273 |
tau(2) = tau(2) - dpair(1) * fpair(2) |
1274 |
tau(3) = tau(3) - dpair(1) * fpair(3) |
1275 |
tau(4) = tau(4) - dpair(2) * fpair(1) |
1276 |
tau(5) = tau(5) - dpair(2) * fpair(2) |
1277 |
tau(6) = tau(6) - dpair(2) * fpair(3) |
1278 |
tau(7) = tau(7) - dpair(3) * fpair(1) |
1279 |
tau(8) = tau(8) - dpair(3) * fpair(2) |
1280 |
tau(9) = tau(9) - dpair(3) * fpair(3) |
1281 |
|
1282 |
end subroutine add_stress_tensor |
1283 |
|
1284 |
end module doForces |