ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 883
Committed: Thu Dec 18 20:46:45 2003 UTC (20 years, 6 months ago) by chuckv
File size: 29234 byte(s)
Log Message:
Added functions for simple profiling in fortran.

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