ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 375
Committed: Fri Mar 21 15:07:14 2003 UTC (22 years, 5 months ago) by chuckv
File size: 20550 byte(s)
Log Message:
more bug fixes...

File Contents

# User Rev Content
1 chuckv 306 !! do_Forces.F90
2     !! module do_Forces
3 chuckv 292 !! Calculates Long Range forces.
4 chuckv 306
5 chuckv 292 !! @author Charles F. Vardeman II
6     !! @author Matthew Meineke
7 chuckv 375 !! @version $Id: do_Forces.F90,v 1.30 2003-03-21 15:07:14 chuckv Exp $, $Date: 2003-03-21 15:07:14 $, $Name: not supported by cvs2svn $, $Revision: 1.30 $
8 chuckv 292
9     module do_Forces
10 gezelter 360 use force_globals
11 chuckv 292 use simulation
12     use definitions
13 gezelter 328 use atype_module
14 gezelter 325 use neighborLists
15 gezelter 330 use lj
16     use sticky_pair
17 gezelter 309 use dipole_dipole
18 gezelter 332 use reaction_field
19 gezelter 364 use gb_pair
20 chuckv 292 #ifdef IS_MPI
21     use mpiSimulation
22     #endif
23 gezelter 356
24 chuckv 292 implicit none
25     PRIVATE
26    
27 gezelter 356 #define __FORTRAN90
28     #include "fForceField.h"
29    
30 gezelter 329 logical, save :: do_forces_initialized = .false.
31     logical, save :: FF_uses_LJ
32     logical, save :: FF_uses_sticky
33     logical, save :: FF_uses_dipoles
34     logical, save :: FF_uses_RF
35     logical, save :: FF_uses_GB
36     logical, save :: FF_uses_EAM
37 chuckv 292
38 gezelter 329 public :: init_FF
39 gezelter 330 public :: do_force_loop
40 gezelter 329
41 chuckv 292 contains
42    
43 gezelter 358 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
44 gezelter 356
45 gezelter 358 integer, intent(in) :: LJMIXPOLICY
46 mmeineke 359 logical, intent(in) :: use_RF_c
47 gezelter 356
48 chuckv 331 integer, intent(out) :: thisStat
49 gezelter 332 integer :: my_status, nMatches
50 chuckv 368 integer, pointer :: MatchList(:) => null()
51 gezelter 360 real(kind=dp) :: rcut, rrf, rt, dielect
52 gezelter 329
53 gezelter 332 !! assume things are copacetic, unless they aren't
54 chuckv 331 thisStat = 0
55 gezelter 356
56     !! Fortran's version of a cast:
57 gezelter 358 FF_uses_RF = use_RF_c
58 gezelter 332
59     !! init_FF is called *after* all of the atom types have been
60     !! defined in atype_module using the new_atype subroutine.
61     !!
62     !! this will scan through the known atypes and figure out what
63     !! interactions are used by the force field.
64 chuckv 368
65 gezelter 332 FF_uses_LJ = .false.
66     FF_uses_sticky = .false.
67     FF_uses_dipoles = .false.
68     FF_uses_GB = .false.
69     FF_uses_EAM = .false.
70    
71     call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
72     if (nMatches .gt. 0) FF_uses_LJ = .true.
73 gezelter 358
74 gezelter 332 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
75     if (nMatches .gt. 0) FF_uses_dipoles = .true.
76 gezelter 358
77 gezelter 332 call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
78     MatchList)
79     if (nMatches .gt. 0) FF_uses_Sticky = .true.
80    
81     call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
82     if (nMatches .gt. 0) FF_uses_GB = .true.
83 gezelter 358
84 gezelter 332 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
85     if (nMatches .gt. 0) FF_uses_EAM = .true.
86 gezelter 358
87 gezelter 356 !! check to make sure the FF_uses_RF setting makes sense
88 gezelter 358
89 gezelter 360 if (FF_uses_dipoles) then
90     rrf = getRrf()
91     rt = getRt()
92     call initialize_dipole(rrf, rt)
93     if (FF_uses_RF) then
94     dielect = getDielect()
95     call initialize_rf(rrf, rt, dielect)
96     endif
97     else
98     if (FF_uses_RF) then
99 gezelter 332 write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
100     thisStat = -1
101     return
102     endif
103     endif
104 gezelter 356
105     if (FF_uses_LJ) then
106 gezelter 358
107 gezelter 360 call getRcut(rcut)
108    
109 gezelter 358 select case (LJMIXPOLICY)
110     case (LB_MIXING_RULE)
111 gezelter 360 call init_lj_FF(LB_MIXING_RULE, rcut, my_status)
112 gezelter 358 case (EXPLICIT_MIXING_RULE)
113 gezelter 360 call init_lj_FF(EXPLICIT_MIXING_RULE, rcut, my_status)
114 gezelter 356 case default
115     write(default_error,*) 'unknown LJ Mixing Policy!'
116     thisStat = -1
117     return
118 gezelter 358 end select
119 gezelter 356 if (my_status /= 0) then
120     thisStat = -1
121     return
122     end if
123     endif
124    
125     if (FF_uses_sticky) then
126     call check_sticky_FF(my_status)
127     if (my_status /= 0) then
128     thisStat = -1
129     return
130     end if
131     endif
132 gezelter 332
133 gezelter 364 if (FF_uses_GB) then
134     call check_gb_pair_FF(my_status)
135     if (my_status .ne. 0) then
136     thisStat = -1
137     return
138     endif
139     endif
140    
141     if (FF_uses_GB .and. FF_uses_LJ) then
142     endif
143    
144    
145 gezelter 329 do_forces_initialized = .true.
146 chuckv 368
147 gezelter 329 end subroutine init_FF
148 gezelter 358
149 gezelter 329
150 gezelter 317 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
151 gezelter 358 !------------------------------------------------------------->
152 gezelter 325 subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, &
153 gezelter 329 error)
154 gezelter 317 !! Position array provided by C, dimensioned by getNlocal
155 chuckv 292 real ( kind = dp ), dimension(3,getNlocal()) :: q
156 gezelter 317 !! Rotation Matrix for each long range particle in simulation.
157 gezelter 325 real( kind = dp), dimension(9,getNlocal()) :: A
158 gezelter 317 !! Unit vectors for dipoles (lab frame)
159 chuckv 292 real( kind = dp ), dimension(3,getNlocal()) :: u_l
160 gezelter 317 !! Force array provided by C, dimensioned by getNlocal
161 chuckv 292 real ( kind = dp ), dimension(3,getNlocal()) :: f
162 gezelter 317 !! Torsion array provided by C, dimensioned by getNlocal
163 gezelter 325 real( kind = dp ), dimension(3,getNlocal()) :: t
164 gezelter 317 !! Stress Tensor
165 gezelter 325 real( kind = dp), dimension(9) :: tau
166     real ( kind = dp ) :: pot
167     logical ( kind = 2) :: do_pot_c, do_stress_c
168 gezelter 329 logical :: do_pot
169     logical :: do_stress
170 gezelter 317 #ifdef IS_MPI
171     real( kind = DP ) :: pot_local
172 gezelter 329 integer :: nrow
173     integer :: ncol
174 gezelter 317 #endif
175 gezelter 330 integer :: nlocal
176 gezelter 329 integer :: natoms
177 gezelter 325 logical :: update_nlist
178     integer :: i, j, jbeg, jend, jnab
179 gezelter 317 integer :: nlist
180 gezelter 325 real( kind = DP ) :: rijsq, rlistsq, rcutsq, rlist, rcut
181 gezelter 330 real(kind=dp),dimension(3) :: d
182     real(kind=dp) :: rfpot, mu_i, virial
183     integer :: me_i
184     logical :: is_dp_i
185 gezelter 317 integer :: neighborListSize
186 gezelter 330 integer :: listerror, error
187 chuckv 331 integer :: localError
188 chuckv 292
189 gezelter 329 !! initialize local variables
190 gezelter 325
191 chuckv 292 #ifdef IS_MPI
192 gezelter 329 nlocal = getNlocal()
193     nrow = getNrow(plan_row)
194     ncol = getNcol(plan_col)
195     #else
196     nlocal = getNlocal()
197     natoms = nlocal
198 chuckv 292 #endif
199 gezelter 329
200 chuckv 331 call getRcut(rcut,rc2=rcutsq)
201 gezelter 317 call getRlist(rlist,rlistsq)
202    
203 chuckv 331 call check_initialization(localError)
204     if ( localError .ne. 0 ) then
205     error = -1
206     return
207     end if
208 gezelter 329 call zero_work_arrays()
209    
210     do_pot = do_pot_c
211     do_stress = do_stress_c
212    
213     ! Gather all information needed by all force loops:
214 gezelter 317
215 gezelter 329 #ifdef IS_MPI
216 chuckv 292
217 gezelter 309 call gather(q,q_Row,plan_row3d)
218     call gather(q,q_Col,plan_col3d)
219 gezelter 329
220     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
221     call gather(u_l,u_l_Row,plan_row3d)
222     call gather(u_l,u_l_Col,plan_col3d)
223    
224     call gather(A,A_Row,plan_row_rotation)
225     call gather(A,A_Col,plan_col_rotation)
226     endif
227 gezelter 317
228 gezelter 329 #endif
229 gezelter 317
230 gezelter 329 if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
231     !! See if we need to update neighbor lists
232     call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
233     !! if_mpi_gather_stuff_for_prepair
234     !! do_prepair_loop_if_needed
235     !! if_mpi_scatter_stuff_from_prepair
236     !! if_mpi_gather_stuff_from_prepair_to_main_loop
237     else
238     !! See if we need to update neighbor lists
239     call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
240     endif
241    
242 gezelter 297 #ifdef IS_MPI
243 chuckv 292
244 gezelter 297 if (update_nlist) then
245    
246 gezelter 329 !! save current configuration, construct neighbor list,
247     !! and calculate forces
248 gezelter 334 call saveNeighborList(q)
249 gezelter 297
250 chuckv 375 neighborListSize = size(list)
251 gezelter 329 nlist = 0
252 gezelter 297
253     do i = 1, nrow
254     point(i) = nlist + 1
255    
256     inner: do j = 1, ncol
257    
258 gezelter 334 if (skipThisPair(i,j)) cycle inner
259 gezelter 329
260 gezelter 317 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
261 gezelter 330
262 gezelter 297 if (rijsq < rlistsq) then
263    
264     nlist = nlist + 1
265    
266     if (nlist > neighborListSize) then
267 gezelter 325 call expandNeighborList(nlocal, listerror)
268 gezelter 297 if (listerror /= 0) then
269 gezelter 329 error = -1
270 gezelter 297 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
271     return
272     end if
273 chuckv 375 neighborListSize = size(list)
274 gezelter 297 endif
275    
276     list(nlist) = j
277 gezelter 317
278 gezelter 297 if (rijsq < rcutsq) then
279 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
280 chuckv 368 u_l, A, f, t,pot)
281 gezelter 297 endif
282     endif
283     enddo inner
284     enddo
285 chuckv 292
286 gezelter 297 point(nrow + 1) = nlist + 1
287    
288 gezelter 329 else !! (of update_check)
289 chuckv 292
290 gezelter 297 ! use the list to find the neighbors
291     do i = 1, nrow
292     JBEG = POINT(i)
293     JEND = POINT(i+1) - 1
294     ! check thiat molecule i has neighbors
295     if (jbeg .le. jend) then
296 gezelter 317
297 gezelter 297 do jnab = jbeg, jend
298     j = list(jnab)
299 gezelter 317
300     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
301 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
302 chuckv 368 u_l, A, f, t,pot)
303 gezelter 317
304 gezelter 297 enddo
305     endif
306     enddo
307     endif
308 gezelter 329
309 gezelter 297 #else
310    
311     if (update_nlist) then
312 chuckv 292
313 gezelter 297 ! save current configuration, contruct neighbor list,
314     ! and calculate forces
315 gezelter 334 call saveNeighborList(q)
316 gezelter 297
317 chuckv 375 neighborListSize = size(list)
318    
319 gezelter 297 nlist = 0
320 gezelter 329
321 gezelter 297 do i = 1, natoms-1
322     point(i) = nlist + 1
323 gezelter 329
324 gezelter 297 inner: do j = i+1, natoms
325 gezelter 329
326 gezelter 334 if (skipThisPair(i,j)) cycle inner
327 gezelter 329
328 gezelter 317 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
329 chuckv 295
330 gezelter 297 if (rijsq < rlistsq) then
331 gezelter 317
332 gezelter 297 nlist = nlist + 1
333 chuckv 375
334 gezelter 297 if (nlist > neighborListSize) then
335 mmeineke 367 call expandNeighborList(natoms, listerror)
336 gezelter 297 if (listerror /= 0) then
337 gezelter 329 error = -1
338 gezelter 297 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
339     return
340     end if
341 chuckv 375 neighborListSize = size(list)
342 gezelter 297 endif
343    
344     list(nlist) = j
345 gezelter 329
346 gezelter 297 if (rijsq < rcutsq) then
347 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
348 chuckv 368 u_l, A, f, t,pot)
349 gezelter 297 endif
350     endif
351 chuckv 292 enddo inner
352 gezelter 297 enddo
353    
354     point(natoms) = nlist + 1
355    
356     else !! (update)
357    
358     ! use the list to find the neighbors
359 gezelter 330 do i = 1, natoms-1
360 gezelter 297 JBEG = POINT(i)
361     JEND = POINT(i+1) - 1
362     ! check thiat molecule i has neighbors
363     if (jbeg .le. jend) then
364    
365     do jnab = jbeg, jend
366     j = list(jnab)
367 gezelter 317
368     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
369 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
370 chuckv 368 u_l, A, f, t,pot)
371 gezelter 317
372 gezelter 297 enddo
373     endif
374     enddo
375     endif
376 gezelter 317
377 gezelter 297 #endif
378 gezelter 317
379 gezelter 329 ! phew, done with main loop.
380 gezelter 317
381 chuckv 292 #ifdef IS_MPI
382 gezelter 329 !!distribute forces
383 gezelter 317
384 gezelter 309 call scatter(f_Row,f,plan_row3d)
385     call scatter(f_Col,f_temp,plan_col3d)
386 chuckv 295 do i = 1,nlocal
387 gezelter 309 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
388 chuckv 295 end do
389 gezelter 329
390     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
391 gezelter 309 call scatter(t_Row,t,plan_row3d)
392     call scatter(t_Col,t_temp,plan_col3d)
393 gezelter 329
394 chuckv 295 do i = 1,nlocal
395 gezelter 309 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
396 chuckv 295 end do
397     endif
398 gezelter 309
399 chuckv 295 if (do_pot) then
400     ! scatter/gather pot_row into the members of my column
401 gezelter 309 call scatter(pot_Row, pot_Temp, plan_row)
402 chuckv 295
403     ! scatter/gather pot_local into all other procs
404     ! add resultant to get total pot
405     do i = 1, nlocal
406 gezelter 309 pot_local = pot_local + pot_Temp(i)
407 chuckv 295 enddo
408    
409 gezelter 309 pot_Temp = 0.0_DP
410    
411     call scatter(pot_Col, pot_Temp, plan_col)
412 chuckv 295 do i = 1, nlocal
413 gezelter 309 pot_local = pot_local + pot_Temp(i)
414 chuckv 295 enddo
415    
416 gezelter 330 endif
417     #endif
418    
419 gezelter 329 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
420    
421     if (FF_uses_RF .and. SimUsesRF()) then
422    
423     #ifdef IS_MPI
424     call scatter(rf_Row,rf,plan_row3d)
425     call scatter(rf_Col,rf_Temp,plan_col3d)
426     do i = 1,nlocal
427     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
428     end do
429     #endif
430    
431     do i = 1, getNlocal()
432    
433 gezelter 330 rfpot = 0.0_DP
434 gezelter 329 #ifdef IS_MPI
435     me_i = atid_row(i)
436     #else
437     me_i = atid(i)
438     #endif
439     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
440     if ( is_DP_i ) then
441     call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
442     !! The reaction field needs to include a self contribution
443     !! to the field:
444     call accumulate_self_rf(i, mu_i, u_l)
445     !! Get the reaction field contribution to the
446     !! potential and torques:
447 gezelter 330 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
448 gezelter 329 #ifdef IS_MPI
449     pot_local = pot_local + rfpot
450     #else
451     pot = pot + rfpot
452 chuckv 375
453 gezelter 329 #endif
454     endif
455     enddo
456     endif
457     endif
458    
459    
460     #ifdef IS_MPI
461    
462     if (do_pot) then
463 gezelter 309 pot = pot_local
464 gezelter 329 !! we assume the c code will do the allreduce to get the total potential
465     !! we could do it right here if we needed to...
466 chuckv 295 endif
467    
468 gezelter 329 if (do_stress) then
469 gezelter 330 call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
470 gezelter 309 mpi_comm_world,mpi_err)
471 gezelter 330 call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
472 gezelter 309 mpi_comm_world,mpi_err)
473     endif
474 chuckv 295
475 gezelter 329 #else
476 chuckv 295
477 gezelter 329 if (do_stress) then
478 gezelter 309 tau = tau_Temp
479     virial = virial_Temp
480 chuckv 295 endif
481    
482 gezelter 329 #endif
483    
484 chuckv 295 end subroutine do_force_loop
485    
486 chuckv 368 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
487 chuckv 295
488 gezelter 330 real( kind = dp ) :: pot
489 gezelter 364 real( kind = dp ), dimension(:,:) :: u_l
490     real (kind=dp), dimension(:,:) :: A
491     real (kind=dp), dimension(:,:) :: f
492     real (kind=dp), dimension(:,:) :: t
493 gezelter 330
494 gezelter 329 logical, intent(inout) :: do_pot, do_stress
495 gezelter 317 integer, intent(in) :: i, j
496 chuckv 331 real ( kind = dp ), intent(inout) :: rijsq
497 gezelter 317 real ( kind = dp ) :: r
498     real ( kind = dp ), intent(inout) :: d(3)
499     logical :: is_LJ_i, is_LJ_j
500     logical :: is_DP_i, is_DP_j
501 gezelter 364 logical :: is_GB_i, is_GB_j
502 gezelter 317 logical :: is_Sticky_i, is_Sticky_j
503     integer :: me_i, me_j
504 chuckv 295
505 gezelter 330 r = sqrt(rijsq)
506 chuckv 368
507 gezelter 302 #ifdef IS_MPI
508    
509 gezelter 317 me_i = atid_row(i)
510     me_j = atid_col(j)
511 chuckv 295
512 gezelter 317 #else
513 chuckv 295
514 gezelter 317 me_i = atid(i)
515     me_j = atid(j)
516 gezelter 302
517 gezelter 317 #endif
518 gezelter 302
519 gezelter 329 if (FF_uses_LJ .and. SimUsesLJ()) then
520     call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
521     call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
522 chuckv 375
523     if ( is_LJ_i .and. is_LJ_j ) &
524     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
525 gezelter 329 endif
526 chuckv 375
527 gezelter 330 if (FF_uses_dipoles .and. SimUsesDipoles()) then
528 gezelter 329 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
529     call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
530    
531     if ( is_DP_i .and. is_DP_j ) then
532    
533 gezelter 360 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
534     do_pot, do_stress)
535 gezelter 329
536     if (FF_uses_RF .and. SimUsesRF()) then
537    
538     call accumulate_rf(i, j, r, u_l)
539     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
540    
541     endif
542    
543 gezelter 297 endif
544     endif
545    
546 gezelter 329 if (FF_uses_Sticky .and. SimUsesSticky()) then
547 gezelter 317
548 gezelter 329 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
549     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
550    
551     if ( is_Sticky_i .and. is_Sticky_j ) then
552     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
553     do_pot, do_stress)
554     endif
555 gezelter 297 endif
556 gezelter 364
557    
558     if (FF_uses_GB .and. SimUsesGB()) then
559    
560     call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
561     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
562    
563     if ( is_GB_i .and. is_GB_j ) then
564     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
565     do_pot, do_stress)
566     endif
567     endif
568 chuckv 375
569 chuckv 295 end subroutine do_pair
570 chuckv 292
571    
572 gezelter 317 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
573    
574 chuckv 295 real (kind = dp), dimension(3) :: q_i
575     real (kind = dp), dimension(3) :: q_j
576     real ( kind = dp ), intent(out) :: r_sq
577     real( kind = dp ) :: d(3)
578 mmeineke 370 real( kind = dp ) :: d_old(3)
579 chuckv 295 d(1:3) = q_i(1:3) - q_j(1:3)
580 mmeineke 370 d_old = d
581 gezelter 317 ! Wrap back into periodic box if necessary
582 gezelter 330 if ( SimUsesPBC() ) then
583 mmeineke 370
584     d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
585     int(abs(d(1:3)/box(1:3)) + 0.5_dp)
586    
587 gezelter 317 endif
588 chuckv 295 r_sq = dot_product(d,d)
589 gezelter 317
590 chuckv 295 end subroutine get_interatomic_vector
591 gezelter 329
592     subroutine check_initialization(error)
593     integer, intent(out) :: error
594 gezelter 332
595 gezelter 329 error = 0
596     ! Make sure we are properly initialized.
597 gezelter 332 if (.not. do_forces_initialized) then
598 gezelter 329 error = -1
599     return
600     endif
601 gezelter 332
602 gezelter 329 #ifdef IS_MPI
603     if (.not. isMPISimSet()) then
604     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
605     error = -1
606     return
607     endif
608     #endif
609 gezelter 332
610 gezelter 329 return
611     end subroutine check_initialization
612    
613 gezelter 317
614 gezelter 325 subroutine zero_work_arrays()
615    
616     #ifdef IS_MPI
617 chuckv 292
618 gezelter 325 q_Row = 0.0_dp
619     q_Col = 0.0_dp
620 chuckv 295
621 gezelter 325 u_l_Row = 0.0_dp
622     u_l_Col = 0.0_dp
623 chuckv 295
624 gezelter 325 A_Row = 0.0_dp
625     A_Col = 0.0_dp
626 chuckv 295
627 gezelter 325 f_Row = 0.0_dp
628     f_Col = 0.0_dp
629     f_Temp = 0.0_dp
630    
631     t_Row = 0.0_dp
632     t_Col = 0.0_dp
633     t_Temp = 0.0_dp
634 chuckv 295
635 gezelter 325 pot_Row = 0.0_dp
636     pot_Col = 0.0_dp
637     pot_Temp = 0.0_dp
638 chuckv 292
639 gezelter 332 rf_Row = 0.0_dp
640     rf_Col = 0.0_dp
641     rf_Temp = 0.0_dp
642    
643 chuckv 295 #endif
644 chuckv 292
645 gezelter 332 rf = 0.0_dp
646 gezelter 325 tau_Temp = 0.0_dp
647     virial_Temp = 0.0_dp
648    
649     end subroutine zero_work_arrays
650    
651 gezelter 334 function skipThisPair(atom1, atom2) result(skip_it)
652 gezelter 332
653 gezelter 334 integer, intent(in) :: atom1
654     integer, intent(in), optional :: atom2
655     logical :: skip_it
656 gezelter 332 integer :: unique_id_1, unique_id_2
657 gezelter 334 integer :: i
658 gezelter 332
659 gezelter 334 skip_it = .false.
660 chuckv 375
661 gezelter 334 !! there are a number of reasons to skip a pair or a particle
662     !! mostly we do this to exclude atoms who are involved in short
663     !! range interactions (bonds, bends, torsions), but we also need
664     !! to exclude some overcounted interactions that result from
665     !! the parallel decomposition
666 gezelter 330
667 chuckv 306 #ifdef IS_MPI
668 gezelter 334 !! in MPI, we have to look up the unique IDs for each atom
669 gezelter 332 unique_id_1 = tagRow(atom1)
670 chuckv 306 #else
671 gezelter 334 !! in the normal loop, the atom numbers are unique
672     unique_id_1 = atom1
673 chuckv 306 #endif
674 gezelter 330
675 gezelter 334 !! We were called with only one atom, so just check the global exclude
676     !! list for this atom
677 chuckv 306 if (.not. present(atom2)) then
678 gezelter 330 do i = 1, nExcludes_global
679 gezelter 332 if (excludesGlobal(i) == unique_id_1) then
680 gezelter 334 skip_it = .true.
681 chuckv 306 return
682     end if
683     end do
684 gezelter 334 return
685 chuckv 306 end if
686 gezelter 330
687 gezelter 332 #ifdef IS_MPI
688     unique_id_2 = tagColumn(atom2)
689     #else
690 gezelter 334 unique_id_2 = atom2
691 gezelter 332 #endif
692    
693 gezelter 334 #ifdef IS_MPI
694     !! this situation should only arise in MPI simulations
695 gezelter 332 if (unique_id_1 == unique_id_2) then
696 gezelter 334 skip_it = .true.
697 chuckv 295 return
698     end if
699 gezelter 330
700 gezelter 334 !! this prevents us from doing the pair on multiple processors
701 gezelter 332 if (unique_id_1 < unique_id_2) then
702 gezelter 334 if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
703 chuckv 295 return
704     else
705 gezelter 334 if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
706 chuckv 375 return
707 chuckv 295 endif
708 gezelter 334 #endif
709    
710     !! the rest of these situations can happen in all simulations:
711     do i = 1, nExcludes_global
712     if ((excludesGlobal(i) == unique_id_1) .or. &
713     (excludesGlobal(i) == unique_id_2)) then
714     skip_it = .true.
715     return
716     endif
717     enddo
718 chuckv 375
719 gezelter 330 do i = 1, nExcludes_local
720 gezelter 334 if (excludesLocal(1,i) == unique_id_1) then
721     if (excludesLocal(2,i) == unique_id_2) then
722     skip_it = .true.
723     return
724     endif
725     else
726     if (excludesLocal(1,i) == unique_id_2) then
727     if (excludesLocal(2,i) == unique_id_1) then
728     skip_it = .true.
729     return
730     endif
731     endif
732     endif
733 chuckv 306 end do
734 gezelter 330
735 gezelter 334 return
736     end function skipThisPair
737 chuckv 306
738 gezelter 329 function FF_UsesDirectionalAtoms() result(doesit)
739     logical :: doesit
740     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
741     FF_uses_GB .or. FF_uses_RF
742     end function FF_UsesDirectionalAtoms
743    
744     function FF_RequiresPrepairCalc() result(doesit)
745     logical :: doesit
746     doesit = FF_uses_EAM
747     end function FF_RequiresPrepairCalc
748    
749     function FF_RequiresPostpairCalc() result(doesit)
750     logical :: doesit
751     doesit = FF_uses_RF
752     end function FF_RequiresPostpairCalc
753    
754 chuckv 292 end module do_Forces