ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 841
Committed: Wed Oct 29 17:55:28 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 30145 byte(s)
Log Message:
som efixes to the way rcut is setup, as well as additional debugging comments.

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