ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2261
Committed: Tue Jun 28 13:58:45 2005 UTC (19 years ago) by gezelter
File size: 38404 byte(s)
Log Message:
*** empty log message ***

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