ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 1948
Committed: Fri Jan 14 20:31:16 2005 UTC (19 years, 5 months ago) by gezelter
File size: 37543 byte(s)
Log Message:
separating modules and C/Fortran interface subroutines

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