ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 2211
Committed: Thu Apr 21 14:12:19 2005 UTC (19 years, 2 months ago) by chrisfen
File size: 37511 byte(s)
Log Message:
Shapes is limping along with a array bounds overwrite (I think...). At least the parser loads the forcefield fine...

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