ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1197
Committed: Wed May 26 16:41:23 2004 UTC (21 years, 1 month ago) by gezelter
File size: 33748 byte(s)
Log Message:
Compacted all of the 8 copies of the force loop into one.

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