ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/doForces.F90
Revision: 1634
Committed: Fri Oct 22 21:21:02 2004 UTC (19 years, 8 months ago) by gezelter
File size: 36650 byte(s)
Log Message:
fixey fixey

File Contents

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