ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/doForces.F90
Revision: 2530
Committed: Fri Dec 30 00:18:28 2005 UTC (18 years, 6 months ago) by chuckv
File size: 49665 byte(s)
Log Message:
Sutton-Chen bug fixes. Almost there...

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 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. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
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: doForces.F90,v 1.72 2005-12-30 00:18:28 chuckv Exp $, $Date: 2005-12-30 00:18:28 $, $Name: not supported by cvs2svn $, $Revision: 1.72 $
49
50
51 module doForces
52 use force_globals
53 use simulation
54 use definitions
55 use atype_module
56 use switcheroo
57 use neighborLists
58 use lj
59 use sticky
60 use electrostatic_module
61 use gayberne
62 use shapes
63 use vector_class
64 use eam
65 use suttonchen
66 use status
67 #ifdef IS_MPI
68 use mpiSimulation
69 #endif
70
71 implicit none
72 PRIVATE
73
74 #define __FORTRAN90
75 #include "UseTheForce/fSwitchingFunction.h"
76 #include "UseTheForce/fCutoffPolicy.h"
77 #include "UseTheForce/DarkSide/fInteractionMap.h"
78 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
79
80
81 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
82 INTEGER, PARAMETER:: PAIR_LOOP = 2
83
84 logical, save :: haveNeighborList = .false.
85 logical, save :: haveSIMvariables = .false.
86 logical, save :: haveSaneForceField = .false.
87 logical, save :: haveInteractionHash = .false.
88 logical, save :: haveGtypeCutoffMap = .false.
89 logical, save :: haveDefaultCutoffs = .false.
90 logical, save :: haveSkinThickness = .false.
91 logical, save :: haveElectrostaticSummationMethod = .false.
92 logical, save :: haveCutoffPolicy = .false.
93 logical, save :: VisitCutoffsAfterComputing = .false.
94
95 logical, save :: FF_uses_DirectionalAtoms
96 logical, save :: FF_uses_Dipoles
97 logical, save :: FF_uses_GayBerne
98 logical, save :: FF_uses_EAM
99 logical, save :: FF_uses_SC
100 logical, save :: FF_uses_MEAM
101
102
103 logical, save :: SIM_uses_DirectionalAtoms
104 logical, save :: SIM_uses_EAM
105 logical, save :: SIM_uses_SC
106 logical, save :: SIM_uses_MEAM
107 logical, save :: SIM_requires_postpair_calc
108 logical, save :: SIM_requires_prepair_calc
109 logical, save :: SIM_uses_PBC
110
111 integer, save :: electrostaticSummationMethod
112 integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY
113
114 real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut
115 real(kind=dp), save :: skinThickness
116 logical, save :: defaultDoShift
117
118 public :: init_FF
119 public :: setCutoffs
120 public :: cWasLame
121 public :: setElectrostaticMethod
122 public :: setCutoffPolicy
123 public :: setSkinThickness
124 public :: do_force_loop
125
126 #ifdef PROFILE
127 public :: getforcetime
128 real, save :: forceTime = 0
129 real :: forceTimeInitial, forceTimeFinal
130 integer :: nLoops
131 #endif
132
133 !! Variables for cutoff mapping and interaction mapping
134 ! Bit hash to determine pair-pair interactions.
135 integer, dimension(:,:), allocatable :: InteractionHash
136 real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff
137 real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow
138 real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol
139
140 integer, dimension(:), allocatable, target :: groupToGtypeRow
141 integer, dimension(:), pointer :: groupToGtypeCol => null()
142
143 real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow
144 real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol
145 type ::gtypeCutoffs
146 real(kind=dp) :: rcut
147 real(kind=dp) :: rcutsq
148 real(kind=dp) :: rlistsq
149 end type gtypeCutoffs
150 type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap
151
152 contains
153
154 subroutine createInteractionHash()
155 integer :: nAtypes
156 integer :: i
157 integer :: j
158 integer :: iHash
159 !! Test Types
160 logical :: i_is_LJ
161 logical :: i_is_Elect
162 logical :: i_is_Sticky
163 logical :: i_is_StickyP
164 logical :: i_is_GB
165 logical :: i_is_EAM
166 logical :: i_is_Shape
167 logical :: i_is_SC
168 logical :: i_is_MEAM
169 logical :: j_is_LJ
170 logical :: j_is_Elect
171 logical :: j_is_Sticky
172 logical :: j_is_StickyP
173 logical :: j_is_GB
174 logical :: j_is_EAM
175 logical :: j_is_Shape
176 logical :: j_is_SC
177 logical :: j_is_MEAM
178 real(kind=dp) :: myRcut
179
180 if (.not. associated(atypes)) then
181 call handleError("doForces", "atypes was not present before call of createInteractionHash!")
182 return
183 endif
184
185 nAtypes = getSize(atypes)
186
187 if (nAtypes == 0) then
188 call handleError("doForces", "nAtypes was zero during call of createInteractionHash!")
189 return
190 end if
191
192 if (.not. allocated(InteractionHash)) then
193 allocate(InteractionHash(nAtypes,nAtypes))
194 else
195 deallocate(InteractionHash)
196 allocate(InteractionHash(nAtypes,nAtypes))
197 endif
198
199 if (.not. allocated(atypeMaxCutoff)) then
200 allocate(atypeMaxCutoff(nAtypes))
201 else
202 deallocate(atypeMaxCutoff)
203 allocate(atypeMaxCutoff(nAtypes))
204 endif
205
206 do i = 1, nAtypes
207 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
208 call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
209 call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
210 call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
211 call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
212 call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
213 call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
214 call getElementProperty(atypes, i, "is_SC", i_is_SC)
215 call getElementProperty(atypes, i, "is_MEAM", i_is_MEAM)
216
217 do j = i, nAtypes
218
219 iHash = 0
220 myRcut = 0.0_dp
221
222 call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ)
223 call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect)
224 call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky)
225 call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP)
226 call getElementProperty(atypes, j, "is_GayBerne", j_is_GB)
227 call getElementProperty(atypes, j, "is_EAM", j_is_EAM)
228 call getElementProperty(atypes, j, "is_Shape", j_is_Shape)
229 call getElementProperty(atypes, j, "is_SC", j_is_SC)
230 call getElementProperty(atypes, j, "is_MEAM", j_is_MEAM)
231
232 if (i_is_LJ .and. j_is_LJ) then
233 iHash = ior(iHash, LJ_PAIR)
234 endif
235
236 if (i_is_Elect .and. j_is_Elect) then
237 iHash = ior(iHash, ELECTROSTATIC_PAIR)
238 endif
239
240 if (i_is_Sticky .and. j_is_Sticky) then
241 iHash = ior(iHash, STICKY_PAIR)
242 endif
243
244 if (i_is_StickyP .and. j_is_StickyP) then
245 iHash = ior(iHash, STICKYPOWER_PAIR)
246 endif
247
248 if (i_is_EAM .and. j_is_EAM) then
249 iHash = ior(iHash, EAM_PAIR)
250 endif
251
252 if (i_is_SC .and. j_is_SC) then
253 iHash = ior(iHash, SC_PAIR)
254 endif
255
256 if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR)
257 if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ)
258 if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ)
259
260 if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR)
261 if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ)
262 if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ)
263
264
265 InteractionHash(i,j) = iHash
266 InteractionHash(j,i) = iHash
267
268 end do
269
270 end do
271
272 haveInteractionHash = .true.
273 end subroutine createInteractionHash
274
275 subroutine createGtypeCutoffMap()
276
277 logical :: i_is_LJ
278 logical :: i_is_Elect
279 logical :: i_is_Sticky
280 logical :: i_is_StickyP
281 logical :: i_is_GB
282 logical :: i_is_EAM
283 logical :: i_is_Shape
284 logical :: i_is_SC
285 logical :: GtypeFound
286
287 integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend
288 integer :: n_in_i, me_i, ia, g, atom1, ja, n_in_j,me_j
289 integer :: nGroupsInRow
290 integer :: nGroupsInCol
291 integer :: nGroupTypesRow,nGroupTypesCol
292 real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol
293 real(kind=dp) :: biggestAtypeCutoff
294
295 if (.not. haveInteractionHash) then
296 call createInteractionHash()
297 endif
298 #ifdef IS_MPI
299 nGroupsInRow = getNgroupsInRow(plan_group_row)
300 nGroupsInCol = getNgroupsInCol(plan_group_col)
301 #endif
302 nAtypes = getSize(atypes)
303 ! Set all of the initial cutoffs to zero.
304 atypeMaxCutoff = 0.0_dp
305 do i = 1, nAtypes
306 if (SimHasAtype(i)) then
307 call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ)
308 call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect)
309 call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky)
310 call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP)
311 call getElementProperty(atypes, i, "is_GayBerne", i_is_GB)
312 call getElementProperty(atypes, i, "is_EAM", i_is_EAM)
313 call getElementProperty(atypes, i, "is_Shape", i_is_Shape)
314 call getElementProperty(atypes, i, "is_SC", i_is_SC)
315
316 if (haveDefaultCutoffs) then
317 atypeMaxCutoff(i) = defaultRcut
318 else
319 if (i_is_LJ) then
320 thisRcut = getSigma(i) * 2.5_dp
321 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
322 endif
323 if (i_is_Elect) then
324 thisRcut = defaultRcut
325 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
326 endif
327 if (i_is_Sticky) then
328 thisRcut = getStickyCut(i)
329 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
330 endif
331 if (i_is_StickyP) then
332 thisRcut = getStickyPowerCut(i)
333 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
334 endif
335 if (i_is_GB) then
336 thisRcut = getGayBerneCut(i)
337 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
338 endif
339 if (i_is_EAM) then
340 thisRcut = getEAMCut(i)
341 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
342 endif
343 if (i_is_Shape) then
344 thisRcut = getShapeCut(i)
345 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
346 endif
347 if (i_is_SC) then
348 thisRcut = getSCCut(i)
349 if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut
350 endif
351 endif
352
353 if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then
354 biggestAtypeCutoff = atypeMaxCutoff(i)
355 endif
356
357 endif
358 enddo
359
360 istart = 1
361 jstart = 1
362 #ifdef IS_MPI
363 iend = nGroupsInRow
364 jend = nGroupsInCol
365 #else
366 iend = nGroups
367 jend = nGroups
368 #endif
369
370 !! allocate the groupToGtype and gtypeMaxCutoff here.
371 if(.not.allocated(groupToGtypeRow)) then
372 ! allocate(groupToGtype(iend))
373 allocate(groupToGtypeRow(iend))
374 else
375 deallocate(groupToGtypeRow)
376 allocate(groupToGtypeRow(iend))
377 endif
378 if(.not.allocated(groupMaxCutoffRow)) then
379 allocate(groupMaxCutoffRow(iend))
380 else
381 deallocate(groupMaxCutoffRow)
382 allocate(groupMaxCutoffRow(iend))
383 end if
384
385 if(.not.allocated(gtypeMaxCutoffRow)) then
386 allocate(gtypeMaxCutoffRow(iend))
387 else
388 deallocate(gtypeMaxCutoffRow)
389 allocate(gtypeMaxCutoffRow(iend))
390 endif
391
392
393 #ifdef IS_MPI
394 ! We only allocate new storage if we are in MPI because Ncol /= Nrow
395 if(.not.associated(groupToGtypeCol)) then
396 allocate(groupToGtypeCol(jend))
397 else
398 deallocate(groupToGtypeCol)
399 allocate(groupToGtypeCol(jend))
400 end if
401
402 if(.not.associated(groupToGtypeCol)) then
403 allocate(groupToGtypeCol(jend))
404 else
405 deallocate(groupToGtypeCol)
406 allocate(groupToGtypeCol(jend))
407 end if
408 if(.not.associated(gtypeMaxCutoffCol)) then
409 allocate(gtypeMaxCutoffCol(jend))
410 else
411 deallocate(gtypeMaxCutoffCol)
412 allocate(gtypeMaxCutoffCol(jend))
413 end if
414
415 groupMaxCutoffCol = 0.0_dp
416 gtypeMaxCutoffCol = 0.0_dp
417
418 #endif
419 groupMaxCutoffRow = 0.0_dp
420 gtypeMaxCutoffRow = 0.0_dp
421
422
423 !! first we do a single loop over the cutoff groups to find the
424 !! largest cutoff for any atypes present in this group. We also
425 !! create gtypes at this point.
426
427 tol = 1.0d-6
428 nGroupTypesRow = 0
429
430 do i = istart, iend
431 n_in_i = groupStartRow(i+1) - groupStartRow(i)
432 groupMaxCutoffRow(i) = 0.0_dp
433 do ia = groupStartRow(i), groupStartRow(i+1)-1
434 atom1 = groupListRow(ia)
435 #ifdef IS_MPI
436 me_i = atid_row(atom1)
437 #else
438 me_i = atid(atom1)
439 #endif
440 if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then
441 groupMaxCutoffRow(i)=atypeMaxCutoff(me_i)
442 endif
443 enddo
444 if (nGroupTypesRow.eq.0) then
445 nGroupTypesRow = nGroupTypesRow + 1
446 gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
447 groupToGtypeRow(i) = nGroupTypesRow
448 else
449 GtypeFound = .false.
450 do g = 1, nGroupTypesRow
451 if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then
452 groupToGtypeRow(i) = g
453 GtypeFound = .true.
454 endif
455 enddo
456 if (.not.GtypeFound) then
457 nGroupTypesRow = nGroupTypesRow + 1
458 gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i)
459 groupToGtypeRow(i) = nGroupTypesRow
460 endif
461 endif
462 enddo
463
464 #ifdef IS_MPI
465 do j = jstart, jend
466 n_in_j = groupStartCol(j+1) - groupStartCol(j)
467 groupMaxCutoffCol(j) = 0.0_dp
468 do ja = groupStartCol(j), groupStartCol(j+1)-1
469 atom1 = groupListCol(ja)
470
471 me_j = atid_col(atom1)
472
473 if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then
474 groupMaxCutoffCol(j)=atypeMaxCutoff(me_j)
475 endif
476 enddo
477
478 if (nGroupTypesCol.eq.0) then
479 nGroupTypesCol = nGroupTypesCol + 1
480 gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
481 groupToGtypeCol(j) = nGroupTypesCol
482 else
483 GtypeFound = .false.
484 do g = 1, nGroupTypesCol
485 if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then
486 groupToGtypeCol(j) = g
487 GtypeFound = .true.
488 endif
489 enddo
490 if (.not.GtypeFound) then
491 nGroupTypesCol = nGroupTypesCol + 1
492 gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j)
493 groupToGtypeCol(j) = nGroupTypesCol
494 endif
495 endif
496 enddo
497
498 #else
499 ! Set pointers to information we just found
500 nGroupTypesCol = nGroupTypesRow
501 groupToGtypeCol => groupToGtypeRow
502 gtypeMaxCutoffCol => gtypeMaxCutoffRow
503 groupMaxCutoffCol => groupMaxCutoffRow
504 #endif
505
506 !! allocate the gtypeCutoffMap here.
507 allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol))
508 !! then we do a double loop over all the group TYPES to find the cutoff
509 !! map between groups of two types
510 tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol))
511
512 do i = 1, nGroupTypesRow
513 do j = 1, nGroupTypesCol
514
515 select case(cutoffPolicy)
516 case(TRADITIONAL_CUTOFF_POLICY)
517 thisRcut = tradRcut
518 case(MIX_CUTOFF_POLICY)
519 thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j))
520 case(MAX_CUTOFF_POLICY)
521 thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j))
522 case default
523 call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy")
524 return
525 end select
526 gtypeCutoffMap(i,j)%rcut = thisRcut
527
528 if (thisRcut.gt.largestRcut) largestRcut = thisRcut
529
530 gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut
531
532 if (.not.haveSkinThickness) then
533 skinThickness = 1.0_dp
534 endif
535
536 gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2
537
538 ! sanity check
539
540 if (haveDefaultCutoffs) then
541 if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then
542 call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff")
543 endif
544 endif
545 enddo
546 enddo
547
548 if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow)
549 if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow)
550 if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff)
551 #ifdef IS_MPI
552 if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol)
553 if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol)
554 #endif
555 groupMaxCutoffCol => null()
556 gtypeMaxCutoffCol => null()
557
558 haveGtypeCutoffMap = .true.
559 end subroutine createGtypeCutoffMap
560
561 subroutine setCutoffs(defRcut, defRsw)
562
563 real(kind=dp),intent(in) :: defRcut, defRsw
564 character(len = statusMsgSize) :: errMsg
565 integer :: localError
566
567 defaultRcut = defRcut
568 defaultRsw = defRsw
569
570 defaultDoShift = .false.
571 if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then
572
573 write(errMsg, *) &
574 'cutoffRadius and switchingRadius are set to the same', newline &
575 // tab, 'value. OOPSE will use shifted ', newline &
576 // tab, 'potentials instead of switching functions.'
577
578 call handleInfo("setCutoffs", errMsg)
579
580 defaultDoShift = .true.
581
582 endif
583
584 localError = 0
585 call setLJDefaultCutoff( defaultRcut, defaultDoShift )
586 call setElectrostaticCutoffRadius( defaultRcut, defaultRsw )
587 call setCutoffEAM( defaultRcut, localError)
588 if (localError /= 0) then
589 write(errMsg, *) 'An error has occured in setting the EAM cutoff'
590 call handleError("setCutoffs", errMsg)
591 end if
592 call set_switch(GROUP_SWITCH, defaultRsw, defaultRcut)
593
594 haveDefaultCutoffs = .true.
595 haveGtypeCutoffMap = .false.
596 end subroutine setCutoffs
597
598 subroutine cWasLame()
599
600 VisitCutoffsAfterComputing = .true.
601 return
602
603 end subroutine cWasLame
604
605 subroutine setCutoffPolicy(cutPolicy)
606
607 integer, intent(in) :: cutPolicy
608
609 cutoffPolicy = cutPolicy
610 haveCutoffPolicy = .true.
611 haveGtypeCutoffMap = .false.
612
613 end subroutine setCutoffPolicy
614
615 subroutine setElectrostaticMethod( thisESM )
616
617 integer, intent(in) :: thisESM
618
619 electrostaticSummationMethod = thisESM
620 haveElectrostaticSummationMethod = .true.
621
622 end subroutine setElectrostaticMethod
623
624 subroutine setSkinThickness( thisSkin )
625
626 real(kind=dp), intent(in) :: thisSkin
627
628 skinThickness = thisSkin
629 haveSkinThickness = .true.
630 haveGtypeCutoffMap = .false.
631
632 end subroutine setSkinThickness
633
634 subroutine setSimVariables()
635 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
636 SIM_uses_EAM = SimUsesEAM()
637 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
638 SIM_requires_prepair_calc = SimRequiresPrepairCalc()
639 SIM_uses_PBC = SimUsesPBC()
640
641 haveSIMvariables = .true.
642
643 return
644 end subroutine setSimVariables
645
646 subroutine doReadyCheck(error)
647 integer, intent(out) :: error
648
649 integer :: myStatus
650
651 error = 0
652
653 if (.not. haveInteractionHash) then
654 call createInteractionHash()
655 endif
656
657 if (.not. haveGtypeCutoffMap) then
658 call createGtypeCutoffMap()
659 endif
660
661
662 if (VisitCutoffsAfterComputing) then
663 call set_switch(GROUP_SWITCH, largestRcut, largestRcut)
664 endif
665
666
667 if (.not. haveSIMvariables) then
668 call setSimVariables()
669 endif
670
671 ! if (.not. haveRlist) then
672 ! write(default_error, *) 'rList has not been set in doForces!'
673 ! error = -1
674 ! return
675 ! endif
676
677 if (.not. haveNeighborList) then
678 write(default_error, *) 'neighbor list has not been initialized in doForces!'
679 error = -1
680 return
681 end if
682
683 if (.not. haveSaneForceField) then
684 write(default_error, *) 'Force Field is not sane in doForces!'
685 error = -1
686 return
687 end if
688
689 #ifdef IS_MPI
690 if (.not. isMPISimSet()) then
691 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
692 error = -1
693 return
694 endif
695 #endif
696 return
697 end subroutine doReadyCheck
698
699
700 subroutine init_FF(thisStat)
701
702 integer, intent(out) :: thisStat
703 integer :: my_status, nMatches
704 integer, pointer :: MatchList(:) => null()
705
706 !! assume things are copacetic, unless they aren't
707 thisStat = 0
708
709 !! init_FF is called *after* all of the atom types have been
710 !! defined in atype_module using the new_atype subroutine.
711 !!
712 !! this will scan through the known atypes and figure out what
713 !! interactions are used by the force field.
714
715 FF_uses_DirectionalAtoms = .false.
716 FF_uses_Dipoles = .false.
717 FF_uses_GayBerne = .false.
718 FF_uses_EAM = .false.
719
720 call getMatchingElementList(atypes, "is_Directional", .true., &
721 nMatches, MatchList)
722 if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
723
724 call getMatchingElementList(atypes, "is_Dipole", .true., &
725 nMatches, MatchList)
726 if (nMatches .gt. 0) FF_uses_Dipoles = .true.
727
728 call getMatchingElementList(atypes, "is_GayBerne", .true., &
729 nMatches, MatchList)
730 if (nMatches .gt. 0) FF_uses_GayBerne = .true.
731
732 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
733 if (nMatches .gt. 0) FF_uses_EAM = .true.
734
735
736 haveSaneForceField = .true.
737
738 if (FF_uses_EAM) then
739 call init_EAM_FF(my_status)
740 if (my_status /= 0) then
741 write(default_error, *) "init_EAM_FF returned a bad status"
742 thisStat = -1
743 haveSaneForceField = .false.
744 return
745 end if
746 endif
747
748 if (.not. haveNeighborList) then
749 !! Create neighbor lists
750 call expandNeighborList(nLocal, my_status)
751 if (my_Status /= 0) then
752 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
753 thisStat = -1
754 return
755 endif
756 haveNeighborList = .true.
757 endif
758
759 end subroutine init_FF
760
761
762 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
763 !------------------------------------------------------------->
764 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
765 do_pot_c, do_stress_c, error)
766 !! Position array provided by C, dimensioned by getNlocal
767 real ( kind = dp ), dimension(3, nLocal) :: q
768 !! molecular center-of-mass position array
769 real ( kind = dp ), dimension(3, nGroups) :: q_group
770 !! Rotation Matrix for each long range particle in simulation.
771 real( kind = dp), dimension(9, nLocal) :: A
772 !! Unit vectors for dipoles (lab frame)
773 real( kind = dp ), dimension(9,nLocal) :: eFrame
774 !! Force array provided by C, dimensioned by getNlocal
775 real ( kind = dp ), dimension(3,nLocal) :: f
776 !! Torsion array provided by C, dimensioned by getNlocal
777 real( kind = dp ), dimension(3,nLocal) :: t
778
779 !! Stress Tensor
780 real( kind = dp), dimension(9) :: tau
781 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
782 logical ( kind = 2) :: do_pot_c, do_stress_c
783 logical :: do_pot
784 logical :: do_stress
785 logical :: in_switching_region
786 #ifdef IS_MPI
787 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
788 integer :: nAtomsInRow
789 integer :: nAtomsInCol
790 integer :: nprocs
791 integer :: nGroupsInRow
792 integer :: nGroupsInCol
793 #endif
794 integer :: natoms
795 logical :: update_nlist
796 integer :: i, j, jstart, jend, jnab
797 integer :: istart, iend
798 integer :: ia, jb, atom1, atom2
799 integer :: nlist
800 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
801 real( kind = DP ) :: sw, dswdr, swderiv, mf
802 real( kind = DP ) :: rVal
803 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij
804 real(kind=dp) :: rfpot, mu_i, virial
805 real(kind=dp):: rCut
806 integer :: me_i, me_j, n_in_i, n_in_j
807 logical :: is_dp_i
808 integer :: neighborListSize
809 integer :: listerror, error
810 integer :: localError
811 integer :: propPack_i, propPack_j
812 integer :: loopStart, loopEnd, loop
813 integer :: iHash
814 integer :: i1
815
816
817 !! initialize local variables
818
819 #ifdef IS_MPI
820 pot_local = 0.0_dp
821 nAtomsInRow = getNatomsInRow(plan_atom_row)
822 nAtomsInCol = getNatomsInCol(plan_atom_col)
823 nGroupsInRow = getNgroupsInRow(plan_group_row)
824 nGroupsInCol = getNgroupsInCol(plan_group_col)
825 #else
826 natoms = nlocal
827 #endif
828
829 call doReadyCheck(localError)
830 if ( localError .ne. 0 ) then
831 call handleError("do_force_loop", "Not Initialized")
832 error = -1
833 return
834 end if
835 call zero_work_arrays()
836
837 do_pot = do_pot_c
838 do_stress = do_stress_c
839
840 ! Gather all information needed by all force loops:
841
842 #ifdef IS_MPI
843
844 call gather(q, q_Row, plan_atom_row_3d)
845 call gather(q, q_Col, plan_atom_col_3d)
846
847 call gather(q_group, q_group_Row, plan_group_row_3d)
848 call gather(q_group, q_group_Col, plan_group_col_3d)
849
850 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
851 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
852 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
853
854 call gather(A, A_Row, plan_atom_row_rotation)
855 call gather(A, A_Col, plan_atom_col_rotation)
856 endif
857
858 #endif
859
860 !! Begin force loop timing:
861 #ifdef PROFILE
862 call cpu_time(forceTimeInitial)
863 nloops = nloops + 1
864 #endif
865
866 loopEnd = PAIR_LOOP
867 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
868 loopStart = PREPAIR_LOOP
869 else
870 loopStart = PAIR_LOOP
871 endif
872
873 do loop = loopStart, loopEnd
874
875 ! See if we need to update neighbor lists
876 ! (but only on the first time through):
877 if (loop .eq. loopStart) then
878 #ifdef IS_MPI
879 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
880 update_nlist)
881 #else
882 call checkNeighborList(nGroups, q_group, skinThickness, &
883 update_nlist)
884 #endif
885 endif
886
887 if (update_nlist) then
888 !! save current configuration and construct neighbor list
889 #ifdef IS_MPI
890 call saveNeighborList(nGroupsInRow, q_group_row)
891 #else
892 call saveNeighborList(nGroups, q_group)
893 #endif
894 neighborListSize = size(list)
895 nlist = 0
896 endif
897
898 istart = 1
899 #ifdef IS_MPI
900 iend = nGroupsInRow
901 #else
902 iend = nGroups - 1
903 #endif
904 outer: do i = istart, iend
905
906 if (update_nlist) point(i) = nlist + 1
907
908 n_in_i = groupStartRow(i+1) - groupStartRow(i)
909
910 if (update_nlist) then
911 #ifdef IS_MPI
912 jstart = 1
913 jend = nGroupsInCol
914 #else
915 jstart = i+1
916 jend = nGroups
917 #endif
918 else
919 jstart = point(i)
920 jend = point(i+1) - 1
921 ! make sure group i has neighbors
922 if (jstart .gt. jend) cycle outer
923 endif
924
925 do jnab = jstart, jend
926 if (update_nlist) then
927 j = jnab
928 else
929 j = list(jnab)
930 endif
931
932 #ifdef IS_MPI
933 me_j = atid_col(j)
934 call get_interatomic_vector(q_group_Row(:,i), &
935 q_group_Col(:,j), d_grp, rgrpsq)
936 #else
937 me_j = atid(j)
938 call get_interatomic_vector(q_group(:,i), &
939 q_group(:,j), d_grp, rgrpsq)
940 #endif
941
942 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then
943 if (update_nlist) then
944 nlist = nlist + 1
945
946 if (nlist > neighborListSize) then
947 #ifdef IS_MPI
948 call expandNeighborList(nGroupsInRow, listerror)
949 #else
950 call expandNeighborList(nGroups, listerror)
951 #endif
952 if (listerror /= 0) then
953 error = -1
954 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
955 return
956 end if
957 neighborListSize = size(list)
958 endif
959
960 list(nlist) = j
961 endif
962
963
964
965 if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then
966
967 rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut
968 if (loop .eq. PAIR_LOOP) then
969 vij = 0.0d0
970 fij(1:3) = 0.0d0
971 endif
972
973 call get_switch(rgrpsq, sw, dswdr, rgrp, &
974 group_switch, in_switching_region)
975
976 n_in_j = groupStartCol(j+1) - groupStartCol(j)
977
978 do ia = groupStartRow(i), groupStartRow(i+1)-1
979
980 atom1 = groupListRow(ia)
981
982 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
983
984 atom2 = groupListCol(jb)
985
986 if (skipThisPair(atom1, atom2)) cycle inner
987
988 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
989 d_atm(1:3) = d_grp(1:3)
990 ratmsq = rgrpsq
991 else
992 #ifdef IS_MPI
993 call get_interatomic_vector(q_Row(:,atom1), &
994 q_Col(:,atom2), d_atm, ratmsq)
995 #else
996 call get_interatomic_vector(q(:,atom1), &
997 q(:,atom2), d_atm, ratmsq)
998 #endif
999 endif
1000
1001 if (loop .eq. PREPAIR_LOOP) then
1002 #ifdef IS_MPI
1003 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1004 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1005 eFrame, A, f, t, pot_local)
1006 #else
1007 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
1008 rgrpsq, d_grp, rCut, do_pot, do_stress, &
1009 eFrame, A, f, t, pot)
1010 #endif
1011 else
1012 #ifdef IS_MPI
1013 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1014 do_pot, eFrame, A, f, t, pot_local, vpair, &
1015 fpair, d_grp, rgrp, rCut)
1016 #else
1017 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
1018 do_pot, eFrame, A, f, t, pot, vpair, fpair, &
1019 d_grp, rgrp, rCut)
1020 #endif
1021 vij = vij + vpair
1022 fij(1:3) = fij(1:3) + fpair(1:3)
1023 endif
1024 enddo inner
1025 enddo
1026
1027 if (loop .eq. PAIR_LOOP) then
1028 if (in_switching_region) then
1029 swderiv = vij*dswdr/rgrp
1030 fij(1) = fij(1) + swderiv*d_grp(1)
1031 fij(2) = fij(2) + swderiv*d_grp(2)
1032 fij(3) = fij(3) + swderiv*d_grp(3)
1033
1034 do ia=groupStartRow(i), groupStartRow(i+1)-1
1035 atom1=groupListRow(ia)
1036 mf = mfactRow(atom1)
1037 #ifdef IS_MPI
1038 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
1039 f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
1040 f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
1041 #else
1042 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
1043 f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
1044 f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
1045 #endif
1046 enddo
1047
1048 do jb=groupStartCol(j), groupStartCol(j+1)-1
1049 atom2=groupListCol(jb)
1050 mf = mfactCol(atom2)
1051 #ifdef IS_MPI
1052 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
1053 f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
1054 f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
1055 #else
1056 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
1057 f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
1058 f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
1059 #endif
1060 enddo
1061 endif
1062
1063 if (do_stress) call add_stress_tensor(d_grp, fij)
1064 endif
1065 endif
1066 endif
1067 enddo
1068
1069 enddo outer
1070
1071 if (update_nlist) then
1072 #ifdef IS_MPI
1073 point(nGroupsInRow + 1) = nlist + 1
1074 #else
1075 point(nGroups) = nlist + 1
1076 #endif
1077 if (loop .eq. PREPAIR_LOOP) then
1078 ! we just did the neighbor list update on the first
1079 ! pass, so we don't need to do it
1080 ! again on the second pass
1081 update_nlist = .false.
1082 endif
1083 endif
1084
1085 if (loop .eq. PREPAIR_LOOP) then
1086 call do_preforce(nlocal, pot)
1087 endif
1088
1089 enddo
1090
1091 !! Do timing
1092 #ifdef PROFILE
1093 call cpu_time(forceTimeFinal)
1094 forceTime = forceTime + forceTimeFinal - forceTimeInitial
1095 #endif
1096
1097 #ifdef IS_MPI
1098 !!distribute forces
1099
1100 f_temp = 0.0_dp
1101 call scatter(f_Row,f_temp,plan_atom_row_3d)
1102 do i = 1,nlocal
1103 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1104 end do
1105
1106 f_temp = 0.0_dp
1107 call scatter(f_Col,f_temp,plan_atom_col_3d)
1108 do i = 1,nlocal
1109 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1110 end do
1111
1112 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
1113 t_temp = 0.0_dp
1114 call scatter(t_Row,t_temp,plan_atom_row_3d)
1115 do i = 1,nlocal
1116 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1117 end do
1118 t_temp = 0.0_dp
1119 call scatter(t_Col,t_temp,plan_atom_col_3d)
1120
1121 do i = 1,nlocal
1122 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1123 end do
1124 endif
1125
1126 if (do_pot) then
1127 ! scatter/gather pot_row into the members of my column
1128 do i = 1,LR_POT_TYPES
1129 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
1130 end do
1131 ! scatter/gather pot_local into all other procs
1132 ! add resultant to get total pot
1133 do i = 1, nlocal
1134 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
1135 + pot_Temp(1:LR_POT_TYPES,i)
1136 enddo
1137
1138 pot_Temp = 0.0_DP
1139 do i = 1,LR_POT_TYPES
1140 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
1141 end do
1142 do i = 1, nlocal
1143 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
1144 + pot_Temp(1:LR_POT_TYPES,i)
1145 enddo
1146
1147 endif
1148 #endif
1149
1150 if (SIM_requires_postpair_calc) then
1151 do i = 1, nlocal
1152
1153 ! we loop only over the local atoms, so we don't need row and column
1154 ! lookups for the types
1155
1156 me_i = atid(i)
1157
1158 ! is the atom electrostatic? See if it would have an
1159 ! electrostatic interaction with itself
1160 iHash = InteractionHash(me_i,me_i)
1161
1162 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1163 #ifdef IS_MPI
1164 call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), &
1165 t, do_pot)
1166 #else
1167 call self_self(i, eFrame, pot(ELECTROSTATIC_POT), &
1168 t, do_pot)
1169 #endif
1170 endif
1171
1172
1173 if (electrostaticSummationMethod.eq.REACTION_FIELD) then
1174
1175 ! loop over the excludes to accumulate RF stuff we've
1176 ! left out of the normal pair loop
1177
1178 do i1 = 1, nSkipsForAtom(i)
1179 j = skipsForAtom(i, i1)
1180
1181 ! prevent overcounting of the skips
1182 if (i.lt.j) then
1183 call get_interatomic_vector(q(:,i), &
1184 q(:,j), d_atm, ratmsq)
1185 rVal = dsqrt(ratmsq)
1186 call get_switch(ratmsq, sw, dswdr, rVal, group_switch, &
1187 in_switching_region)
1188 #ifdef IS_MPI
1189 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1190 vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot)
1191 #else
1192 call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, &
1193 vpair, pot(ELECTROSTATIC_POT), f, t, do_pot)
1194 #endif
1195 endif
1196 enddo
1197 endif
1198 enddo
1199 endif
1200
1201 #ifdef IS_MPI
1202
1203 if (do_pot) then
1204 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision,mpi_sum, &
1205 mpi_comm_world,mpi_err)
1206 endif
1207
1208 if (do_stress) then
1209 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
1210 mpi_comm_world,mpi_err)
1211 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
1212 mpi_comm_world,mpi_err)
1213 endif
1214
1215 #else
1216
1217 if (do_stress) then
1218 tau = tau_Temp
1219 virial = virial_Temp
1220 endif
1221
1222 #endif
1223
1224 end subroutine do_force_loop
1225
1226 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
1227 eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp, rCut)
1228
1229 real( kind = dp ) :: vpair, sw
1230 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1231 real( kind = dp ), dimension(3) :: fpair
1232 real( kind = dp ), dimension(nLocal) :: mfact
1233 real( kind = dp ), dimension(9,nLocal) :: eFrame
1234 real( kind = dp ), dimension(9,nLocal) :: A
1235 real( kind = dp ), dimension(3,nLocal) :: f
1236 real( kind = dp ), dimension(3,nLocal) :: t
1237
1238 logical, intent(inout) :: do_pot
1239 integer, intent(in) :: i, j
1240 real ( kind = dp ), intent(inout) :: rijsq
1241 real ( kind = dp ), intent(inout) :: r_grp
1242 real ( kind = dp ), intent(inout) :: d(3)
1243 real ( kind = dp ), intent(inout) :: d_grp(3)
1244 real ( kind = dp ), intent(inout) :: rCut
1245 real ( kind = dp ) :: r
1246 integer :: me_i, me_j
1247
1248 integer :: iHash
1249
1250 r = sqrt(rijsq)
1251 vpair = 0.0d0
1252 fpair(1:3) = 0.0d0
1253
1254 #ifdef IS_MPI
1255 me_i = atid_row(i)
1256 me_j = atid_col(j)
1257 #else
1258 me_i = atid(i)
1259 me_j = atid(j)
1260 #endif
1261
1262 iHash = InteractionHash(me_i, me_j)
1263
1264 if ( iand(iHash, LJ_PAIR).ne.0 ) then
1265 call do_lj_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1266 pot(VDW_POT), f, do_pot)
1267 endif
1268
1269 if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then
1270 call doElectrostaticPair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1271 pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot)
1272 endif
1273
1274 if ( iand(iHash, STICKY_PAIR).ne.0 ) then
1275 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1276 pot(HB_POT), A, f, t, do_pot)
1277 endif
1278
1279 if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then
1280 call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1281 pot(HB_POT), A, f, t, do_pot)
1282 endif
1283
1284 if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then
1285 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1286 pot(VDW_POT), A, f, t, do_pot)
1287 endif
1288
1289 if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then
1290 call do_gb_lj_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1291 pot(VDW_POT), A, f, t, do_pot)
1292 endif
1293
1294 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1295 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1296 pot(METALLIC_POT), f, do_pot)
1297 endif
1298
1299 if ( iand(iHash, SHAPE_PAIR).ne.0 ) then
1300 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1301 pot(VDW_POT), A, f, t, do_pot)
1302 endif
1303
1304 if ( iand(iHash, SHAPE_LJ).ne.0 ) then
1305 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
1306 pot(VDW_POT), A, f, t, do_pot)
1307 endif
1308
1309 if ( iand(iHash, SC_PAIR).ne.0 ) then
1310 call do_SC_pair(i, j, d, r, rijsq, rcut, sw, vpair, fpair, &
1311 pot(METALLIC_POT), f, do_pot)
1312 endif
1313
1314
1315
1316 end subroutine do_pair
1317
1318 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
1319 do_pot, do_stress, eFrame, A, f, t, pot)
1320
1321 real( kind = dp ) :: sw
1322 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
1323 real( kind = dp ), dimension(9,nLocal) :: eFrame
1324 real (kind=dp), dimension(9,nLocal) :: A
1325 real (kind=dp), dimension(3,nLocal) :: f
1326 real (kind=dp), dimension(3,nLocal) :: t
1327
1328 logical, intent(inout) :: do_pot, do_stress
1329 integer, intent(in) :: i, j
1330 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
1331 real ( kind = dp ) :: r, rc
1332 real ( kind = dp ), intent(inout) :: d(3), dc(3)
1333
1334 integer :: me_i, me_j, iHash
1335
1336 r = sqrt(rijsq)
1337
1338 #ifdef IS_MPI
1339 me_i = atid_row(i)
1340 me_j = atid_col(j)
1341 #else
1342 me_i = atid(i)
1343 me_j = atid(j)
1344 #endif
1345
1346 iHash = InteractionHash(me_i, me_j)
1347
1348 if ( iand(iHash, EAM_PAIR).ne.0 ) then
1349 call calc_EAM_prepair_rho(i, j, d, r, rijsq)
1350 endif
1351
1352 if ( iand(iHash, SC_PAIR).ne.0 ) then
1353 call calc_SC_prepair_rho(i, j, d, r, rijsq, rcut )
1354 endif
1355
1356 end subroutine do_prepair
1357
1358
1359 subroutine do_preforce(nlocal,pot)
1360 integer :: nlocal
1361 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
1362
1363 if (FF_uses_EAM .and. SIM_uses_EAM) then
1364 call calc_EAM_preforce_Frho(nlocal,pot(METALLIC_POT))
1365 endif
1366 if (FF_uses_SC .and. SIM_uses_SC) then
1367 call calc_SC_preforce_Frho(nlocal,pot(METALLIC_POT))
1368 endif
1369
1370
1371 end subroutine do_preforce
1372
1373
1374 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1375
1376 real (kind = dp), dimension(3) :: q_i
1377 real (kind = dp), dimension(3) :: q_j
1378 real ( kind = dp ), intent(out) :: r_sq
1379 real( kind = dp ) :: d(3), scaled(3)
1380 integer i
1381
1382 d(1:3) = q_j(1:3) - q_i(1:3)
1383
1384 ! Wrap back into periodic box if necessary
1385 if ( SIM_uses_PBC ) then
1386
1387 if( .not.boxIsOrthorhombic ) then
1388 ! calc the scaled coordinates.
1389
1390 scaled = matmul(HmatInv, d)
1391
1392 ! wrap the scaled coordinates
1393
1394 scaled = scaled - anint(scaled)
1395
1396
1397 ! calc the wrapped real coordinates from the wrapped scaled
1398 ! coordinates
1399
1400 d = matmul(Hmat,scaled)
1401
1402 else
1403 ! calc the scaled coordinates.
1404
1405 do i = 1, 3
1406 scaled(i) = d(i) * HmatInv(i,i)
1407
1408 ! wrap the scaled coordinates
1409
1410 scaled(i) = scaled(i) - anint(scaled(i))
1411
1412 ! calc the wrapped real coordinates from the wrapped scaled
1413 ! coordinates
1414
1415 d(i) = scaled(i)*Hmat(i,i)
1416 enddo
1417 endif
1418
1419 endif
1420
1421 r_sq = dot_product(d,d)
1422
1423 end subroutine get_interatomic_vector
1424
1425 subroutine zero_work_arrays()
1426
1427 #ifdef IS_MPI
1428
1429 q_Row = 0.0_dp
1430 q_Col = 0.0_dp
1431
1432 q_group_Row = 0.0_dp
1433 q_group_Col = 0.0_dp
1434
1435 eFrame_Row = 0.0_dp
1436 eFrame_Col = 0.0_dp
1437
1438 A_Row = 0.0_dp
1439 A_Col = 0.0_dp
1440
1441 f_Row = 0.0_dp
1442 f_Col = 0.0_dp
1443 f_Temp = 0.0_dp
1444
1445 t_Row = 0.0_dp
1446 t_Col = 0.0_dp
1447 t_Temp = 0.0_dp
1448
1449 pot_Row = 0.0_dp
1450 pot_Col = 0.0_dp
1451 pot_Temp = 0.0_dp
1452
1453 #endif
1454
1455 if (FF_uses_EAM .and. SIM_uses_EAM) then
1456 call clean_EAM()
1457 endif
1458
1459 tau_Temp = 0.0_dp
1460 virial_Temp = 0.0_dp
1461 end subroutine zero_work_arrays
1462
1463 function skipThisPair(atom1, atom2) result(skip_it)
1464 integer, intent(in) :: atom1
1465 integer, intent(in), optional :: atom2
1466 logical :: skip_it
1467 integer :: unique_id_1, unique_id_2
1468 integer :: me_i,me_j
1469 integer :: i
1470
1471 skip_it = .false.
1472
1473 !! there are a number of reasons to skip a pair or a particle
1474 !! mostly we do this to exclude atoms who are involved in short
1475 !! range interactions (bonds, bends, torsions), but we also need
1476 !! to exclude some overcounted interactions that result from
1477 !! the parallel decomposition
1478
1479 #ifdef IS_MPI
1480 !! in MPI, we have to look up the unique IDs for each atom
1481 unique_id_1 = AtomRowToGlobal(atom1)
1482 #else
1483 !! in the normal loop, the atom numbers are unique
1484 unique_id_1 = atom1
1485 #endif
1486
1487 !! We were called with only one atom, so just check the global exclude
1488 !! list for this atom
1489 if (.not. present(atom2)) then
1490 do i = 1, nExcludes_global
1491 if (excludesGlobal(i) == unique_id_1) then
1492 skip_it = .true.
1493 return
1494 end if
1495 end do
1496 return
1497 end if
1498
1499 #ifdef IS_MPI
1500 unique_id_2 = AtomColToGlobal(atom2)
1501 #else
1502 unique_id_2 = atom2
1503 #endif
1504
1505 #ifdef IS_MPI
1506 !! this situation should only arise in MPI simulations
1507 if (unique_id_1 == unique_id_2) then
1508 skip_it = .true.
1509 return
1510 end if
1511
1512 !! this prevents us from doing the pair on multiple processors
1513 if (unique_id_1 < unique_id_2) then
1514 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1515 skip_it = .true.
1516 return
1517 endif
1518 else
1519 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1520 skip_it = .true.
1521 return
1522 endif
1523 endif
1524 #endif
1525
1526 !! the rest of these situations can happen in all simulations:
1527 do i = 1, nExcludes_global
1528 if ((excludesGlobal(i) == unique_id_1) .or. &
1529 (excludesGlobal(i) == unique_id_2)) then
1530 skip_it = .true.
1531 return
1532 endif
1533 enddo
1534
1535 do i = 1, nSkipsForAtom(atom1)
1536 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1537 skip_it = .true.
1538 return
1539 endif
1540 end do
1541
1542 return
1543 end function skipThisPair
1544
1545 function FF_UsesDirectionalAtoms() result(doesit)
1546 logical :: doesit
1547 doesit = FF_uses_DirectionalAtoms
1548 end function FF_UsesDirectionalAtoms
1549
1550 function FF_RequiresPrepairCalc() result(doesit)
1551 logical :: doesit
1552 doesit = FF_uses_EAM .or. FF_uses_SC &
1553 .or. FF_uses_MEAM
1554 end function FF_RequiresPrepairCalc
1555
1556 #ifdef PROFILE
1557 function getforcetime() result(totalforcetime)
1558 real(kind=dp) :: totalforcetime
1559 totalforcetime = forcetime
1560 end function getforcetime
1561 #endif
1562
1563 !! This cleans componets of force arrays belonging only to fortran
1564
1565 subroutine add_stress_tensor(dpair, fpair)
1566
1567 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1568
1569 ! because the d vector is the rj - ri vector, and
1570 ! because fx, fy, fz are the force on atom i, we need a
1571 ! negative sign here:
1572
1573 tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1574 tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1575 tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1576 tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1577 tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1578 tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1579 tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1580 tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1581 tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1582
1583 virial_Temp = virial_Temp + &
1584 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1585
1586 end subroutine add_stress_tensor
1587
1588 end module doForces