ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 694
Committed: Wed Aug 13 21:20:20 2003 UTC (20 years, 10 months ago) by chuckv
File size: 29923 byte(s)
Log Message:
Added some profiling code -DPROFILE.

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 chuckv 694 !! @version $Id: do_Forces.F90,v 1.29 2003-08-13 21:20:20 chuckv Exp $, $Date: 2003-08-13 21:20:20 $, $Name: not supported by cvs2svn $, $Revision: 1.29 $
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     case (LB_MIXING_RULE)
133 mmeineke 626 call init_lj_FF(LB_MIXING_RULE, my_status)
134 mmeineke 377 case (EXPLICIT_MIXING_RULE)
135 mmeineke 626 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     call init_EAM_FF(my_status)
158     if (my_status /= 0) then
159     thisStat = -1
160     return
161     end if
162     endif
163    
164    
165 mmeineke 377
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 chuckv 480 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 chuckv 657
186 mmeineke 377
187 mmeineke 626 havePolicies = .true.
188     if( haveRlist ) do_forces_initialized = .true.
189 chuckv 657
190 mmeineke 377 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 chuckv 439 #ifdef IS_MPI
214 chuckv 441 real( kind = DP ) :: pot_local
215 mmeineke 377 integer :: nrow
216     integer :: ncol
217 chuckv 694 integer :: nprocs
218 mmeineke 377 #endif
219     integer :: nlocal
220     integer :: natoms
221     logical :: update_nlist
222     integer :: i, j, jbeg, jend, jnab
223     integer :: nlist
224 mmeineke 626 real( kind = DP ) :: rijsq
225 mmeineke 377 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 mmeineke 626
233     real(kind=dp) :: listSkin = 1.0
234 mmeineke 619
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    
448     !! save current configuration, construct neighbor list,
449     !! and calculate forces
450 mmeineke 459 call saveNeighborList(nlocal, q)
451 mmeineke 377
452     neighborListSize = size(list)
453     nlist = 0
454    
455     do i = 1, nrow
456     point(i) = nlist + 1
457    
458     inner: do j = 1, ncol
459    
460     if (skipThisPair(i,j)) cycle inner
461    
462     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
463    
464 mmeineke 626 if (rijsq < rlistsq) then
465 mmeineke 377
466     nlist = nlist + 1
467    
468     if (nlist > neighborListSize) then
469     call expandNeighborList(nlocal, listerror)
470     if (listerror /= 0) then
471     error = -1
472     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
473     return
474     end if
475     neighborListSize = size(list)
476     endif
477    
478     list(nlist) = j
479    
480 mmeineke 626 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
481     u_l, A, f, t, pot_local)
482    
483 mmeineke 377 endif
484     enddo inner
485     enddo
486    
487     point(nrow + 1) = nlist + 1
488    
489     else !! (of update_check)
490    
491     ! use the list to find the neighbors
492     do i = 1, nrow
493     JBEG = POINT(i)
494     JEND = POINT(i+1) - 1
495     ! check thiat molecule i has neighbors
496     if (jbeg .le. jend) then
497    
498     do jnab = jbeg, jend
499     j = list(jnab)
500    
501     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
502     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
503 chuckv 441 u_l, A, f, t, pot_local)
504 mmeineke 377
505     enddo
506     endif
507     enddo
508     endif
509    
510     #else
511    
512     if (update_nlist) then
513    
514     ! save current configuration, contruct neighbor list,
515     ! and calculate forces
516 mmeineke 459 call saveNeighborList(natoms, q)
517 mmeineke 377
518     neighborListSize = size(list)
519    
520     nlist = 0
521    
522     do i = 1, natoms-1
523     point(i) = nlist + 1
524    
525     inner: do j = i+1, natoms
526    
527 chuckv 388 if (skipThisPair(i,j)) cycle inner
528    
529 mmeineke 377 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
530    
531 chuckv 388
532 mmeineke 626 if (rijsq < rlistsq) then
533 mmeineke 377
534     nlist = nlist + 1
535    
536     if (nlist > neighborListSize) then
537     call expandNeighborList(natoms, listerror)
538     if (listerror /= 0) then
539     error = -1
540     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
541     return
542     end if
543     neighborListSize = size(list)
544     endif
545    
546     list(nlist) = j
547    
548 mmeineke 626 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
549 chuckv 441 u_l, A, f, t, pot)
550 mmeineke 626
551 mmeineke 377 endif
552     enddo inner
553     enddo
554    
555     point(natoms) = nlist + 1
556    
557     else !! (update)
558    
559     ! use the list to find the neighbors
560     do i = 1, natoms-1
561     JBEG = POINT(i)
562     JEND = POINT(i+1) - 1
563     ! check thiat molecule i has neighbors
564     if (jbeg .le. jend) then
565    
566     do jnab = jbeg, jend
567     j = list(jnab)
568    
569     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
570     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
571 chuckv 441 u_l, A, f, t, pot)
572 mmeineke 377
573     enddo
574     endif
575     enddo
576     endif
577    
578     #endif
579    
580     ! phew, done with main loop.
581 chuckv 694
582     !! Do timing
583     #ifdef PROFILE
584     call cpu_time(forceTimeFinal)
585     forceTime = forceTime + forceTimeFinal - forceTimeInitial
586     #endif
587    
588    
589 mmeineke 377 #ifdef IS_MPI
590     !!distribute forces
591 chuckv 438
592     f_temp = 0.0_dp
593     call scatter(f_Row,f_temp,plan_row3d)
594     do i = 1,nlocal
595     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
596     end do
597    
598     f_temp = 0.0_dp
599 mmeineke 377 call scatter(f_Col,f_temp,plan_col3d)
600     do i = 1,nlocal
601     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
602     end do
603    
604     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
605 chuckv 438 t_temp = 0.0_dp
606     call scatter(t_Row,t_temp,plan_row3d)
607     do i = 1,nlocal
608     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
609     end do
610     t_temp = 0.0_dp
611 mmeineke 377 call scatter(t_Col,t_temp,plan_col3d)
612    
613     do i = 1,nlocal
614     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
615     end do
616     endif
617    
618     if (do_pot) then
619     ! scatter/gather pot_row into the members of my column
620     call scatter(pot_Row, pot_Temp, plan_row)
621 chuckv 439
622 mmeineke 377 ! scatter/gather pot_local into all other procs
623     ! add resultant to get total pot
624     do i = 1, nlocal
625     pot_local = pot_local + pot_Temp(i)
626     enddo
627 chuckv 439
628     pot_Temp = 0.0_DP
629 mmeineke 377
630     call scatter(pot_Col, pot_Temp, plan_col)
631     do i = 1, nlocal
632     pot_local = pot_local + pot_Temp(i)
633     enddo
634 chuckv 439
635 mmeineke 377 endif
636     #endif
637    
638     if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
639    
640     if (FF_uses_RF .and. SimUsesRF()) then
641    
642     #ifdef IS_MPI
643     call scatter(rf_Row,rf,plan_row3d)
644     call scatter(rf_Col,rf_Temp,plan_col3d)
645     do i = 1,nlocal
646     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
647     end do
648     #endif
649    
650     do i = 1, getNlocal()
651    
652     rfpot = 0.0_DP
653     #ifdef IS_MPI
654     me_i = atid_row(i)
655     #else
656     me_i = atid(i)
657     #endif
658     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
659     if ( is_DP_i ) then
660     call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
661     !! The reaction field needs to include a self contribution
662     !! to the field:
663     call accumulate_self_rf(i, mu_i, u_l)
664     !! Get the reaction field contribution to the
665     !! potential and torques:
666     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
667     #ifdef IS_MPI
668     pot_local = pot_local + rfpot
669     #else
670     pot = pot + rfpot
671    
672     #endif
673     endif
674     enddo
675     endif
676     endif
677    
678    
679     #ifdef IS_MPI
680    
681     if (do_pot) then
682 chuckv 441 pot = pot + pot_local
683 mmeineke 377 !! we assume the c code will do the allreduce to get the total potential
684     !! we could do it right here if we needed to...
685     endif
686    
687     if (do_stress) then
688 gezelter 490 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
689 mmeineke 377 mpi_comm_world,mpi_err)
690 chuckv 470 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
691 mmeineke 377 mpi_comm_world,mpi_err)
692     endif
693    
694     #else
695    
696     if (do_stress) then
697     tau = tau_Temp
698     virial = virial_Temp
699     endif
700    
701     #endif
702 chuckv 694
703     #ifdef PROFILE
704     if (do_pot) then
705    
706     #ifdef IS_MPI
707    
708    
709     call printCommTime()
710    
711     call mpi_allreduce(forceTime,globalForceTime,1,MPI_DOUBLE_PRECISION, &
712     mpi_sum,mpi_comm_world,mpi_err)
713    
714     call mpi_allreduce(forceTime,maxForceTime,1,MPI_DOUBLE_PRECISION, &
715     MPI_MAX,mpi_comm_world,mpi_err)
716    
717     call mpi_comm_size( MPI_COMM_WORLD, nprocs,mpi_err)
718    
719     if (getMyNode() == 0) then
720     write(*,*) "Total processor time spent in force calculations is: ", globalForceTime
721     write(*,*) "Total Time spent in force loop per processor is: ", globalforceTime/nprocs
722     write(*,*) "Maximum force time on any processor is: ", maxForceTime
723     end if
724     #else
725     write(*,*) "Time spent in force loop is: ", forceTime
726     #endif
727    
728    
729     endif
730    
731     #endif
732    
733    
734    
735 mmeineke 377 end subroutine do_force_loop
736    
737 chuckv 441 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
738 mmeineke 377
739     real( kind = dp ) :: pot
740 chuckv 460 real( kind = dp ), dimension(3,getNlocal()) :: u_l
741     real (kind=dp), dimension(9,getNlocal()) :: A
742     real (kind=dp), dimension(3,getNlocal()) :: f
743     real (kind=dp), dimension(3,getNlocal()) :: t
744 mmeineke 377
745     logical, intent(inout) :: do_pot, do_stress
746     integer, intent(in) :: i, j
747     real ( kind = dp ), intent(inout) :: rijsq
748     real ( kind = dp ) :: r
749     real ( kind = dp ), intent(inout) :: d(3)
750     logical :: is_LJ_i, is_LJ_j
751     logical :: is_DP_i, is_DP_j
752     logical :: is_GB_i, is_GB_j
753 chuckv 648 logical :: is_EAM_i,is_EAM_j
754 mmeineke 377 logical :: is_Sticky_i, is_Sticky_j
755     integer :: me_i, me_j
756    
757     r = sqrt(rijsq)
758    
759     #ifdef IS_MPI
760 gezelter 490 if (tagRow(i) .eq. tagColumn(j)) then
761     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
762     endif
763 mmeineke 377
764     me_i = atid_row(i)
765     me_j = atid_col(j)
766    
767     #else
768    
769     me_i = atid(i)
770     me_j = atid(j)
771    
772     #endif
773    
774     if (FF_uses_LJ .and. SimUsesLJ()) then
775     call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
776     call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
777    
778     if ( is_LJ_i .and. is_LJ_j ) &
779     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
780     endif
781    
782     if (FF_uses_dipoles .and. SimUsesDipoles()) then
783     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
784     call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
785    
786     if ( is_DP_i .and. is_DP_j ) then
787    
788 gezelter 462 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
789 mmeineke 377 do_pot, do_stress)
790     if (FF_uses_RF .and. SimUsesRF()) then
791     call accumulate_rf(i, j, r, u_l)
792     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
793     endif
794    
795     endif
796     endif
797    
798     if (FF_uses_Sticky .and. SimUsesSticky()) then
799    
800     call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
801     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
802 chuckv 388
803 mmeineke 377 if ( is_Sticky_i .and. is_Sticky_j ) then
804     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
805     do_pot, do_stress)
806     endif
807     endif
808    
809    
810     if (FF_uses_GB .and. SimUsesGB()) then
811    
812     call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
813     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
814    
815     if ( is_GB_i .and. is_GB_j ) then
816     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
817     do_pot, do_stress)
818     endif
819     endif
820    
821 mmeineke 597
822 chuckv 648
823     if (FF_uses_EAM .and. SimUsesEAM()) then
824     call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
825     call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
826    
827     if ( is_EAM_i .and. is_EAM_j ) &
828     call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
829     endif
830    
831 mmeineke 597
832 chuckv 648
833    
834 mmeineke 377 end subroutine do_pair
835    
836    
837 chuckv 631
838 chuckv 648 subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
839 chuckv 631 real( kind = dp ) :: pot
840     real( kind = dp ), dimension(3,getNlocal()) :: u_l
841     real (kind=dp), dimension(9,getNlocal()) :: A
842     real (kind=dp), dimension(3,getNlocal()) :: f
843     real (kind=dp), dimension(3,getNlocal()) :: t
844    
845     logical, intent(inout) :: do_pot, do_stress
846     integer, intent(in) :: i, j
847     real ( kind = dp ), intent(inout) :: rijsq
848     real ( kind = dp ) :: r
849     real ( kind = dp ), intent(inout) :: d(3)
850    
851     logical :: is_EAM_i, is_EAM_j
852    
853     integer :: me_i, me_j
854    
855     r = sqrt(rijsq)
856    
857 chuckv 669
858 chuckv 631 #ifdef IS_MPI
859     if (tagRow(i) .eq. tagColumn(j)) then
860     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
861     endif
862    
863     me_i = atid_row(i)
864     me_j = atid_col(j)
865    
866     #else
867    
868     me_i = atid(i)
869     me_j = atid(j)
870    
871     #endif
872 chuckv 673
873 chuckv 631 if (FF_uses_EAM .and. SimUsesEAM()) then
874     call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
875     call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
876    
877 chuckv 648 if ( is_EAM_i .and. is_EAM_j ) &
878     call calc_EAM_prepair_rho(i, j, d, r, rijsq )
879 chuckv 631 endif
880 chuckv 648
881 chuckv 673 end subroutine do_prepair
882 chuckv 648
883    
884    
885 chuckv 673
886 chuckv 648 subroutine do_preforce(nlocal,pot)
887     integer :: nlocal
888     real( kind = dp ) :: pot
889    
890 chuckv 669 if (FF_uses_EAM .and. SimUsesEAM()) then
891     call calc_EAM_preforce_Frho(nlocal,pot)
892     endif
893 chuckv 648
894    
895 chuckv 631 end subroutine do_preforce
896    
897    
898 mmeineke 377 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
899    
900     real (kind = dp), dimension(3) :: q_i
901     real (kind = dp), dimension(3) :: q_j
902     real ( kind = dp ), intent(out) :: r_sq
903 gezelter 571 real( kind = dp ) :: d(3), scaled(3)
904     integer i
905    
906 chuckv 482 d(1:3) = q_j(1:3) - q_i(1:3)
907 gezelter 571
908 mmeineke 377 ! Wrap back into periodic box if necessary
909     if ( SimUsesPBC() ) then
910 mmeineke 393
911 gezelter 571 if( .not.boxIsOrthorhombic ) then
912     ! calc the scaled coordinates.
913    
914 mmeineke 572 scaled = matmul(HmatInv, d)
915 gezelter 571
916     ! wrap the scaled coordinates
917    
918 mmeineke 572 scaled = scaled - anint(scaled)
919    
920 gezelter 571
921     ! calc the wrapped real coordinates from the wrapped scaled
922     ! coordinates
923    
924 mmeineke 572 d = matmul(Hmat,scaled)
925 gezelter 571
926     else
927     ! calc the scaled coordinates.
928    
929     do i = 1, 3
930     scaled(i) = d(i) * HmatInv(i,i)
931    
932     ! wrap the scaled coordinates
933    
934     scaled(i) = scaled(i) - anint(scaled(i))
935    
936     ! calc the wrapped real coordinates from the wrapped scaled
937     ! coordinates
938    
939     d(i) = scaled(i)*Hmat(i,i)
940     enddo
941     endif
942 mmeineke 393
943 mmeineke 377 endif
944 gezelter 571
945 mmeineke 377 r_sq = dot_product(d,d)
946 gezelter 571
947 mmeineke 377 end subroutine get_interatomic_vector
948 gezelter 571
949 mmeineke 377 subroutine check_initialization(error)
950     integer, intent(out) :: error
951    
952     error = 0
953     ! Make sure we are properly initialized.
954     if (.not. do_forces_initialized) then
955 chuckv 657 write(*,*) "Forces not initialized"
956 mmeineke 377 error = -1
957     return
958     endif
959    
960     #ifdef IS_MPI
961     if (.not. isMPISimSet()) then
962     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
963     error = -1
964     return
965     endif
966     #endif
967    
968     return
969     end subroutine check_initialization
970    
971    
972     subroutine zero_work_arrays()
973    
974     #ifdef IS_MPI
975    
976     q_Row = 0.0_dp
977     q_Col = 0.0_dp
978    
979     u_l_Row = 0.0_dp
980     u_l_Col = 0.0_dp
981    
982     A_Row = 0.0_dp
983     A_Col = 0.0_dp
984    
985     f_Row = 0.0_dp
986     f_Col = 0.0_dp
987     f_Temp = 0.0_dp
988    
989     t_Row = 0.0_dp
990     t_Col = 0.0_dp
991     t_Temp = 0.0_dp
992    
993     pot_Row = 0.0_dp
994     pot_Col = 0.0_dp
995     pot_Temp = 0.0_dp
996    
997     rf_Row = 0.0_dp
998     rf_Col = 0.0_dp
999     rf_Temp = 0.0_dp
1000    
1001     #endif
1002    
1003 chuckv 673
1004     if (FF_uses_EAM .and. SimUsesEAM()) then
1005     call clean_EAM()
1006     endif
1007    
1008    
1009    
1010    
1011    
1012 mmeineke 377 rf = 0.0_dp
1013     tau_Temp = 0.0_dp
1014     virial_Temp = 0.0_dp
1015     end subroutine zero_work_arrays
1016    
1017     function skipThisPair(atom1, atom2) result(skip_it)
1018     integer, intent(in) :: atom1
1019     integer, intent(in), optional :: atom2
1020     logical :: skip_it
1021     integer :: unique_id_1, unique_id_2
1022 chuckv 388 integer :: me_i,me_j
1023 mmeineke 377 integer :: i
1024    
1025     skip_it = .false.
1026    
1027     !! there are a number of reasons to skip a pair or a particle
1028     !! mostly we do this to exclude atoms who are involved in short
1029     !! range interactions (bonds, bends, torsions), but we also need
1030     !! to exclude some overcounted interactions that result from
1031     !! the parallel decomposition
1032    
1033     #ifdef IS_MPI
1034     !! in MPI, we have to look up the unique IDs for each atom
1035     unique_id_1 = tagRow(atom1)
1036     #else
1037     !! in the normal loop, the atom numbers are unique
1038     unique_id_1 = atom1
1039     #endif
1040 chuckv 388
1041 mmeineke 377 !! We were called with only one atom, so just check the global exclude
1042     !! list for this atom
1043     if (.not. present(atom2)) then
1044     do i = 1, nExcludes_global
1045     if (excludesGlobal(i) == unique_id_1) then
1046     skip_it = .true.
1047     return
1048     end if
1049     end do
1050     return
1051     end if
1052    
1053     #ifdef IS_MPI
1054     unique_id_2 = tagColumn(atom2)
1055     #else
1056     unique_id_2 = atom2
1057     #endif
1058 chuckv 441
1059 mmeineke 377 #ifdef IS_MPI
1060     !! this situation should only arise in MPI simulations
1061     if (unique_id_1 == unique_id_2) then
1062     skip_it = .true.
1063     return
1064     end if
1065    
1066     !! this prevents us from doing the pair on multiple processors
1067     if (unique_id_1 < unique_id_2) then
1068 chuckv 441 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1069     skip_it = .true.
1070     return
1071     endif
1072 mmeineke 377 else
1073 chuckv 441 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1074     skip_it = .true.
1075     return
1076     endif
1077 mmeineke 377 endif
1078     #endif
1079 chuckv 441
1080 mmeineke 377 !! the rest of these situations can happen in all simulations:
1081     do i = 1, nExcludes_global
1082     if ((excludesGlobal(i) == unique_id_1) .or. &
1083     (excludesGlobal(i) == unique_id_2)) then
1084     skip_it = .true.
1085     return
1086     endif
1087     enddo
1088 chuckv 441
1089 mmeineke 377 do i = 1, nExcludes_local
1090     if (excludesLocal(1,i) == unique_id_1) then
1091     if (excludesLocal(2,i) == unique_id_2) then
1092     skip_it = .true.
1093     return
1094     endif
1095     else
1096     if (excludesLocal(1,i) == unique_id_2) then
1097     if (excludesLocal(2,i) == unique_id_1) then
1098     skip_it = .true.
1099     return
1100     endif
1101     endif
1102     endif
1103     end do
1104    
1105     return
1106     end function skipThisPair
1107    
1108     function FF_UsesDirectionalAtoms() result(doesit)
1109     logical :: doesit
1110     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1111     FF_uses_GB .or. FF_uses_RF
1112     end function FF_UsesDirectionalAtoms
1113    
1114     function FF_RequiresPrepairCalc() result(doesit)
1115     logical :: doesit
1116     doesit = FF_uses_EAM
1117     end function FF_RequiresPrepairCalc
1118    
1119     function FF_RequiresPostpairCalc() result(doesit)
1120     logical :: doesit
1121     doesit = FF_uses_RF
1122     end function FF_RequiresPostpairCalc
1123    
1124 chuckv 673 !! This cleans componets of force arrays belonging only to fortran
1125    
1126 mmeineke 377 end module do_Forces