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 (21 years, 5 months ago) by chuckv
File size: 20549 byte(s)
Log Message:
*** empty log message ***

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.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
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 nlocal = getNlocal()
193 nrow = getNrow(plan_row)
194 ncol = getNcol(plan_col)
195 #else
196 nlocal = getNlocal()
197 natoms = nlocal
198 #endif
199
200 call getRcut(rcut,rc2=rcutsq)
201 call getRlist(rlist,rlistsq)
202
203 call check_initialization(localError)
204 if ( localError .ne. 0 ) then
205 error = -1
206 return
207 end if
208 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
215 #ifdef IS_MPI
216
217 call gather(q,q_Row,plan_row3d)
218 call gather(q,q_Col,plan_col3d)
219
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
228 #endif
229
230 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 #ifdef IS_MPI
243
244 if (update_nlist) then
245
246 !! save current configuration, construct neighbor list,
247 !! and calculate forces
248 call saveNeighborList(q)
249
250 neighborListSize = getNeighborListSize()
251 nlist = 0
252
253 do i = 1, nrow
254 point(i) = nlist + 1
255
256 inner: do j = 1, ncol
257
258 if (skipThisPair(i,j)) cycle inner
259
260 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
261
262 if (rijsq < rlistsq) then
263
264 nlist = nlist + 1
265
266 if (nlist > neighborListSize) then
267 call expandNeighborList(nlocal, listerror)
268 if (listerror /= 0) then
269 error = -1
270 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
271 return
272 end if
273 endif
274
275 list(nlist) = j
276
277 if (rijsq < rcutsq) then
278 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
279 u_l, A, f, t,pot)
280 endif
281 endif
282 enddo inner
283 enddo
284
285 point(nrow + 1) = nlist + 1
286
287 else !! (of update_check)
288
289 ! 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
296 do jnab = jbeg, jend
297 j = list(jnab)
298
299 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
300 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
301 u_l, A, f, t,pot)
302
303 enddo
304 endif
305 enddo
306 endif
307
308 #else
309
310 if (update_nlist) then
311
312 ! save current configuration, contruct neighbor list,
313 ! and calculate forces
314 call saveNeighborList(q)
315
316 neighborListSize = getNeighborListSize()
317 nlist = 0
318
319 do i = 1, natoms-1
320 point(i) = nlist + 1
321
322 inner: do j = i+1, natoms
323
324 if (skipThisPair(i,j)) cycle inner
325
326 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
327
328 if (rijsq < rlistsq) then
329
330 nlist = nlist + 1
331
332 if (nlist > neighborListSize) then
333 write(*,*) "do_Forces: Calling expand nlist w/ ",nlist , size(list)
334 call expandNeighborList(natoms, listerror)
335 if (listerror /= 0) then
336 error = -1
337 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
338 return
339 end if
340 endif
341
342 list(nlist) = j
343
344 if (rijsq < rcutsq) then
345 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
346 u_l, A, f, t,pot)
347 endif
348 endif
349 enddo inner
350 enddo
351
352 point(natoms) = nlist + 1
353
354 else !! (update)
355
356 ! use the list to find the neighbors
357 do i = 1, natoms-1
358 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
366 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
367 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
368 u_l, A, f, t,pot)
369
370 enddo
371 endif
372 enddo
373 endif
374
375 #endif
376
377 ! phew, done with main loop.
378
379 #ifdef IS_MPI
380 !!distribute forces
381
382 call scatter(f_Row,f,plan_row3d)
383 call scatter(f_Col,f_temp,plan_col3d)
384 do i = 1,nlocal
385 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
386 end do
387
388 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
389 call scatter(t_Row,t,plan_row3d)
390 call scatter(t_Col,t_temp,plan_col3d)
391
392 do i = 1,nlocal
393 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
394 end do
395 endif
396
397 if (do_pot) then
398 ! scatter/gather pot_row into the members of my column
399 call scatter(pot_Row, pot_Temp, plan_row)
400
401 ! scatter/gather pot_local into all other procs
402 ! add resultant to get total pot
403 do i = 1, nlocal
404 pot_local = pot_local + pot_Temp(i)
405 enddo
406
407 pot_Temp = 0.0_DP
408
409 call scatter(pot_Col, pot_Temp, plan_col)
410 do i = 1, nlocal
411 pot_local = pot_local + pot_Temp(i)
412 enddo
413
414 endif
415 #endif
416
417 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 rfpot = 0.0_DP
432 #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 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
446 #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 pot = pot_local
461 !! 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 endif
464
465 if (do_stress) then
466 call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
467 mpi_comm_world,mpi_err)
468 call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
469 mpi_comm_world,mpi_err)
470 endif
471
472 #else
473
474 if (do_stress) then
475 tau = tau_Temp
476 virial = virial_Temp
477 endif
478
479 #endif
480
481 end subroutine do_force_loop
482
483 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
484
485 real( kind = dp ) :: pot
486 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
491 logical, intent(inout) :: do_pot, do_stress
492 integer, intent(in) :: i, j
493 real ( kind = dp ), intent(inout) :: rijsq
494 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 logical :: is_GB_i, is_GB_j
499 logical :: is_Sticky_i, is_Sticky_j
500 integer :: me_i, me_j
501
502 r = sqrt(rijsq)
503
504 #ifdef IS_MPI
505
506 me_i = atid_row(i)
507 me_j = atid_col(j)
508
509 #else
510
511 me_i = atid(i)
512 me_j = atid(j)
513
514 #endif
515
516 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 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 endif
522
523 if (FF_uses_dipoles .and. SimUsesDipoles()) then
524 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 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
530 do_pot, do_stress)
531
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 endif
540 endif
541
542 if (FF_uses_Sticky .and. SimUsesSticky()) then
543
544 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 endif
552
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
565 end subroutine do_pair
566
567
568 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
569
570 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 real( kind = dp ) :: d_old(3)
575 d(1:3) = q_i(1:3) - q_j(1:3)
576 d_old = d
577 ! Wrap back into periodic box if necessary
578 if ( SimUsesPBC() ) then
579
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 endif
584 r_sq = dot_product(d,d)
585
586 end subroutine get_interatomic_vector
587
588 subroutine check_initialization(error)
589 integer, intent(out) :: error
590
591 error = 0
592 ! Make sure we are properly initialized.
593 if (.not. do_forces_initialized) then
594 error = -1
595 return
596 endif
597
598 #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
606 return
607 end subroutine check_initialization
608
609
610 subroutine zero_work_arrays()
611
612 #ifdef IS_MPI
613
614 q_Row = 0.0_dp
615 q_Col = 0.0_dp
616
617 u_l_Row = 0.0_dp
618 u_l_Col = 0.0_dp
619
620 A_Row = 0.0_dp
621 A_Col = 0.0_dp
622
623 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
631 pot_Row = 0.0_dp
632 pot_Col = 0.0_dp
633 pot_Temp = 0.0_dp
634
635 rf_Row = 0.0_dp
636 rf_Col = 0.0_dp
637 rf_Temp = 0.0_dp
638
639 #endif
640
641 rf = 0.0_dp
642 tau_Temp = 0.0_dp
643 virial_Temp = 0.0_dp
644
645 end subroutine zero_work_arrays
646
647 function skipThisPair(atom1, atom2) result(skip_it)
648
649 integer, intent(in) :: atom1
650 integer, intent(in), optional :: atom2
651 logical :: skip_it
652 integer :: unique_id_1, unique_id_2
653 integer :: i
654
655 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
663 #ifdef IS_MPI
664 !! in MPI, we have to look up the unique IDs for each atom
665 unique_id_1 = tagRow(atom1)
666 #else
667 !! in the normal loop, the atom numbers are unique
668 unique_id_1 = atom1
669 #endif
670
671 !! We were called with only one atom, so just check the global exclude
672 !! list for this atom
673 if (.not. present(atom2)) then
674 do i = 1, nExcludes_global
675 if (excludesGlobal(i) == unique_id_1) then
676 skip_it = .true.
677 return
678 end if
679 end do
680 return
681 end if
682
683 #ifdef IS_MPI
684 unique_id_2 = tagColumn(atom2)
685 #else
686 unique_id_2 = atom2
687 #endif
688
689 #ifdef IS_MPI
690 !! this situation should only arise in MPI simulations
691 if (unique_id_1 == unique_id_2) then
692 skip_it = .true.
693 return
694 end if
695
696 !! this prevents us from doing the pair on multiple processors
697 if (unique_id_1 < unique_id_2) then
698 if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
699 return
700 else
701 if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
702 endif
703 #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
714 do i = 1, nExcludes_local
715 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 end do
729
730 return
731 end function skipThisPair
732
733 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 end module do_Forces