ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/do_Forces.F90
Revision: 1334
Committed: Fri Jul 16 18:58:03 2004 UTC (19 years, 11 months ago) by gezelter
File size: 33783 byte(s)
Log Message:
Initial import of OOPSE-1.0 source tree

File Contents

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