ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 2274
Committed: Wed Aug 17 15:26:42 2005 UTC (18 years, 10 months ago) by gezelter
File size: 38984 byte(s)
Log Message:
added fCutoffPolicy.h

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