ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2756
Committed: Wed May 17 15:37:15 2006 UTC (18 years, 1 month ago) by gezelter
File size: 50940 byte(s)
Log Message:
Getting fortran side prepped for single precision...

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