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