ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 368
Committed: Thu Mar 20 00:02:39 2003 UTC (21 years, 5 months ago) by chuckv
File size: 20500 byte(s)
Log Message:
Fixed bugs. Single version now runs w/o segfault. Still a conservation of energy bug.
do_Forces.F90: Fixed pot not being passed to do_pair.
neighborLists.F90: Fixed bugs in checkNeighborLists
atype_module.F90: Fixed bug with allocating atypes on each new_atype call.Now checks to see if atypes is null, then calles initialize(16).
vector_class.F90: Fixed some bugs with how MatchList was being allocated.

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