ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 891
Committed: Fri Dec 19 20:36:35 2003 UTC (20 years, 6 months ago) by mmeineke
File size: 29249 byte(s)
Log Message:
More profiling fixes.

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