ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 370
Committed: Thu Mar 20 17:10:43 2003 UTC (21 years, 4 months ago) by mmeineke
File size: 20464 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.28 2003-03-20 17:10:43 mmeineke Exp $, $Date: 2003-03-20 17:10:43 $, $Name: not supported by cvs2svn $, $Revision: 1.28 $
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 call expandNeighborList(natoms, listerror)
334 if (listerror /= 0) then
335 error = -1
336 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
337 return
338 end if
339 endif
340
341 list(nlist) = j
342
343 if (rijsq < rcutsq) then
344 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
345 u_l, A, f, t,pot)
346 endif
347 endif
348 enddo inner
349 enddo
350
351 point(natoms) = nlist + 1
352
353 else !! (update)
354
355 ! use the list to find the neighbors
356 do i = 1, natoms-1
357 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
365 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
366 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
367 u_l, A, f, t,pot)
368
369 enddo
370 endif
371 enddo
372 endif
373
374 #endif
375
376 ! phew, done with main loop.
377
378 #ifdef IS_MPI
379 !!distribute forces
380
381 call scatter(f_Row,f,plan_row3d)
382 call scatter(f_Col,f_temp,plan_col3d)
383 do i = 1,nlocal
384 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
385 end do
386
387 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
388 call scatter(t_Row,t,plan_row3d)
389 call scatter(t_Col,t_temp,plan_col3d)
390
391 do i = 1,nlocal
392 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
393 end do
394 endif
395
396 if (do_pot) then
397 ! scatter/gather pot_row into the members of my column
398 call scatter(pot_Row, pot_Temp, plan_row)
399
400 ! scatter/gather pot_local into all other procs
401 ! add resultant to get total pot
402 do i = 1, nlocal
403 pot_local = pot_local + pot_Temp(i)
404 enddo
405
406 pot_Temp = 0.0_DP
407
408 call scatter(pot_Col, pot_Temp, plan_col)
409 do i = 1, nlocal
410 pot_local = pot_local + pot_Temp(i)
411 enddo
412
413 endif
414 #endif
415
416 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 rfpot = 0.0_DP
431 #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 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
445 #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 pot = pot_local
460 !! 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 endif
463
464 if (do_stress) then
465 call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
466 mpi_comm_world,mpi_err)
467 call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
468 mpi_comm_world,mpi_err)
469 endif
470
471 #else
472
473 if (do_stress) then
474 tau = tau_Temp
475 virial = virial_Temp
476 endif
477
478 #endif
479
480 end subroutine do_force_loop
481
482 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
483
484 real( kind = dp ) :: pot
485 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
490 logical, intent(inout) :: do_pot, do_stress
491 integer, intent(in) :: i, j
492 real ( kind = dp ), intent(inout) :: rijsq
493 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 logical :: is_GB_i, is_GB_j
498 logical :: is_Sticky_i, is_Sticky_j
499 integer :: me_i, me_j
500
501 r = sqrt(rijsq)
502
503 #ifdef IS_MPI
504
505 me_i = atid_row(i)
506 me_j = atid_col(j)
507
508 #else
509
510 me_i = atid(i)
511 me_j = atid(j)
512
513 #endif
514
515 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 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 endif
521
522 if (FF_uses_dipoles .and. SimUsesDipoles()) then
523 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 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
529 do_pot, do_stress)
530
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 endif
539 endif
540
541 if (FF_uses_Sticky .and. SimUsesSticky()) then
542
543 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 endif
551
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
564 end subroutine do_pair
565
566
567 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
568
569 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 real( kind = dp ) :: d_old(3)
574 d(1:3) = q_i(1:3) - q_j(1:3)
575 d_old = d
576 ! Wrap back into periodic box if necessary
577 if ( SimUsesPBC() ) then
578
579 d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
580 int(abs(d(1:3)/box(1:3)) + 0.5_dp)
581
582 endif
583 r_sq = dot_product(d,d)
584
585 end subroutine get_interatomic_vector
586
587 subroutine check_initialization(error)
588 integer, intent(out) :: error
589
590 error = 0
591 ! Make sure we are properly initialized.
592 if (.not. do_forces_initialized) then
593 error = -1
594 return
595 endif
596
597 #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
605 return
606 end subroutine check_initialization
607
608
609 subroutine zero_work_arrays()
610
611 #ifdef IS_MPI
612
613 q_Row = 0.0_dp
614 q_Col = 0.0_dp
615
616 u_l_Row = 0.0_dp
617 u_l_Col = 0.0_dp
618
619 A_Row = 0.0_dp
620 A_Col = 0.0_dp
621
622 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
630 pot_Row = 0.0_dp
631 pot_Col = 0.0_dp
632 pot_Temp = 0.0_dp
633
634 rf_Row = 0.0_dp
635 rf_Col = 0.0_dp
636 rf_Temp = 0.0_dp
637
638 #endif
639
640 rf = 0.0_dp
641 tau_Temp = 0.0_dp
642 virial_Temp = 0.0_dp
643
644 end subroutine zero_work_arrays
645
646 function skipThisPair(atom1, atom2) result(skip_it)
647
648 integer, intent(in) :: atom1
649 integer, intent(in), optional :: atom2
650 logical :: skip_it
651 integer :: unique_id_1, unique_id_2
652 integer :: i
653
654 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
662 #ifdef IS_MPI
663 !! in MPI, we have to look up the unique IDs for each atom
664 unique_id_1 = tagRow(atom1)
665 #else
666 !! in the normal loop, the atom numbers are unique
667 unique_id_1 = atom1
668 #endif
669
670 !! We were called with only one atom, so just check the global exclude
671 !! list for this atom
672 if (.not. present(atom2)) then
673 do i = 1, nExcludes_global
674 if (excludesGlobal(i) == unique_id_1) then
675 skip_it = .true.
676 return
677 end if
678 end do
679 return
680 end if
681
682 #ifdef IS_MPI
683 unique_id_2 = tagColumn(atom2)
684 #else
685 unique_id_2 = atom2
686 #endif
687
688 #ifdef IS_MPI
689 !! this situation should only arise in MPI simulations
690 if (unique_id_1 == unique_id_2) then
691 skip_it = .true.
692 return
693 end if
694
695 !! this prevents us from doing the pair on multiple processors
696 if (unique_id_1 < unique_id_2) then
697 if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
698 return
699 else
700 if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
701 endif
702 #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
713 do i = 1, nExcludes_local
714 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 end do
728
729 return
730 end function skipThisPair
731
732 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 end module do_Forces