ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 872
Committed: Fri Nov 21 19:31:05 2003 UTC (20 years, 7 months ago) by chrisfen
File size: 29955 byte(s)
Log Message:
Fixed a bug in SimInfo ordering of radii

File Contents

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