ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 470
Committed: Mon Apr 7 20:50:46 2003 UTC (21 years, 3 months ago) by chuckv
File size: 20960 byte(s)
Log Message:
Fixed transpose bug in mpi reduce for tau and virial.

File Contents

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