ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2405
Committed: Tue Nov 1 19:24:57 2005 UTC (18 years, 8 months ago) by chrisfen
File size: 46497 byte(s)
Log Message:
removed some testing code...

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