ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 845
Committed: Thu Oct 30 18:59:20 2003 UTC (20 years, 8 months ago) by gezelter
File size: 29970 byte(s)
Log Message:
bug fixes for rList problems

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 gezelter 845 !! @version $Id: do_Forces.F90,v 1.36 2003-10-30 18:59:20 gezelter Exp $, $Date: 2003-10-30 18:59:20 $, $Name: not supported by cvs2svn $, $Revision: 1.36 $
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    
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 mmeineke 377 end subroutine do_force_loop
734    
735 chuckv 441 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
736 mmeineke 377
737     real( kind = dp ) :: pot
738 chuckv 460 real( kind = dp ), dimension(3,getNlocal()) :: u_l
739     real (kind=dp), dimension(9,getNlocal()) :: A
740     real (kind=dp), dimension(3,getNlocal()) :: f
741     real (kind=dp), dimension(3,getNlocal()) :: t
742 mmeineke 377
743     logical, intent(inout) :: do_pot, do_stress
744     integer, intent(in) :: i, j
745     real ( kind = dp ), intent(inout) :: rijsq
746     real ( kind = dp ) :: r
747     real ( kind = dp ), intent(inout) :: d(3)
748     logical :: is_LJ_i, is_LJ_j
749     logical :: is_DP_i, is_DP_j
750     logical :: is_GB_i, is_GB_j
751 chuckv 648 logical :: is_EAM_i,is_EAM_j
752 mmeineke 377 logical :: is_Sticky_i, is_Sticky_j
753     integer :: me_i, me_j
754    
755     r = sqrt(rijsq)
756    
757     #ifdef IS_MPI
758 gezelter 490 if (tagRow(i) .eq. tagColumn(j)) then
759     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
760     endif
761 mmeineke 377
762     me_i = atid_row(i)
763     me_j = atid_col(j)
764    
765     #else
766    
767     me_i = atid(i)
768     me_j = atid(j)
769    
770     #endif
771    
772     if (FF_uses_LJ .and. SimUsesLJ()) then
773     call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
774     call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
775    
776     if ( is_LJ_i .and. is_LJ_j ) &
777     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
778     endif
779    
780     if (FF_uses_dipoles .and. SimUsesDipoles()) then
781     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
782     call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
783    
784     if ( is_DP_i .and. is_DP_j ) then
785 gezelter 462 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
786 mmeineke 377 do_pot, do_stress)
787     if (FF_uses_RF .and. SimUsesRF()) then
788     call accumulate_rf(i, j, r, u_l)
789     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
790     endif
791    
792     endif
793     endif
794    
795     if (FF_uses_Sticky .and. SimUsesSticky()) then
796    
797     call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
798     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
799 chuckv 388
800 mmeineke 377 if ( is_Sticky_i .and. is_Sticky_j ) then
801     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
802     do_pot, do_stress)
803     endif
804     endif
805    
806    
807     if (FF_uses_GB .and. SimUsesGB()) then
808    
809 tim 726
810 mmeineke 377 call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
811     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
812    
813     if ( is_GB_i .and. is_GB_j ) then
814     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
815     do_pot, do_stress)
816     endif
817     endif
818    
819 mmeineke 597
820 chuckv 648
821     if (FF_uses_EAM .and. SimUsesEAM()) then
822     call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
823     call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
824    
825     if ( is_EAM_i .and. is_EAM_j ) &
826     call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
827     endif
828    
829 mmeineke 597
830 chuckv 648
831    
832 mmeineke 377 end subroutine do_pair
833    
834    
835 chuckv 631
836 chuckv 648 subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
837 chuckv 631 real( kind = dp ) :: pot
838     real( kind = dp ), dimension(3,getNlocal()) :: u_l
839     real (kind=dp), dimension(9,getNlocal()) :: A
840     real (kind=dp), dimension(3,getNlocal()) :: f
841     real (kind=dp), dimension(3,getNlocal()) :: t
842    
843     logical, intent(inout) :: do_pot, do_stress
844     integer, intent(in) :: i, j
845     real ( kind = dp ), intent(inout) :: rijsq
846     real ( kind = dp ) :: r
847     real ( kind = dp ), intent(inout) :: d(3)
848    
849     logical :: is_EAM_i, is_EAM_j
850    
851     integer :: me_i, me_j
852    
853     r = sqrt(rijsq)
854    
855 chuckv 669
856 chuckv 631 #ifdef IS_MPI
857     if (tagRow(i) .eq. tagColumn(j)) then
858     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
859     endif
860    
861     me_i = atid_row(i)
862     me_j = atid_col(j)
863    
864     #else
865    
866     me_i = atid(i)
867     me_j = atid(j)
868    
869     #endif
870 chuckv 673
871 chuckv 631 if (FF_uses_EAM .and. SimUsesEAM()) then
872     call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
873     call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
874    
875 chuckv 648 if ( is_EAM_i .and. is_EAM_j ) &
876     call calc_EAM_prepair_rho(i, j, d, r, rijsq )
877 chuckv 631 endif
878 chuckv 648
879 chuckv 673 end subroutine do_prepair
880 chuckv 648
881    
882    
883 chuckv 673
884 chuckv 648 subroutine do_preforce(nlocal,pot)
885     integer :: nlocal
886     real( kind = dp ) :: pot
887    
888 chuckv 669 if (FF_uses_EAM .and. SimUsesEAM()) then
889     call calc_EAM_preforce_Frho(nlocal,pot)
890     endif
891 chuckv 648
892    
893 chuckv 631 end subroutine do_preforce
894    
895    
896 mmeineke 377 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
897    
898     real (kind = dp), dimension(3) :: q_i
899     real (kind = dp), dimension(3) :: q_j
900     real ( kind = dp ), intent(out) :: r_sq
901 gezelter 571 real( kind = dp ) :: d(3), scaled(3)
902     integer i
903    
904 chuckv 482 d(1:3) = q_j(1:3) - q_i(1:3)
905 gezelter 571
906 mmeineke 377 ! Wrap back into periodic box if necessary
907     if ( SimUsesPBC() ) then
908 mmeineke 393
909 gezelter 571 if( .not.boxIsOrthorhombic ) then
910     ! calc the scaled coordinates.
911    
912 mmeineke 572 scaled = matmul(HmatInv, d)
913 gezelter 571
914     ! wrap the scaled coordinates
915    
916 mmeineke 572 scaled = scaled - anint(scaled)
917    
918 gezelter 571
919     ! calc the wrapped real coordinates from the wrapped scaled
920     ! coordinates
921    
922 mmeineke 572 d = matmul(Hmat,scaled)
923 gezelter 571
924     else
925     ! calc the scaled coordinates.
926    
927     do i = 1, 3
928     scaled(i) = d(i) * HmatInv(i,i)
929    
930     ! wrap the scaled coordinates
931    
932     scaled(i) = scaled(i) - anint(scaled(i))
933    
934     ! calc the wrapped real coordinates from the wrapped scaled
935     ! coordinates
936    
937     d(i) = scaled(i)*Hmat(i,i)
938     enddo
939     endif
940 mmeineke 393
941 mmeineke 377 endif
942 gezelter 571
943 mmeineke 377 r_sq = dot_product(d,d)
944 gezelter 571
945 mmeineke 377 end subroutine get_interatomic_vector
946 gezelter 571
947 mmeineke 377 subroutine check_initialization(error)
948     integer, intent(out) :: error
949    
950     error = 0
951     ! Make sure we are properly initialized.
952     if (.not. do_forces_initialized) then
953 chuckv 657 write(*,*) "Forces not initialized"
954 mmeineke 377 error = -1
955     return
956     endif
957    
958     #ifdef IS_MPI
959     if (.not. isMPISimSet()) then
960     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
961     error = -1
962     return
963     endif
964     #endif
965    
966     return
967     end subroutine check_initialization
968    
969    
970     subroutine zero_work_arrays()
971    
972     #ifdef IS_MPI
973    
974     q_Row = 0.0_dp
975     q_Col = 0.0_dp
976    
977     u_l_Row = 0.0_dp
978     u_l_Col = 0.0_dp
979    
980     A_Row = 0.0_dp
981     A_Col = 0.0_dp
982    
983     f_Row = 0.0_dp
984     f_Col = 0.0_dp
985     f_Temp = 0.0_dp
986    
987     t_Row = 0.0_dp
988     t_Col = 0.0_dp
989     t_Temp = 0.0_dp
990    
991     pot_Row = 0.0_dp
992     pot_Col = 0.0_dp
993     pot_Temp = 0.0_dp
994    
995     rf_Row = 0.0_dp
996     rf_Col = 0.0_dp
997     rf_Temp = 0.0_dp
998    
999     #endif
1000    
1001 chuckv 673
1002     if (FF_uses_EAM .and. SimUsesEAM()) then
1003     call clean_EAM()
1004     endif
1005    
1006    
1007    
1008    
1009    
1010 mmeineke 377 rf = 0.0_dp
1011     tau_Temp = 0.0_dp
1012     virial_Temp = 0.0_dp
1013     end subroutine zero_work_arrays
1014    
1015     function skipThisPair(atom1, atom2) result(skip_it)
1016     integer, intent(in) :: atom1
1017     integer, intent(in), optional :: atom2
1018     logical :: skip_it
1019     integer :: unique_id_1, unique_id_2
1020 chuckv 388 integer :: me_i,me_j
1021 mmeineke 377 integer :: i
1022    
1023     skip_it = .false.
1024    
1025     !! there are a number of reasons to skip a pair or a particle
1026     !! mostly we do this to exclude atoms who are involved in short
1027     !! range interactions (bonds, bends, torsions), but we also need
1028     !! to exclude some overcounted interactions that result from
1029     !! the parallel decomposition
1030    
1031     #ifdef IS_MPI
1032     !! in MPI, we have to look up the unique IDs for each atom
1033     unique_id_1 = tagRow(atom1)
1034     #else
1035     !! in the normal loop, the atom numbers are unique
1036     unique_id_1 = atom1
1037     #endif
1038 chuckv 388
1039 mmeineke 377 !! We were called with only one atom, so just check the global exclude
1040     !! list for this atom
1041     if (.not. present(atom2)) then
1042     do i = 1, nExcludes_global
1043     if (excludesGlobal(i) == unique_id_1) then
1044     skip_it = .true.
1045     return
1046     end if
1047     end do
1048     return
1049     end if
1050    
1051     #ifdef IS_MPI
1052     unique_id_2 = tagColumn(atom2)
1053     #else
1054     unique_id_2 = atom2
1055     #endif
1056 chuckv 441
1057 mmeineke 377 #ifdef IS_MPI
1058     !! this situation should only arise in MPI simulations
1059     if (unique_id_1 == unique_id_2) then
1060     skip_it = .true.
1061     return
1062     end if
1063    
1064     !! this prevents us from doing the pair on multiple processors
1065     if (unique_id_1 < unique_id_2) then
1066 chuckv 441 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1067     skip_it = .true.
1068     return
1069     endif
1070 mmeineke 377 else
1071 chuckv 441 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1072     skip_it = .true.
1073     return
1074     endif
1075 mmeineke 377 endif
1076     #endif
1077 chuckv 441
1078 mmeineke 377 !! the rest of these situations can happen in all simulations:
1079     do i = 1, nExcludes_global
1080     if ((excludesGlobal(i) == unique_id_1) .or. &
1081     (excludesGlobal(i) == unique_id_2)) then
1082     skip_it = .true.
1083     return
1084     endif
1085     enddo
1086 chuckv 441
1087 mmeineke 377 do i = 1, nExcludes_local
1088     if (excludesLocal(1,i) == unique_id_1) then
1089     if (excludesLocal(2,i) == unique_id_2) then
1090     skip_it = .true.
1091     return
1092     endif
1093     else
1094     if (excludesLocal(1,i) == unique_id_2) then
1095     if (excludesLocal(2,i) == unique_id_1) then
1096     skip_it = .true.
1097     return
1098     endif
1099     endif
1100     endif
1101     end do
1102    
1103     return
1104     end function skipThisPair
1105    
1106     function FF_UsesDirectionalAtoms() result(doesit)
1107     logical :: doesit
1108     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1109     FF_uses_GB .or. FF_uses_RF
1110     end function FF_UsesDirectionalAtoms
1111    
1112     function FF_RequiresPrepairCalc() result(doesit)
1113     logical :: doesit
1114     doesit = FF_uses_EAM
1115     end function FF_RequiresPrepairCalc
1116    
1117     function FF_RequiresPostpairCalc() result(doesit)
1118     logical :: doesit
1119     doesit = FF_uses_RF
1120     end function FF_RequiresPostpairCalc
1121    
1122 chuckv 673 !! This cleans componets of force arrays belonging only to fortran
1123    
1124 mmeineke 377 end module do_Forces