ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 374
Committed: Thu Mar 20 19:50:42 2003 UTC (22 years, 5 months ago) by chuckv
File size: 20549 byte(s)
Log Message:
*** empty log message ***

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 374 !! @version $Id: do_Forces.F90,v 1.29 2003-03-20 19:50:42 chuckv Exp $, $Date: 2003-03-20 19:50:42 $, $Name: not supported by cvs2svn $, $Revision: 1.29 $
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     neighborListSize = getNeighborListSize()
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     endif
274    
275     list(nlist) = j
276 gezelter 317
277 gezelter 297 if (rijsq < rcutsq) then
278 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
279 chuckv 368 u_l, A, f, t,pot)
280 gezelter 297 endif
281     endif
282     enddo inner
283     enddo
284 chuckv 292
285 gezelter 297 point(nrow + 1) = nlist + 1
286    
287 gezelter 329 else !! (of update_check)
288 chuckv 292
289 gezelter 297 ! use the list to find the neighbors
290     do i = 1, nrow
291     JBEG = POINT(i)
292     JEND = POINT(i+1) - 1
293     ! check thiat molecule i has neighbors
294     if (jbeg .le. jend) then
295 gezelter 317
296 gezelter 297 do jnab = jbeg, jend
297     j = list(jnab)
298 gezelter 317
299     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
300 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
301 chuckv 368 u_l, A, f, t,pot)
302 gezelter 317
303 gezelter 297 enddo
304     endif
305     enddo
306     endif
307 gezelter 329
308 gezelter 297 #else
309    
310     if (update_nlist) then
311 chuckv 292
312 gezelter 297 ! save current configuration, contruct neighbor list,
313     ! and calculate forces
314 gezelter 334 call saveNeighborList(q)
315 gezelter 297
316     neighborListSize = getNeighborListSize()
317     nlist = 0
318 gezelter 329
319 gezelter 297 do i = 1, natoms-1
320     point(i) = nlist + 1
321 gezelter 329
322 gezelter 297 inner: do j = i+1, natoms
323 gezelter 329
324 gezelter 334 if (skipThisPair(i,j)) cycle inner
325 gezelter 329
326 gezelter 317 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
327 chuckv 295
328 gezelter 297 if (rijsq < rlistsq) then
329 gezelter 317
330 gezelter 297 nlist = nlist + 1
331    
332     if (nlist > neighborListSize) then
333 chuckv 374 write(*,*) "do_Forces: Calling expand nlist w/ ",nlist , size(list)
334 mmeineke 367 call expandNeighborList(natoms, listerror)
335 gezelter 297 if (listerror /= 0) then
336 gezelter 329 error = -1
337 gezelter 297 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
338     return
339     end if
340     endif
341    
342     list(nlist) = j
343 gezelter 329
344 gezelter 297 if (rijsq < rcutsq) then
345 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
346 chuckv 368 u_l, A, f, t,pot)
347 gezelter 297 endif
348     endif
349 chuckv 292 enddo inner
350 gezelter 297 enddo
351    
352     point(natoms) = nlist + 1
353    
354     else !! (update)
355    
356     ! use the list to find the neighbors
357 gezelter 330 do i = 1, natoms-1
358 gezelter 297 JBEG = POINT(i)
359     JEND = POINT(i+1) - 1
360     ! check thiat molecule i has neighbors
361     if (jbeg .le. jend) then
362    
363     do jnab = jbeg, jend
364     j = list(jnab)
365 gezelter 317
366     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
367 gezelter 356 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
368 chuckv 368 u_l, A, f, t,pot)
369 gezelter 317
370 gezelter 297 enddo
371     endif
372     enddo
373     endif
374 gezelter 317
375 gezelter 297 #endif
376 gezelter 317
377 gezelter 329 ! phew, done with main loop.
378 gezelter 317
379 chuckv 292 #ifdef IS_MPI
380 gezelter 329 !!distribute forces
381 gezelter 317
382 gezelter 309 call scatter(f_Row,f,plan_row3d)
383     call scatter(f_Col,f_temp,plan_col3d)
384 chuckv 295 do i = 1,nlocal
385 gezelter 309 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
386 chuckv 295 end do
387 gezelter 329
388     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
389 gezelter 309 call scatter(t_Row,t,plan_row3d)
390     call scatter(t_Col,t_temp,plan_col3d)
391 gezelter 329
392 chuckv 295 do i = 1,nlocal
393 gezelter 309 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
394 chuckv 295 end do
395     endif
396 gezelter 309
397 chuckv 295 if (do_pot) then
398     ! scatter/gather pot_row into the members of my column
399 gezelter 309 call scatter(pot_Row, pot_Temp, plan_row)
400 chuckv 295
401     ! scatter/gather pot_local into all other procs
402     ! add resultant to get total pot
403     do i = 1, nlocal
404 gezelter 309 pot_local = pot_local + pot_Temp(i)
405 chuckv 295 enddo
406    
407 gezelter 309 pot_Temp = 0.0_DP
408    
409     call scatter(pot_Col, pot_Temp, plan_col)
410 chuckv 295 do i = 1, nlocal
411 gezelter 309 pot_local = pot_local + pot_Temp(i)
412 chuckv 295 enddo
413    
414 gezelter 330 endif
415     #endif
416    
417 gezelter 329 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
418    
419     if (FF_uses_RF .and. SimUsesRF()) then
420    
421     #ifdef IS_MPI
422     call scatter(rf_Row,rf,plan_row3d)
423     call scatter(rf_Col,rf_Temp,plan_col3d)
424     do i = 1,nlocal
425     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
426     end do
427     #endif
428    
429     do i = 1, getNlocal()
430    
431 gezelter 330 rfpot = 0.0_DP
432 gezelter 329 #ifdef IS_MPI
433     me_i = atid_row(i)
434     #else
435     me_i = atid(i)
436     #endif
437     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
438     if ( is_DP_i ) then
439     call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
440     !! The reaction field needs to include a self contribution
441     !! to the field:
442     call accumulate_self_rf(i, mu_i, u_l)
443     !! Get the reaction field contribution to the
444     !! potential and torques:
445 gezelter 330 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
446 gezelter 329 #ifdef IS_MPI
447     pot_local = pot_local + rfpot
448     #else
449     pot = pot + rfpot
450     #endif
451     endif
452     enddo
453     endif
454     endif
455    
456    
457     #ifdef IS_MPI
458    
459     if (do_pot) then
460 gezelter 309 pot = pot_local
461 gezelter 329 !! we assume the c code will do the allreduce to get the total potential
462     !! we could do it right here if we needed to...
463 chuckv 295 endif
464    
465 gezelter 329 if (do_stress) then
466 gezelter 330 call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
467 gezelter 309 mpi_comm_world,mpi_err)
468 gezelter 330 call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
469 gezelter 309 mpi_comm_world,mpi_err)
470     endif
471 chuckv 295
472 gezelter 329 #else
473 chuckv 295
474 gezelter 329 if (do_stress) then
475 gezelter 309 tau = tau_Temp
476     virial = virial_Temp
477 chuckv 295 endif
478    
479 gezelter 329 #endif
480    
481 chuckv 295 end subroutine do_force_loop
482    
483 chuckv 368 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
484 chuckv 295
485 gezelter 330 real( kind = dp ) :: pot
486 gezelter 364 real( kind = dp ), dimension(:,:) :: u_l
487     real (kind=dp), dimension(:,:) :: A
488     real (kind=dp), dimension(:,:) :: f
489     real (kind=dp), dimension(:,:) :: t
490 gezelter 330
491 gezelter 329 logical, intent(inout) :: do_pot, do_stress
492 gezelter 317 integer, intent(in) :: i, j
493 chuckv 331 real ( kind = dp ), intent(inout) :: rijsq
494 gezelter 317 real ( kind = dp ) :: r
495     real ( kind = dp ), intent(inout) :: d(3)
496     logical :: is_LJ_i, is_LJ_j
497     logical :: is_DP_i, is_DP_j
498 gezelter 364 logical :: is_GB_i, is_GB_j
499 gezelter 317 logical :: is_Sticky_i, is_Sticky_j
500     integer :: me_i, me_j
501 chuckv 295
502 gezelter 330 r = sqrt(rijsq)
503 chuckv 368
504 gezelter 302 #ifdef IS_MPI
505    
506 gezelter 317 me_i = atid_row(i)
507     me_j = atid_col(j)
508 chuckv 295
509 gezelter 317 #else
510 chuckv 295
511 gezelter 317 me_i = atid(i)
512     me_j = atid(j)
513 gezelter 302
514 gezelter 317 #endif
515 gezelter 302
516 gezelter 329 if (FF_uses_LJ .and. SimUsesLJ()) then
517     call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
518     call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
519 chuckv 368 if ( is_LJ_i .and. is_LJ_j ) &
520     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
521 gezelter 329 endif
522    
523 gezelter 330 if (FF_uses_dipoles .and. SimUsesDipoles()) then
524 gezelter 329 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
525     call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
526    
527     if ( is_DP_i .and. is_DP_j ) then
528    
529 gezelter 360 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
530     do_pot, do_stress)
531 gezelter 329
532     if (FF_uses_RF .and. SimUsesRF()) then
533    
534     call accumulate_rf(i, j, r, u_l)
535     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
536    
537     endif
538    
539 gezelter 297 endif
540     endif
541    
542 gezelter 329 if (FF_uses_Sticky .and. SimUsesSticky()) then
543 gezelter 317
544 gezelter 329 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
545     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
546    
547     if ( is_Sticky_i .and. is_Sticky_j ) then
548     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
549     do_pot, do_stress)
550     endif
551 gezelter 297 endif
552 gezelter 364
553    
554     if (FF_uses_GB .and. SimUsesGB()) then
555    
556     call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
557     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
558    
559     if ( is_GB_i .and. is_GB_j ) then
560     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
561     do_pot, do_stress)
562     endif
563     endif
564 gezelter 309
565 chuckv 295 end subroutine do_pair
566 chuckv 292
567    
568 gezelter 317 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
569    
570 chuckv 295 real (kind = dp), dimension(3) :: q_i
571     real (kind = dp), dimension(3) :: q_j
572     real ( kind = dp ), intent(out) :: r_sq
573     real( kind = dp ) :: d(3)
574 mmeineke 370 real( kind = dp ) :: d_old(3)
575 chuckv 295 d(1:3) = q_i(1:3) - q_j(1:3)
576 mmeineke 370 d_old = d
577 gezelter 317 ! Wrap back into periodic box if necessary
578 gezelter 330 if ( SimUsesPBC() ) then
579 mmeineke 370
580     d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
581     int(abs(d(1:3)/box(1:3)) + 0.5_dp)
582    
583 gezelter 317 endif
584 chuckv 295 r_sq = dot_product(d,d)
585 gezelter 317
586 chuckv 295 end subroutine get_interatomic_vector
587 gezelter 329
588     subroutine check_initialization(error)
589     integer, intent(out) :: error
590 gezelter 332
591 gezelter 329 error = 0
592     ! Make sure we are properly initialized.
593 gezelter 332 if (.not. do_forces_initialized) then
594 gezelter 329 error = -1
595     return
596     endif
597 gezelter 332
598 gezelter 329 #ifdef IS_MPI
599     if (.not. isMPISimSet()) then
600     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
601     error = -1
602     return
603     endif
604     #endif
605 gezelter 332
606 gezelter 329 return
607     end subroutine check_initialization
608    
609 gezelter 317
610 gezelter 325 subroutine zero_work_arrays()
611    
612     #ifdef IS_MPI
613 chuckv 292
614 gezelter 325 q_Row = 0.0_dp
615     q_Col = 0.0_dp
616 chuckv 295
617 gezelter 325 u_l_Row = 0.0_dp
618     u_l_Col = 0.0_dp
619 chuckv 295
620 gezelter 325 A_Row = 0.0_dp
621     A_Col = 0.0_dp
622 chuckv 295
623 gezelter 325 f_Row = 0.0_dp
624     f_Col = 0.0_dp
625     f_Temp = 0.0_dp
626    
627     t_Row = 0.0_dp
628     t_Col = 0.0_dp
629     t_Temp = 0.0_dp
630 chuckv 295
631 gezelter 325 pot_Row = 0.0_dp
632     pot_Col = 0.0_dp
633     pot_Temp = 0.0_dp
634 chuckv 292
635 gezelter 332 rf_Row = 0.0_dp
636     rf_Col = 0.0_dp
637     rf_Temp = 0.0_dp
638    
639 chuckv 295 #endif
640 chuckv 292
641 gezelter 332 rf = 0.0_dp
642 gezelter 325 tau_Temp = 0.0_dp
643     virial_Temp = 0.0_dp
644    
645     end subroutine zero_work_arrays
646    
647 gezelter 334 function skipThisPair(atom1, atom2) result(skip_it)
648 gezelter 332
649 gezelter 334 integer, intent(in) :: atom1
650     integer, intent(in), optional :: atom2
651     logical :: skip_it
652 gezelter 332 integer :: unique_id_1, unique_id_2
653 gezelter 334 integer :: i
654 gezelter 332
655 gezelter 334 skip_it = .false.
656    
657     !! there are a number of reasons to skip a pair or a particle
658     !! mostly we do this to exclude atoms who are involved in short
659     !! range interactions (bonds, bends, torsions), but we also need
660     !! to exclude some overcounted interactions that result from
661     !! the parallel decomposition
662 gezelter 330
663 chuckv 306 #ifdef IS_MPI
664 gezelter 334 !! in MPI, we have to look up the unique IDs for each atom
665 gezelter 332 unique_id_1 = tagRow(atom1)
666 chuckv 306 #else
667 gezelter 334 !! in the normal loop, the atom numbers are unique
668     unique_id_1 = atom1
669 chuckv 306 #endif
670 gezelter 330
671 gezelter 334 !! We were called with only one atom, so just check the global exclude
672     !! list for this atom
673 chuckv 306 if (.not. present(atom2)) then
674 gezelter 330 do i = 1, nExcludes_global
675 gezelter 332 if (excludesGlobal(i) == unique_id_1) then
676 gezelter 334 skip_it = .true.
677 chuckv 306 return
678     end if
679     end do
680 gezelter 334 return
681 chuckv 306 end if
682 gezelter 330
683 gezelter 332 #ifdef IS_MPI
684     unique_id_2 = tagColumn(atom2)
685     #else
686 gezelter 334 unique_id_2 = atom2
687 gezelter 332 #endif
688    
689 gezelter 334 #ifdef IS_MPI
690     !! this situation should only arise in MPI simulations
691 gezelter 332 if (unique_id_1 == unique_id_2) then
692 gezelter 334 skip_it = .true.
693 chuckv 295 return
694     end if
695 gezelter 330
696 gezelter 334 !! this prevents us from doing the pair on multiple processors
697 gezelter 332 if (unique_id_1 < unique_id_2) then
698 gezelter 334 if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
699 chuckv 295 return
700     else
701 gezelter 334 if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
702 chuckv 295 endif
703 gezelter 334 #endif
704    
705     !! the rest of these situations can happen in all simulations:
706     do i = 1, nExcludes_global
707     if ((excludesGlobal(i) == unique_id_1) .or. &
708     (excludesGlobal(i) == unique_id_2)) then
709     skip_it = .true.
710     return
711     endif
712     enddo
713 gezelter 332
714 gezelter 330 do i = 1, nExcludes_local
715 gezelter 334 if (excludesLocal(1,i) == unique_id_1) then
716     if (excludesLocal(2,i) == unique_id_2) then
717     skip_it = .true.
718     return
719     endif
720     else
721     if (excludesLocal(1,i) == unique_id_2) then
722     if (excludesLocal(2,i) == unique_id_1) then
723     skip_it = .true.
724     return
725     endif
726     endif
727     endif
728 chuckv 306 end do
729 gezelter 330
730 gezelter 334 return
731     end function skipThisPair
732 chuckv 306
733 gezelter 329 function FF_UsesDirectionalAtoms() result(doesit)
734     logical :: doesit
735     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
736     FF_uses_GB .or. FF_uses_RF
737     end function FF_UsesDirectionalAtoms
738    
739     function FF_RequiresPrepairCalc() result(doesit)
740     logical :: doesit
741     doesit = FF_uses_EAM
742     end function FF_RequiresPrepairCalc
743    
744     function FF_RequiresPostpairCalc() result(doesit)
745     logical :: doesit
746     doesit = FF_uses_RF
747     end function FF_RequiresPostpairCalc
748    
749 chuckv 292 end module do_Forces