ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2461
Committed: Mon Nov 21 22:59:02 2005 UTC (18 years, 7 months ago) by gezelter
File size: 49328 byte(s)
Log Message:
Cutoff mixing fixes have been made.

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