ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2432
Committed: Tue Nov 15 16:01:06 2005 UTC (18 years, 7 months ago) by chuckv
File size: 48120 byte(s)
Log Message:
Added Sutton-Chen to force loop...

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