ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 1927
Committed: Wed Jan 12 17:14:35 2005 UTC (21 years ago) by tim
File size: 38883 byte(s)
Log Message:
change license

File Contents

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