ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1206
Committed: Thu May 27 19:51:18 2004 UTC (20 years, 1 month ago) by tim
File size: 34235 byte(s)
Log Message:
Bug fix for SkipList

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.65 2004-05-27 19:51:18 tim Exp $, $Date: 2004-05-27 19:51:18 $, $Name: not supported by cvs2svn $, $Revision: 1.65 $
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 call checkNeighborList(nGroups, q_group, listSkin, update_nlist)
485 endif
486
487 if (update_nlist) then
488 !! save current configuration and construct neighbor list
489 #ifdef IS_MPI
490 call saveNeighborList(nGroupsInRow, q_group)
491 #else
492 call saveNeighborList(nGroups, q_group)
493 #endif
494 neighborListSize = size(list)
495 nlist = 0
496 endif
497
498 istart = 1
499 #ifdef IS_MPI
500 iend = nGroupsInRow
501 #else
502 iend = nGroups - 1
503 #endif
504 outer: do i = istart, iend
505
506 if (update_nlist) point(i) = nlist + 1
507
508 n_in_i = groupStartRow(i+1) - groupStartRow(i)
509
510 if (update_nlist) then
511 #ifdef IS_MPI
512 jstart = 1
513 jend = nGroupsInCol
514 #else
515 jstart = i+1
516 jend = nGroups
517 #endif
518 else
519 jstart = point(i)
520 jend = point(i+1) - 1
521 ! make sure group i has neighbors
522 if (jstart .gt. jend) cycle outer
523 endif
524
525 do jnab = jstart, jend
526 if (update_nlist) then
527 j = jnab
528 else
529 j = list(jnab)
530 endif
531
532 #ifdef IS_MPI
533 call get_interatomic_vector(q_group_Row(:,i), &
534 q_group_Col(:,j), d_grp, rgrpsq)
535 #else
536 call get_interatomic_vector(q_group(:,i), &
537 q_group(:,j), d_grp, rgrpsq)
538 #endif
539
540 if (rgrpsq < rlistsq) then
541 if (update_nlist) then
542 nlist = nlist + 1
543
544 if (nlist > neighborListSize) then
545 #ifdef IS_MPI
546 call expandNeighborList(nGroupsInRow, listerror)
547 #else
548 call expandNeighborList(nGroups, listerror)
549 #endif
550 if (listerror /= 0) then
551 error = -1
552 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
553 return
554 end if
555 neighborListSize = size(list)
556 endif
557
558 list(nlist) = j
559 endif
560
561 if (loop .eq. PAIR_LOOP) then
562 vij = 0.0d0
563 fij(1:3) = 0.0d0
564 endif
565
566 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
567 in_switching_region)
568
569 n_in_j = groupStartCol(j+1) - groupStartCol(j)
570
571 do ia = groupStartRow(i), groupStartRow(i+1)-1
572
573 atom1 = groupListRow(ia)
574
575 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
576
577 atom2 = groupListCol(jb)
578
579 if (skipThisPair(atom1, atom2)) cycle inner
580
581 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
582 d_atm(1:3) = d_grp(1:3)
583 ratmsq = rgrpsq
584 else
585 #ifdef IS_MPI
586 call get_interatomic_vector(q_Row(:,atom1), &
587 q_Col(:,atom2), d_atm, ratmsq)
588 #else
589 call get_interatomic_vector(q(:,atom1), &
590 q(:,atom2), d_atm, ratmsq)
591 #endif
592 endif
593 if (loop .eq. PREPAIR_LOOP) then
594 #ifdef IS_MPI
595 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
596 rgrpsq, d_grp, do_pot, do_stress, &
597 u_l, A, f, t, pot_local)
598 #else
599 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
600 rgrpsq, d_grp, do_pot, do_stress, &
601 u_l, A, f, t, pot)
602 #endif
603 else
604 #ifdef IS_MPI
605 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
606 do_pot, &
607 u_l, A, f, t, pot_local, vpair, fpair)
608 #else
609 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
610 do_pot, &
611 u_l, A, f, t, pot, vpair, fpair)
612 #endif
613 vij = vij + vpair
614 fij(1:3) = fij(1:3) + fpair(1:3)
615 endif
616 enddo inner
617 enddo
618
619 if (loop .eq. PAIR_LOOP) then
620 if (in_switching_region) then
621 swderiv = vij*dswdr/rgrp
622 fij(1) = fij(1) + swderiv*d_grp(1)
623 fij(2) = fij(2) + swderiv*d_grp(2)
624 fij(3) = fij(3) + swderiv*d_grp(3)
625
626 do ia=groupStartRow(i), groupStartRow(i+1)-1
627 atom1=groupListRow(ia)
628 mf = mfactRow(atom1)
629 #ifdef IS_MPI
630 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
631 f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
632 f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
633 #else
634 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
635 f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
636 f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
637 #endif
638 enddo
639
640 do jb=groupStartCol(j), groupStartCol(j+1)-1
641 atom2=groupListCol(jb)
642 mf = mfactCol(atom2)
643 #ifdef IS_MPI
644 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
645 f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
646 f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
647 #else
648 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
649 f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
650 f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
651 #endif
652 enddo
653 endif
654
655 if (do_stress) call add_stress_tensor(d_grp, fij)
656 endif
657 end if
658 enddo
659 enddo outer
660
661 if (update_nlist) then
662 #ifdef IS_MPI
663 point(nGroupsInRow + 1) = nlist + 1
664 #else
665 point(nGroups) = nlist + 1
666 #endif
667 if (loop .eq. PREPAIR_LOOP) then
668 ! we just did the neighbor list update on the first
669 ! pass, so we don't need to do it
670 ! again on the second pass
671 update_nlist = .false.
672 endif
673 endif
674
675 if (loop .eq. PREPAIR_LOOP) then
676 call do_preforce(nlocal, pot)
677 endif
678
679 enddo
680
681 !! Do timing
682 #ifdef PROFILE
683 call cpu_time(forceTimeFinal)
684 forceTime = forceTime + forceTimeFinal - forceTimeInitial
685 #endif
686
687 #ifdef IS_MPI
688 !!distribute forces
689
690 f_temp = 0.0_dp
691 call scatter(f_Row,f_temp,plan_atom_row_3d)
692 do i = 1,nlocal
693 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
694 end do
695
696 f_temp = 0.0_dp
697 call scatter(f_Col,f_temp,plan_atom_col_3d)
698 do i = 1,nlocal
699 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
700 end do
701
702 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
703 t_temp = 0.0_dp
704 call scatter(t_Row,t_temp,plan_atom_row_3d)
705 do i = 1,nlocal
706 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
707 end do
708 t_temp = 0.0_dp
709 call scatter(t_Col,t_temp,plan_atom_col_3d)
710
711 do i = 1,nlocal
712 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
713 end do
714 endif
715
716 if (do_pot) then
717 ! scatter/gather pot_row into the members of my column
718 call scatter(pot_Row, pot_Temp, plan_atom_row)
719
720 ! scatter/gather pot_local into all other procs
721 ! add resultant to get total pot
722 do i = 1, nlocal
723 pot_local = pot_local + pot_Temp(i)
724 enddo
725
726 pot_Temp = 0.0_DP
727
728 call scatter(pot_Col, pot_Temp, plan_atom_col)
729 do i = 1, nlocal
730 pot_local = pot_local + pot_Temp(i)
731 enddo
732
733 endif
734 #endif
735
736 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
737
738 if (FF_uses_RF .and. SIM_uses_RF) then
739
740 #ifdef IS_MPI
741 call scatter(rf_Row,rf,plan_atom_row_3d)
742 call scatter(rf_Col,rf_Temp,plan_atom_col_3d)
743 do i = 1,nlocal
744 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
745 end do
746 #endif
747
748 do i = 1, nLocal
749
750 rfpot = 0.0_DP
751 #ifdef IS_MPI
752 me_i = atid_row(i)
753 #else
754 me_i = atid(i)
755 #endif
756
757 if (PropertyMap(me_i)%is_DP) then
758
759 mu_i = PropertyMap(me_i)%dipole_moment
760
761 !! The reaction field needs to include a self contribution
762 !! to the field:
763 call accumulate_self_rf(i, mu_i, u_l)
764 !! Get the reaction field contribution to the
765 !! potential and torques:
766 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
767 #ifdef IS_MPI
768 pot_local = pot_local + rfpot
769 #else
770 pot = pot + rfpot
771
772 #endif
773 endif
774 enddo
775 endif
776 endif
777
778
779 #ifdef IS_MPI
780
781 if (do_pot) then
782 pot = pot + pot_local
783 !! we assume the c code will do the allreduce to get the total potential
784 !! we could do it right here if we needed to...
785 endif
786
787 if (do_stress) then
788 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
789 mpi_comm_world,mpi_err)
790 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
791 mpi_comm_world,mpi_err)
792 endif
793
794 #else
795
796 if (do_stress) then
797 tau = tau_Temp
798 virial = virial_Temp
799 endif
800
801 #endif
802
803 end subroutine do_force_loop
804
805 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
806 u_l, A, f, t, pot, vpair, fpair)
807
808 real( kind = dp ) :: pot, vpair, sw
809 real( kind = dp ), dimension(3) :: fpair
810 real( kind = dp ), dimension(nLocal) :: mfact
811 real( kind = dp ), dimension(3,nLocal) :: u_l
812 real( kind = dp ), dimension(9,nLocal) :: A
813 real( kind = dp ), dimension(3,nLocal) :: f
814 real( kind = dp ), dimension(3,nLocal) :: t
815
816 logical, intent(inout) :: do_pot
817 integer, intent(in) :: i, j
818 real ( kind = dp ), intent(inout) :: rijsq
819 real ( kind = dp ) :: r
820 real ( kind = dp ), intent(inout) :: d(3)
821 integer :: me_i, me_j
822
823 r = sqrt(rijsq)
824 vpair = 0.0d0
825 fpair(1:3) = 0.0d0
826
827 #ifdef IS_MPI
828 if (AtomRowToGlobal(i) .eq. AtomColToGlobal(j)) then
829 write(0,*) 'do_pair is doing', i , j, AtomRowToGlobal(i), AtomColToGlobal(j)
830 endif
831 me_i = atid_row(i)
832 me_j = atid_col(j)
833 #else
834 me_i = atid(i)
835 me_j = atid(j)
836 #endif
837
838 if (FF_uses_LJ .and. SIM_uses_LJ) then
839
840 if ( PropertyMap(me_i)%is_LJ .and. PropertyMap(me_j)%is_LJ ) then
841 !write(*,*) 'calling lj with'
842 !write(*,*) i, j, r, rijsq
843 !write(*,'(3es12.3)') d(1), d(2), d(3)
844 !write(*,'(3es12.3)') sw, vpair, pot
845 !write(*,*)
846
847 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
848 endif
849
850 endif
851
852 if (FF_uses_charges .and. SIM_uses_charges) then
853
854 if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
855 call do_charge_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
856 endif
857
858 endif
859
860 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
861
862 if ( PropertyMap(me_i)%is_DP .and. PropertyMap(me_j)%is_DP) then
863 call do_dipole_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, u_l, f, t, &
864 do_pot)
865 if (FF_uses_RF .and. SIM_uses_RF) then
866 call accumulate_rf(i, j, r, u_l, sw)
867 call rf_correct_forces(i, j, d, r, u_l, sw, f, fpair)
868 endif
869 endif
870
871 endif
872
873 if (FF_uses_Sticky .and. SIM_uses_sticky) then
874
875 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
876 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, A, f, t, &
877 do_pot)
878 endif
879
880 endif
881
882
883 if (FF_uses_GB .and. SIM_uses_GB) then
884
885 if ( PropertyMap(me_i)%is_GB .and. PropertyMap(me_j)%is_GB) then
886 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, u_l, f, t, &
887 do_pot)
888 endif
889
890 endif
891
892 if (FF_uses_EAM .and. SIM_uses_EAM) then
893
894 if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
895 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
896 do_pot)
897 endif
898
899 endif
900
901 end subroutine do_pair
902
903 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
904 do_pot, do_stress, u_l, A, f, t, pot)
905
906 real( kind = dp ) :: pot, sw
907 real( kind = dp ), dimension(3,nLocal) :: u_l
908 real (kind=dp), dimension(9,nLocal) :: A
909 real (kind=dp), dimension(3,nLocal) :: f
910 real (kind=dp), dimension(3,nLocal) :: t
911
912 logical, intent(inout) :: do_pot, do_stress
913 integer, intent(in) :: i, j
914 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
915 real ( kind = dp ) :: r, rc
916 real ( kind = dp ), intent(inout) :: d(3), dc(3)
917
918 logical :: is_EAM_i, is_EAM_j
919
920 integer :: me_i, me_j
921
922
923 r = sqrt(rijsq)
924 if (SIM_uses_molecular_cutoffs) then
925 rc = sqrt(rcijsq)
926 else
927 rc = r
928 endif
929
930
931 #ifdef IS_MPI
932 if (AtomRowToGlobal(i) .eq. AtomColToGlobal(j)) then
933 write(0,*) 'do_prepair is doing', i , j, AtomRowToGlobal(i), AtomColToGlobal(j)
934 endif
935
936 me_i = atid_row(i)
937 me_j = atid_col(j)
938
939 #else
940
941 me_i = atid(i)
942 me_j = atid(j)
943
944 #endif
945
946 if (FF_uses_EAM .and. SIM_uses_EAM) then
947
948 if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
949 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
950
951 endif
952
953 end subroutine do_prepair
954
955
956 subroutine do_preforce(nlocal,pot)
957 integer :: nlocal
958 real( kind = dp ) :: pot
959
960 if (FF_uses_EAM .and. SIM_uses_EAM) then
961 call calc_EAM_preforce_Frho(nlocal,pot)
962 endif
963
964
965 end subroutine do_preforce
966
967
968 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
969
970 real (kind = dp), dimension(3) :: q_i
971 real (kind = dp), dimension(3) :: q_j
972 real ( kind = dp ), intent(out) :: r_sq
973 real( kind = dp ) :: d(3), scaled(3)
974 integer i
975
976 d(1:3) = q_j(1:3) - q_i(1:3)
977
978 ! Wrap back into periodic box if necessary
979 if ( SIM_uses_PBC ) then
980
981 if( .not.boxIsOrthorhombic ) then
982 ! calc the scaled coordinates.
983
984 scaled = matmul(HmatInv, d)
985
986 ! wrap the scaled coordinates
987
988 scaled = scaled - anint(scaled)
989
990
991 ! calc the wrapped real coordinates from the wrapped scaled
992 ! coordinates
993
994 d = matmul(Hmat,scaled)
995
996 else
997 ! calc the scaled coordinates.
998
999 do i = 1, 3
1000 scaled(i) = d(i) * HmatInv(i,i)
1001
1002 ! wrap the scaled coordinates
1003
1004 scaled(i) = scaled(i) - anint(scaled(i))
1005
1006 ! calc the wrapped real coordinates from the wrapped scaled
1007 ! coordinates
1008
1009 d(i) = scaled(i)*Hmat(i,i)
1010 enddo
1011 endif
1012
1013 endif
1014
1015 r_sq = dot_product(d,d)
1016
1017 end subroutine get_interatomic_vector
1018
1019 subroutine zero_work_arrays()
1020
1021 #ifdef IS_MPI
1022
1023 q_Row = 0.0_dp
1024 q_Col = 0.0_dp
1025
1026 q_group_Row = 0.0_dp
1027 q_group_Col = 0.0_dp
1028
1029 u_l_Row = 0.0_dp
1030 u_l_Col = 0.0_dp
1031
1032 A_Row = 0.0_dp
1033 A_Col = 0.0_dp
1034
1035 f_Row = 0.0_dp
1036 f_Col = 0.0_dp
1037 f_Temp = 0.0_dp
1038
1039 t_Row = 0.0_dp
1040 t_Col = 0.0_dp
1041 t_Temp = 0.0_dp
1042
1043 pot_Row = 0.0_dp
1044 pot_Col = 0.0_dp
1045 pot_Temp = 0.0_dp
1046
1047 rf_Row = 0.0_dp
1048 rf_Col = 0.0_dp
1049 rf_Temp = 0.0_dp
1050
1051 #endif
1052
1053 if (FF_uses_EAM .and. SIM_uses_EAM) then
1054 call clean_EAM()
1055 endif
1056
1057 rf = 0.0_dp
1058 tau_Temp = 0.0_dp
1059 virial_Temp = 0.0_dp
1060 end subroutine zero_work_arrays
1061
1062 function skipThisPair(atom1, atom2) result(skip_it)
1063 integer, intent(in) :: atom1
1064 integer, intent(in), optional :: atom2
1065 logical :: skip_it
1066 integer :: unique_id_1, unique_id_2
1067 integer :: me_i,me_j
1068 integer :: i
1069
1070 skip_it = .false.
1071
1072 !! there are a number of reasons to skip a pair or a particle
1073 !! mostly we do this to exclude atoms who are involved in short
1074 !! range interactions (bonds, bends, torsions), but we also need
1075 !! to exclude some overcounted interactions that result from
1076 !! the parallel decomposition
1077
1078 #ifdef IS_MPI
1079 !! in MPI, we have to look up the unique IDs for each atom
1080 unique_id_1 = AtomRowToGlobal(atom1)
1081 #else
1082 !! in the normal loop, the atom numbers are unique
1083 unique_id_1 = atom1
1084 #endif
1085
1086 !! We were called with only one atom, so just check the global exclude
1087 !! list for this atom
1088 if (.not. present(atom2)) then
1089 do i = 1, nExcludes_global
1090 if (excludesGlobal(i) == unique_id_1) then
1091 skip_it = .true.
1092 return
1093 end if
1094 end do
1095 return
1096 end if
1097
1098 #ifdef IS_MPI
1099 unique_id_2 = AtomColToGlobal(atom2)
1100 #else
1101 unique_id_2 = atom2
1102 #endif
1103
1104 #ifdef IS_MPI
1105 !! this situation should only arise in MPI simulations
1106 if (unique_id_1 == unique_id_2) then
1107 skip_it = .true.
1108 return
1109 end if
1110
1111 !! this prevents us from doing the pair on multiple processors
1112 if (unique_id_1 < unique_id_2) then
1113 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1114 skip_it = .true.
1115 return
1116 endif
1117 else
1118 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1119 skip_it = .true.
1120 return
1121 endif
1122 endif
1123 #endif
1124
1125 !! the rest of these situations can happen in all simulations:
1126 do i = 1, nExcludes_global
1127 if ((excludesGlobal(i) == unique_id_1) .or. &
1128 (excludesGlobal(i) == unique_id_2)) then
1129 skip_it = .true.
1130 return
1131 endif
1132 enddo
1133
1134 do i = 1, nSkipsForAtom(atom1)
1135 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1136 skip_it = .true.
1137 return
1138 endif
1139 end do
1140
1141 return
1142 end function skipThisPair
1143
1144 function FF_UsesDirectionalAtoms() result(doesit)
1145 logical :: doesit
1146 doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1147 FF_uses_GB .or. FF_uses_RF
1148 end function FF_UsesDirectionalAtoms
1149
1150 function FF_RequiresPrepairCalc() result(doesit)
1151 logical :: doesit
1152 doesit = FF_uses_EAM
1153 end function FF_RequiresPrepairCalc
1154
1155 function FF_RequiresPostpairCalc() result(doesit)
1156 logical :: doesit
1157 doesit = FF_uses_RF
1158 end function FF_RequiresPostpairCalc
1159
1160 #ifdef PROFILE
1161 function getforcetime() result(totalforcetime)
1162 real(kind=dp) :: totalforcetime
1163 totalforcetime = forcetime
1164 end function getforcetime
1165 #endif
1166
1167 !! This cleans componets of force arrays belonging only to fortran
1168
1169 subroutine add_stress_tensor(dpair, fpair)
1170
1171 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1172
1173 ! because the d vector is the rj - ri vector, and
1174 ! because fx, fy, fz are the force on atom i, we need a
1175 ! negative sign here:
1176
1177 tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1178 tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1179 tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1180 tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1181 tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1182 tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1183 tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1184 tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1185 tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1186
1187 !write(*,'(6es12.3)') fpair(1:3), tau_Temp(1), tau_Temp(5), tau_temp(9)
1188 virial_Temp = virial_Temp + &
1189 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1190
1191 end subroutine add_stress_tensor
1192
1193 end module do_Forces
1194