ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1509
Committed: Tue Oct 5 19:46:37 2010 UTC (14 years, 7 months ago) by gezelter
File size: 42250 byte(s)
Log Message:
removing group-selective cutoff cruft

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date