ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2277
Committed: Fri Aug 26 21:30:41 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 39324 byte(s)
Log Message:
added some probably nonfunctional get*cut routines

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