ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 669
Committed: Thu Aug 7 00:47:33 2003 UTC (20 years, 11 months ago) by chuckv
File size: 28347 byte(s)
Log Message:
Bug fixes for eam...

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.27 2003-08-07 00:47:33 chuckv Exp $, $Date: 2003-08-07 00:47:33 $, $Name: not supported by cvs2svn $, $Revision: 1.27 $
8
9 module do_Forces
10 use force_globals
11 use simulation
12 use definitions
13 use atype_module
14 use neighborLists
15 use lj
16 use sticky_pair
17 use dipole_dipole
18 use reaction_field
19 use gb_pair
20 use vector_class
21 use eam
22 use status
23 #ifdef IS_MPI
24 use mpiSimulation
25 #endif
26
27 implicit none
28 PRIVATE
29
30 #define __FORTRAN90
31 #include "fForceField.h"
32
33 logical, save :: do_forces_initialized = .false., haveRlist = .false.
34 logical, save :: havePolicies = .false.
35 logical, save :: FF_uses_LJ
36 logical, save :: FF_uses_sticky
37 logical, save :: FF_uses_dipoles
38 logical, save :: FF_uses_RF
39 logical, save :: FF_uses_GB
40 logical, save :: FF_uses_EAM
41
42 real(kind=dp), save :: rlist, rlistsq
43
44 public :: init_FF
45 public :: do_force_loop
46 public :: setRlistDF
47
48 contains
49
50 subroutine setRlistDF( this_rlist )
51
52 real(kind=dp) :: this_rlist
53
54 rlist = this_rlist
55 rlistsq = rlist * rlist
56
57 haveRlist = .true.
58 if( havePolicies ) do_forces_initialized = .true.
59
60 end subroutine setRlistDF
61
62 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
63
64 integer, intent(in) :: LJMIXPOLICY
65 logical, intent(in) :: use_RF_c
66
67 integer, intent(out) :: thisStat
68 integer :: my_status, nMatches
69 integer, pointer :: MatchList(:) => null()
70 real(kind=dp) :: rcut, rrf, rt, dielect
71
72 !! assume things are copacetic, unless they aren't
73 thisStat = 0
74
75 !! Fortran's version of a cast:
76 FF_uses_RF = use_RF_c
77
78 !! init_FF is called *after* all of the atom types have been
79 !! defined in atype_module using the new_atype subroutine.
80 !!
81 !! this will scan through the known atypes and figure out what
82 !! interactions are used by the force field.
83
84 FF_uses_LJ = .false.
85 FF_uses_sticky = .false.
86 FF_uses_dipoles = .false.
87 FF_uses_GB = .false.
88 FF_uses_EAM = .false.
89
90 call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
91 if (nMatches .gt. 0) FF_uses_LJ = .true.
92
93 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
94 if (nMatches .gt. 0) FF_uses_dipoles = .true.
95
96 call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
97 MatchList)
98 if (nMatches .gt. 0) FF_uses_Sticky = .true.
99
100 call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
101 if (nMatches .gt. 0) FF_uses_GB = .true.
102
103 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
104 if (nMatches .gt. 0) FF_uses_EAM = .true.
105
106 !! check to make sure the FF_uses_RF setting makes sense
107
108 if (FF_uses_dipoles) then
109 if (FF_uses_RF) then
110 dielect = getDielect()
111 call initialize_rf(dielect)
112 endif
113 else
114 if (FF_uses_RF) then
115 write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
116 thisStat = -1
117 return
118 endif
119 endif
120
121 if (FF_uses_LJ) then
122
123 select case (LJMIXPOLICY)
124 case (LB_MIXING_RULE)
125 call init_lj_FF(LB_MIXING_RULE, my_status)
126 case (EXPLICIT_MIXING_RULE)
127 call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
128 case default
129 write(default_error,*) 'unknown LJ Mixing Policy!'
130 thisStat = -1
131 return
132 end select
133 if (my_status /= 0) then
134 thisStat = -1
135 return
136 end if
137 endif
138
139 if (FF_uses_sticky) then
140 call check_sticky_FF(my_status)
141 if (my_status /= 0) then
142 thisStat = -1
143 return
144 end if
145 endif
146
147
148 if (FF_uses_EAM) then
149 call init_EAM_FF(my_status)
150 if (my_status /= 0) then
151 thisStat = -1
152 return
153 end if
154 endif
155
156
157
158 if (FF_uses_GB) then
159 call check_gb_pair_FF(my_status)
160 if (my_status .ne. 0) then
161 thisStat = -1
162 return
163 endif
164 endif
165
166 if (FF_uses_GB .and. FF_uses_LJ) then
167 endif
168 if (.not. do_forces_initialized) then
169 !! Create neighbor lists
170 call expandNeighborList(getNlocal(), my_status)
171 if (my_Status /= 0) then
172 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
173 thisStat = -1
174 return
175 endif
176 endif
177
178
179 havePolicies = .true.
180 if( haveRlist ) do_forces_initialized = .true.
181
182 end subroutine init_FF
183
184
185 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
186 !------------------------------------------------------------->
187 subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, &
188 error)
189 !! Position array provided by C, dimensioned by getNlocal
190 real ( kind = dp ), dimension(3,getNlocal()) :: q
191 !! Rotation Matrix for each long range particle in simulation.
192 real( kind = dp), dimension(9,getNlocal()) :: A
193 !! Unit vectors for dipoles (lab frame)
194 real( kind = dp ), dimension(3,getNlocal()) :: u_l
195 !! Force array provided by C, dimensioned by getNlocal
196 real ( kind = dp ), dimension(3,getNlocal()) :: f
197 !! Torsion array provided by C, dimensioned by getNlocal
198 real( kind = dp ), dimension(3,getNlocal()) :: t
199 !! Stress Tensor
200 real( kind = dp), dimension(9) :: tau
201 real ( kind = dp ) :: pot
202 logical ( kind = 2) :: do_pot_c, do_stress_c
203 logical :: do_pot
204 logical :: do_stress
205 #ifdef IS_MPI
206 real( kind = DP ) :: pot_local
207 integer :: nrow
208 integer :: ncol
209 #endif
210 integer :: nlocal
211 integer :: natoms
212 logical :: update_nlist
213 integer :: i, j, jbeg, jend, jnab
214 integer :: nlist
215 real( kind = DP ) :: rijsq
216 real(kind=dp),dimension(3) :: d
217 real(kind=dp) :: rfpot, mu_i, virial
218 integer :: me_i
219 logical :: is_dp_i
220 integer :: neighborListSize
221 integer :: listerror, error
222 integer :: localError
223
224 real(kind=dp) :: listSkin = 1.0
225
226
227 !! initialize local variables
228
229 #ifdef IS_MPI
230 pot_local = 0.0_dp
231 nlocal = getNlocal()
232 nrow = getNrow(plan_row)
233 ncol = getNcol(plan_col)
234 #else
235 nlocal = getNlocal()
236 natoms = nlocal
237 #endif
238
239 call check_initialization(localError)
240 if ( localError .ne. 0 ) then
241 call handleError("do_force_loop","Not Initialized")
242 error = -1
243 return
244 end if
245 call zero_work_arrays()
246
247 do_pot = do_pot_c
248 do_stress = do_stress_c
249
250
251 ! Gather all information needed by all force loops:
252
253 #ifdef IS_MPI
254
255 call gather(q,q_Row,plan_row3d)
256 call gather(q,q_Col,plan_col3d)
257
258 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
259 call gather(u_l,u_l_Row,plan_row3d)
260 call gather(u_l,u_l_Col,plan_col3d)
261
262 call gather(A,A_Row,plan_row_rotation)
263 call gather(A,A_Col,plan_col_rotation)
264 endif
265
266 #endif
267
268 if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
269 !! See if we need to update neighbor lists
270 call checkNeighborList(nlocal, q, listSkin, update_nlist)
271 !! if_mpi_gather_stuff_for_prepair
272 !! do_prepair_loop_if_needed
273 !! if_mpi_scatter_stuff_from_prepair
274 !! if_mpi_gather_stuff_from_prepair_to_main_loop
275
276 !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
277 #ifdef IS_MPI
278
279 if (update_nlist) then
280
281 !! save current configuration, construct neighbor list,
282 !! and calculate forces
283 call saveNeighborList(nlocal, q)
284
285 neighborListSize = size(list)
286 nlist = 0
287
288 do i = 1, nrow
289 point(i) = nlist + 1
290
291 prepair_inner: do j = 1, ncol
292
293 if (skipThisPair(i,j)) cycle prepair_inner
294
295 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
296
297 if (rijsq < rlistsq) then
298
299 nlist = nlist + 1
300
301 if (nlist > neighborListSize) then
302 call expandNeighborList(nlocal, listerror)
303 if (listerror /= 0) then
304 error = -1
305 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
306 return
307 end if
308 neighborListSize = size(list)
309 endif
310
311 list(nlist) = j
312 call do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot_local)
313 endif
314 enddo prepair_inner
315 enddo
316
317 point(nrow + 1) = nlist + 1
318
319 else !! (of update_check)
320
321 ! use the list to find the neighbors
322 do i = 1, nrow
323 JBEG = POINT(i)
324 JEND = POINT(i+1) - 1
325 ! check thiat molecule i has neighbors
326 if (jbeg .le. jend) then
327
328 do jnab = jbeg, jend
329 j = list(jnab)
330
331 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
332 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
333 u_l, A, f, t, pot_local)
334
335 enddo
336 endif
337 enddo
338 endif
339
340 #else
341
342 if (update_nlist) then
343
344 ! save current configuration, contruct neighbor list,
345 ! and calculate forces
346 call saveNeighborList(natoms, q)
347
348 neighborListSize = size(list)
349
350 nlist = 0
351
352 do i = 1, natoms-1
353 point(i) = nlist + 1
354
355 prepair_inner: do j = i+1, natoms
356
357 if (skipThisPair(i,j)) cycle prepair_inner
358
359 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
360
361
362 if (rijsq < rlistsq) then
363
364 nlist = nlist + 1
365
366 if (nlist > neighborListSize) then
367 call expandNeighborList(natoms, listerror)
368 if (listerror /= 0) then
369 error = -1
370 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
371 return
372 end if
373 neighborListSize = size(list)
374 endif
375
376 list(nlist) = j
377
378 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
379 u_l, A, f, t, pot)
380
381 endif
382 enddo prepair_inner
383 enddo
384
385 point(natoms) = nlist + 1
386
387 else !! (update)
388
389 ! use the list to find the neighbors
390 do i = 1, natoms-1
391 JBEG = POINT(i)
392 JEND = POINT(i+1) - 1
393 ! check thiat molecule i has neighbors
394 if (jbeg .le. jend) then
395
396 do jnab = jbeg, jend
397 j = list(jnab)
398
399 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
400 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
401 u_l, A, f, t, pot)
402
403 enddo
404 endif
405 enddo
406 endif
407 #endif
408 !! Do rest of preforce calculations
409 call do_preforce(nlocal,pot)
410 else
411 !! See if we need to update neighbor lists for non pre-pair
412 call checkNeighborList(nlocal, q, listSkin, update_nlist)
413 endif
414
415
416
417
418
419 !---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>>>>>>>>>>>>>>>>
420
421
422
423
424
425 #ifdef IS_MPI
426
427 if (update_nlist) then
428
429 !! save current configuration, construct neighbor list,
430 !! and calculate forces
431 call saveNeighborList(nlocal, q)
432
433 neighborListSize = size(list)
434 nlist = 0
435
436 do i = 1, nrow
437 point(i) = nlist + 1
438
439 inner: do j = 1, ncol
440
441 if (skipThisPair(i,j)) cycle inner
442
443 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
444
445 if (rijsq < rlistsq) then
446
447 nlist = nlist + 1
448
449 if (nlist > neighborListSize) then
450 call expandNeighborList(nlocal, listerror)
451 if (listerror /= 0) then
452 error = -1
453 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
454 return
455 end if
456 neighborListSize = size(list)
457 endif
458
459 list(nlist) = j
460
461 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
462 u_l, A, f, t, pot_local)
463
464 endif
465 enddo inner
466 enddo
467
468 point(nrow + 1) = nlist + 1
469
470 else !! (of update_check)
471
472 ! use the list to find the neighbors
473 do i = 1, nrow
474 JBEG = POINT(i)
475 JEND = POINT(i+1) - 1
476 ! check thiat molecule i has neighbors
477 if (jbeg .le. jend) then
478
479 do jnab = jbeg, jend
480 j = list(jnab)
481
482 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
483 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
484 u_l, A, f, t, pot_local)
485
486 enddo
487 endif
488 enddo
489 endif
490
491 #else
492
493 if (update_nlist) then
494
495 ! save current configuration, contruct neighbor list,
496 ! and calculate forces
497 call saveNeighborList(natoms, q)
498
499 neighborListSize = size(list)
500
501 nlist = 0
502
503 do i = 1, natoms-1
504 point(i) = nlist + 1
505
506 inner: do j = i+1, natoms
507
508 if (skipThisPair(i,j)) cycle inner
509
510 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
511
512
513 if (rijsq < rlistsq) then
514
515 nlist = nlist + 1
516
517 if (nlist > neighborListSize) then
518 call expandNeighborList(natoms, listerror)
519 if (listerror /= 0) then
520 error = -1
521 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
522 return
523 end if
524 neighborListSize = size(list)
525 endif
526
527 list(nlist) = j
528
529 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
530 u_l, A, f, t, pot)
531
532 endif
533 enddo inner
534 enddo
535
536 point(natoms) = nlist + 1
537
538 else !! (update)
539
540 ! use the list to find the neighbors
541 do i = 1, natoms-1
542 JBEG = POINT(i)
543 JEND = POINT(i+1) - 1
544 ! check thiat molecule i has neighbors
545 if (jbeg .le. jend) then
546
547 do jnab = jbeg, jend
548 j = list(jnab)
549
550 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
551 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
552 u_l, A, f, t, pot)
553
554 enddo
555 endif
556 enddo
557 endif
558
559 #endif
560
561 ! phew, done with main loop.
562
563 #ifdef IS_MPI
564 !!distribute forces
565
566 f_temp = 0.0_dp
567 call scatter(f_Row,f_temp,plan_row3d)
568 do i = 1,nlocal
569 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
570 end do
571
572 f_temp = 0.0_dp
573 call scatter(f_Col,f_temp,plan_col3d)
574 do i = 1,nlocal
575 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
576 end do
577
578 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
579 t_temp = 0.0_dp
580 call scatter(t_Row,t_temp,plan_row3d)
581 do i = 1,nlocal
582 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
583 end do
584 t_temp = 0.0_dp
585 call scatter(t_Col,t_temp,plan_col3d)
586
587 do i = 1,nlocal
588 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
589 end do
590 endif
591
592 if (do_pot) then
593 ! scatter/gather pot_row into the members of my column
594 call scatter(pot_Row, pot_Temp, plan_row)
595
596 ! scatter/gather pot_local into all other procs
597 ! add resultant to get total pot
598 do i = 1, nlocal
599 pot_local = pot_local + pot_Temp(i)
600 enddo
601
602 pot_Temp = 0.0_DP
603
604 call scatter(pot_Col, pot_Temp, plan_col)
605 do i = 1, nlocal
606 pot_local = pot_local + pot_Temp(i)
607 enddo
608
609 endif
610 #endif
611
612 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
613
614 if (FF_uses_RF .and. SimUsesRF()) then
615
616 #ifdef IS_MPI
617 call scatter(rf_Row,rf,plan_row3d)
618 call scatter(rf_Col,rf_Temp,plan_col3d)
619 do i = 1,nlocal
620 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
621 end do
622 #endif
623
624 do i = 1, getNlocal()
625
626 rfpot = 0.0_DP
627 #ifdef IS_MPI
628 me_i = atid_row(i)
629 #else
630 me_i = atid(i)
631 #endif
632 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
633 if ( is_DP_i ) then
634 call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
635 !! The reaction field needs to include a self contribution
636 !! to the field:
637 call accumulate_self_rf(i, mu_i, u_l)
638 !! Get the reaction field contribution to the
639 !! potential and torques:
640 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
641 #ifdef IS_MPI
642 pot_local = pot_local + rfpot
643 #else
644 pot = pot + rfpot
645
646 #endif
647 endif
648 enddo
649 endif
650 endif
651
652
653 #ifdef IS_MPI
654
655 if (do_pot) then
656 pot = pot + pot_local
657 !! we assume the c code will do the allreduce to get the total potential
658 !! we could do it right here if we needed to...
659 endif
660
661 if (do_stress) then
662 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
663 mpi_comm_world,mpi_err)
664 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
665 mpi_comm_world,mpi_err)
666 endif
667
668 #else
669
670 if (do_stress) then
671 tau = tau_Temp
672 virial = virial_Temp
673 endif
674
675 #endif
676
677 end subroutine do_force_loop
678
679 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
680
681 real( kind = dp ) :: pot
682 real( kind = dp ), dimension(3,getNlocal()) :: u_l
683 real (kind=dp), dimension(9,getNlocal()) :: A
684 real (kind=dp), dimension(3,getNlocal()) :: f
685 real (kind=dp), dimension(3,getNlocal()) :: t
686
687 logical, intent(inout) :: do_pot, do_stress
688 integer, intent(in) :: i, j
689 real ( kind = dp ), intent(inout) :: rijsq
690 real ( kind = dp ) :: r
691 real ( kind = dp ), intent(inout) :: d(3)
692 logical :: is_LJ_i, is_LJ_j
693 logical :: is_DP_i, is_DP_j
694 logical :: is_GB_i, is_GB_j
695 logical :: is_EAM_i,is_EAM_j
696 logical :: is_Sticky_i, is_Sticky_j
697 integer :: me_i, me_j
698
699 r = sqrt(rijsq)
700
701 #ifdef IS_MPI
702 if (tagRow(i) .eq. tagColumn(j)) then
703 write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
704 endif
705
706 me_i = atid_row(i)
707 me_j = atid_col(j)
708
709 #else
710
711 me_i = atid(i)
712 me_j = atid(j)
713
714 #endif
715
716 if (FF_uses_LJ .and. SimUsesLJ()) then
717 call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
718 call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
719
720 if ( is_LJ_i .and. is_LJ_j ) &
721 call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
722 endif
723
724 if (FF_uses_dipoles .and. SimUsesDipoles()) then
725 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
726 call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
727
728 if ( is_DP_i .and. is_DP_j ) then
729
730 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
731 do_pot, do_stress)
732 if (FF_uses_RF .and. SimUsesRF()) then
733 call accumulate_rf(i, j, r, u_l)
734 call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
735 endif
736
737 endif
738 endif
739
740 if (FF_uses_Sticky .and. SimUsesSticky()) then
741
742 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
743 call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
744
745 if ( is_Sticky_i .and. is_Sticky_j ) then
746 call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
747 do_pot, do_stress)
748 endif
749 endif
750
751
752 if (FF_uses_GB .and. SimUsesGB()) then
753
754 call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
755 call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
756
757 if ( is_GB_i .and. is_GB_j ) then
758 call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
759 do_pot, do_stress)
760 endif
761 endif
762
763
764
765 if (FF_uses_EAM .and. SimUsesEAM()) then
766 call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
767 call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
768
769 if ( is_EAM_i .and. is_EAM_j ) &
770 call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
771 endif
772
773
774
775
776 end subroutine do_pair
777
778
779
780 subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
781 real( kind = dp ) :: pot
782 real( kind = dp ), dimension(3,getNlocal()) :: u_l
783 real (kind=dp), dimension(9,getNlocal()) :: A
784 real (kind=dp), dimension(3,getNlocal()) :: f
785 real (kind=dp), dimension(3,getNlocal()) :: t
786
787 logical, intent(inout) :: do_pot, do_stress
788 integer, intent(in) :: i, j
789 real ( kind = dp ), intent(inout) :: rijsq
790 real ( kind = dp ) :: r
791 real ( kind = dp ), intent(inout) :: d(3)
792
793 logical :: is_EAM_i, is_EAM_j
794
795 integer :: me_i, me_j
796
797 r = sqrt(rijsq)
798
799
800 #ifdef IS_MPI
801 if (tagRow(i) .eq. tagColumn(j)) then
802 write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
803 endif
804
805 me_i = atid_row(i)
806 me_j = atid_col(j)
807
808 #else
809
810 me_i = atid(i)
811 me_j = atid(j)
812
813 #endif
814
815 if (FF_uses_EAM .and. SimUsesEAM()) then
816 call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
817 call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
818
819 if ( is_EAM_i .and. is_EAM_j ) &
820 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
821 endif
822 end subroutine do_prepair
823
824
825
826
827 subroutine do_preforce(nlocal,pot)
828 integer :: nlocal
829 real( kind = dp ) :: pot
830
831 if (FF_uses_EAM .and. SimUsesEAM()) then
832 call calc_EAM_preforce_Frho(nlocal,pot)
833 endif
834
835
836 end subroutine do_preforce
837
838
839 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
840
841 real (kind = dp), dimension(3) :: q_i
842 real (kind = dp), dimension(3) :: q_j
843 real ( kind = dp ), intent(out) :: r_sq
844 real( kind = dp ) :: d(3), scaled(3)
845 integer i
846
847 d(1:3) = q_j(1:3) - q_i(1:3)
848
849 ! Wrap back into periodic box if necessary
850 if ( SimUsesPBC() ) then
851
852 if( .not.boxIsOrthorhombic ) then
853 ! calc the scaled coordinates.
854
855 scaled = matmul(HmatInv, d)
856
857 ! wrap the scaled coordinates
858
859 scaled = scaled - anint(scaled)
860
861
862 ! calc the wrapped real coordinates from the wrapped scaled
863 ! coordinates
864
865 d = matmul(Hmat,scaled)
866
867 else
868 ! calc the scaled coordinates.
869
870 do i = 1, 3
871 scaled(i) = d(i) * HmatInv(i,i)
872
873 ! wrap the scaled coordinates
874
875 scaled(i) = scaled(i) - anint(scaled(i))
876
877 ! calc the wrapped real coordinates from the wrapped scaled
878 ! coordinates
879
880 d(i) = scaled(i)*Hmat(i,i)
881 enddo
882 endif
883
884 endif
885
886 r_sq = dot_product(d,d)
887
888 end subroutine get_interatomic_vector
889
890 subroutine check_initialization(error)
891 integer, intent(out) :: error
892
893 error = 0
894 ! Make sure we are properly initialized.
895 if (.not. do_forces_initialized) then
896 write(*,*) "Forces not initialized"
897 error = -1
898 return
899 endif
900
901 #ifdef IS_MPI
902 if (.not. isMPISimSet()) then
903 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
904 error = -1
905 return
906 endif
907 #endif
908
909 return
910 end subroutine check_initialization
911
912
913 subroutine zero_work_arrays()
914
915 #ifdef IS_MPI
916
917 q_Row = 0.0_dp
918 q_Col = 0.0_dp
919
920 u_l_Row = 0.0_dp
921 u_l_Col = 0.0_dp
922
923 A_Row = 0.0_dp
924 A_Col = 0.0_dp
925
926 f_Row = 0.0_dp
927 f_Col = 0.0_dp
928 f_Temp = 0.0_dp
929
930 t_Row = 0.0_dp
931 t_Col = 0.0_dp
932 t_Temp = 0.0_dp
933
934 pot_Row = 0.0_dp
935 pot_Col = 0.0_dp
936 pot_Temp = 0.0_dp
937
938 rf_Row = 0.0_dp
939 rf_Col = 0.0_dp
940 rf_Temp = 0.0_dp
941
942 #endif
943
944 rf = 0.0_dp
945 tau_Temp = 0.0_dp
946 virial_Temp = 0.0_dp
947 end subroutine zero_work_arrays
948
949 function skipThisPair(atom1, atom2) result(skip_it)
950 integer, intent(in) :: atom1
951 integer, intent(in), optional :: atom2
952 logical :: skip_it
953 integer :: unique_id_1, unique_id_2
954 integer :: me_i,me_j
955 integer :: i
956
957 skip_it = .false.
958
959 !! there are a number of reasons to skip a pair or a particle
960 !! mostly we do this to exclude atoms who are involved in short
961 !! range interactions (bonds, bends, torsions), but we also need
962 !! to exclude some overcounted interactions that result from
963 !! the parallel decomposition
964
965 #ifdef IS_MPI
966 !! in MPI, we have to look up the unique IDs for each atom
967 unique_id_1 = tagRow(atom1)
968 #else
969 !! in the normal loop, the atom numbers are unique
970 unique_id_1 = atom1
971 #endif
972
973 !! We were called with only one atom, so just check the global exclude
974 !! list for this atom
975 if (.not. present(atom2)) then
976 do i = 1, nExcludes_global
977 if (excludesGlobal(i) == unique_id_1) then
978 skip_it = .true.
979 return
980 end if
981 end do
982 return
983 end if
984
985 #ifdef IS_MPI
986 unique_id_2 = tagColumn(atom2)
987 #else
988 unique_id_2 = atom2
989 #endif
990
991 #ifdef IS_MPI
992 !! this situation should only arise in MPI simulations
993 if (unique_id_1 == unique_id_2) then
994 skip_it = .true.
995 return
996 end if
997
998 !! this prevents us from doing the pair on multiple processors
999 if (unique_id_1 < unique_id_2) then
1000 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1001 skip_it = .true.
1002 return
1003 endif
1004 else
1005 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1006 skip_it = .true.
1007 return
1008 endif
1009 endif
1010 #endif
1011
1012 !! the rest of these situations can happen in all simulations:
1013 do i = 1, nExcludes_global
1014 if ((excludesGlobal(i) == unique_id_1) .or. &
1015 (excludesGlobal(i) == unique_id_2)) then
1016 skip_it = .true.
1017 return
1018 endif
1019 enddo
1020
1021 do i = 1, nExcludes_local
1022 if (excludesLocal(1,i) == unique_id_1) then
1023 if (excludesLocal(2,i) == unique_id_2) then
1024 skip_it = .true.
1025 return
1026 endif
1027 else
1028 if (excludesLocal(1,i) == unique_id_2) then
1029 if (excludesLocal(2,i) == unique_id_1) then
1030 skip_it = .true.
1031 return
1032 endif
1033 endif
1034 endif
1035 end do
1036
1037 return
1038 end function skipThisPair
1039
1040 function FF_UsesDirectionalAtoms() result(doesit)
1041 logical :: doesit
1042 doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1043 FF_uses_GB .or. FF_uses_RF
1044 end function FF_UsesDirectionalAtoms
1045
1046 function FF_RequiresPrepairCalc() result(doesit)
1047 logical :: doesit
1048 doesit = FF_uses_EAM
1049 end function FF_RequiresPrepairCalc
1050
1051 function FF_RequiresPostpairCalc() result(doesit)
1052 logical :: doesit
1053 doesit = FF_uses_RF
1054 end function FF_RequiresPostpairCalc
1055
1056 end module do_Forces