ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 393
Committed: Mon Mar 24 18:33:51 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 20667 byte(s)
Log Message:
little bug fixes here and there

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.3 2003-03-24 18:33:51 mmeineke Exp $, $Date: 2003-03-24 18:33:51 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
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 write(*,*) 'u_l', u_l(1,1), u_l(2,1), u_l(3,1)
201
202
203 call getRcut(rcut,rc2=rcutsq)
204 call getRlist(rlist,rlistsq)
205
206 call check_initialization(localError)
207 if ( localError .ne. 0 ) then
208 error = -1
209 return
210 end if
211 call zero_work_arrays()
212
213 do_pot = do_pot_c
214 do_stress = do_stress_c
215
216 ! Gather all information needed by all force loops:
217
218 #ifdef IS_MPI
219
220 call gather(q,q_Row,plan_row3d)
221 call gather(q,q_Col,plan_col3d)
222
223 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
224 call gather(u_l,u_l_Row,plan_row3d)
225 call gather(u_l,u_l_Col,plan_col3d)
226
227 call gather(A,A_Row,plan_row_rotation)
228 call gather(A,A_Col,plan_col_rotation)
229 endif
230
231 #endif
232
233 if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
234 !! See if we need to update neighbor lists
235 call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
236 !! if_mpi_gather_stuff_for_prepair
237 !! do_prepair_loop_if_needed
238 !! if_mpi_scatter_stuff_from_prepair
239 !! if_mpi_gather_stuff_from_prepair_to_main_loop
240 else
241 !! See if we need to update neighbor lists
242 call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)
243 endif
244
245 #ifdef IS_MPI
246
247 if (update_nlist) then
248
249 !! save current configuration, construct neighbor list,
250 !! and calculate forces
251 call saveNeighborList(q)
252
253 neighborListSize = size(list)
254 nlist = 0
255
256 do i = 1, nrow
257 point(i) = nlist + 1
258
259 inner: do j = 1, ncol
260
261 if (skipThisPair(i,j)) cycle inner
262
263 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
264
265 if (rijsq < rlistsq) then
266
267 nlist = nlist + 1
268
269 if (nlist > neighborListSize) then
270 call expandNeighborList(nlocal, listerror)
271 if (listerror /= 0) then
272 error = -1
273 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
274 return
275 end if
276 neighborListSize = size(list)
277 endif
278
279 list(nlist) = j
280
281 if (rijsq < rcutsq) then
282 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
283 u_l, A, f, t,pot)
284 endif
285 endif
286 enddo inner
287 enddo
288
289 point(nrow + 1) = nlist + 1
290
291 else !! (of update_check)
292
293 ! use the list to find the neighbors
294 do i = 1, nrow
295 JBEG = POINT(i)
296 JEND = POINT(i+1) - 1
297 ! check thiat molecule i has neighbors
298 if (jbeg .le. jend) then
299
300 do jnab = jbeg, jend
301 j = list(jnab)
302
303 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
304 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
305 u_l, A, f, t,pot)
306
307 enddo
308 endif
309 enddo
310 endif
311
312 #else
313
314 if (update_nlist) then
315
316 ! save current configuration, contruct neighbor list,
317 ! and calculate forces
318 call saveNeighborList(q)
319
320 neighborListSize = size(list)
321
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
334 if (rijsq < rlistsq) then
335
336 nlist = nlist + 1
337
338 if (nlist > neighborListSize) then
339 call expandNeighborList(natoms, listerror)
340 if (listerror /= 0) then
341 error = -1
342 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
343 return
344 end if
345 neighborListSize = size(list)
346 endif
347
348 list(nlist) = j
349
350 if (rijsq < rcutsq) then
351 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
352 u_l, A, f, t,pot)
353 endif
354 endif
355 enddo inner
356 enddo
357
358 point(natoms) = nlist + 1
359
360 else !! (update)
361
362 ! use the list to find the neighbors
363 do i = 1, natoms-1
364 JBEG = POINT(i)
365 JEND = POINT(i+1) - 1
366 ! check thiat molecule i has neighbors
367 if (jbeg .le. jend) then
368
369 do jnab = jbeg, jend
370 j = list(jnab)
371
372 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
373 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
374 u_l, A, f, t,pot)
375
376 enddo
377 endif
378 enddo
379 endif
380
381 #endif
382
383 ! phew, done with main loop.
384
385 #ifdef IS_MPI
386 !!distribute forces
387
388 call scatter(f_Row,f,plan_row3d)
389 call scatter(f_Col,f_temp,plan_col3d)
390 do i = 1,nlocal
391 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
392 end do
393
394 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
395 call scatter(t_Row,t,plan_row3d)
396 call scatter(t_Col,t_temp,plan_col3d)
397
398 do i = 1,nlocal
399 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
400 end do
401 endif
402
403 if (do_pot) then
404 ! scatter/gather pot_row into the members of my column
405 call scatter(pot_Row, pot_Temp, plan_row)
406
407 ! scatter/gather pot_local into all other procs
408 ! add resultant to get total pot
409 do i = 1, nlocal
410 pot_local = pot_local + pot_Temp(i)
411 enddo
412
413 pot_Temp = 0.0_DP
414
415 call scatter(pot_Col, pot_Temp, plan_col)
416 do i = 1, nlocal
417 pot_local = pot_local + pot_Temp(i)
418 enddo
419
420 endif
421 #endif
422
423 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
424
425 if (FF_uses_RF .and. SimUsesRF()) then
426
427 #ifdef IS_MPI
428 call scatter(rf_Row,rf,plan_row3d)
429 call scatter(rf_Col,rf_Temp,plan_col3d)
430 do i = 1,nlocal
431 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
432 end do
433 #endif
434
435 do i = 1, getNlocal()
436
437 rfpot = 0.0_DP
438 #ifdef IS_MPI
439 me_i = atid_row(i)
440 #else
441 me_i = atid(i)
442 #endif
443 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
444 if ( is_DP_i ) then
445 call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
446 !! The reaction field needs to include a self contribution
447 !! to the field:
448 call accumulate_self_rf(i, mu_i, u_l)
449 !! Get the reaction field contribution to the
450 !! potential and torques:
451 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
452 #ifdef IS_MPI
453 pot_local = pot_local + rfpot
454 #else
455 pot = pot + rfpot
456
457 #endif
458 endif
459 enddo
460 endif
461 endif
462
463
464 #ifdef IS_MPI
465
466 if (do_pot) then
467 pot = pot_local
468 !! we assume the c code will do the allreduce to get the total potential
469 !! we could do it right here if we needed to...
470 endif
471
472 if (do_stress) then
473 call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
474 mpi_comm_world,mpi_err)
475 call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
476 mpi_comm_world,mpi_err)
477 endif
478
479 #else
480
481 if (do_stress) then
482 tau = tau_Temp
483 virial = virial_Temp
484 endif
485
486 #endif
487
488 end subroutine do_force_loop
489
490 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
491
492 real( kind = dp ) :: pot
493 real( kind = dp ), dimension(:,:) :: u_l
494 real (kind=dp), dimension(:,:) :: A
495 real (kind=dp), dimension(:,:) :: f
496 real (kind=dp), dimension(:,:) :: t
497
498 logical, intent(inout) :: do_pot, do_stress
499 integer, intent(in) :: i, j
500 real ( kind = dp ), intent(inout) :: rijsq
501 real ( kind = dp ) :: r
502 real ( kind = dp ), intent(inout) :: d(3)
503 logical :: is_LJ_i, is_LJ_j
504 logical :: is_DP_i, is_DP_j
505 logical :: is_GB_i, is_GB_j
506 logical :: is_Sticky_i, is_Sticky_j
507 integer :: me_i, me_j
508
509 r = sqrt(rijsq)
510
511 #ifdef IS_MPI
512
513 me_i = atid_row(i)
514 me_j = atid_col(j)
515
516 #else
517
518 me_i = atid(i)
519 me_j = atid(j)
520
521 #endif
522
523 if (FF_uses_LJ .and. SimUsesLJ()) then
524 call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
525 call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
526
527 if ( is_LJ_i .and. is_LJ_j ) &
528 call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
529 endif
530
531 if (FF_uses_dipoles .and. SimUsesDipoles()) then
532 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
533 call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
534
535 if ( is_DP_i .and. is_DP_j ) then
536
537 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
538 do_pot, do_stress)
539 if (FF_uses_RF .and. SimUsesRF()) then
540 call accumulate_rf(i, j, r, u_l)
541 call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
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 real( kind = dp ) :: d_old(3)
580 d(1:3) = q_i(1:3) - q_j(1:3)
581 d_old = d
582 ! Wrap back into periodic box if necessary
583 if ( SimUsesPBC() ) then
584
585 d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
586 int(abs(d(1:3)/box(1:3)) + 0.5_dp)
587
588 ! write(*,'(a11,3es12.3)') 'wrapped to ', d(1), d(2), d(3)
589
590 endif
591 r_sq = dot_product(d,d)
592
593 end subroutine get_interatomic_vector
594
595 subroutine check_initialization(error)
596 integer, intent(out) :: error
597
598 error = 0
599 ! Make sure we are properly initialized.
600 if (.not. do_forces_initialized) then
601 error = -1
602 return
603 endif
604
605 #ifdef IS_MPI
606 if (.not. isMPISimSet()) then
607 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
608 error = -1
609 return
610 endif
611 #endif
612
613 return
614 end subroutine check_initialization
615
616
617 subroutine zero_work_arrays()
618
619 #ifdef IS_MPI
620
621 q_Row = 0.0_dp
622 q_Col = 0.0_dp
623
624 u_l_Row = 0.0_dp
625 u_l_Col = 0.0_dp
626
627 A_Row = 0.0_dp
628 A_Col = 0.0_dp
629
630 f_Row = 0.0_dp
631 f_Col = 0.0_dp
632 f_Temp = 0.0_dp
633
634 t_Row = 0.0_dp
635 t_Col = 0.0_dp
636 t_Temp = 0.0_dp
637
638 pot_Row = 0.0_dp
639 pot_Col = 0.0_dp
640 pot_Temp = 0.0_dp
641
642 rf_Row = 0.0_dp
643 rf_Col = 0.0_dp
644 rf_Temp = 0.0_dp
645
646 #endif
647
648 rf = 0.0_dp
649 tau_Temp = 0.0_dp
650 virial_Temp = 0.0_dp
651
652 end subroutine zero_work_arrays
653
654 function skipThisPair(atom1, atom2) result(skip_it)
655 integer, intent(in) :: atom1
656 integer, intent(in), optional :: atom2
657 logical :: skip_it
658 integer :: unique_id_1, unique_id_2
659 integer :: me_i,me_j
660 integer :: i
661
662 skip_it = .false.
663
664 !! there are a number of reasons to skip a pair or a particle
665 !! mostly we do this to exclude atoms who are involved in short
666 !! range interactions (bonds, bends, torsions), but we also need
667 !! to exclude some overcounted interactions that result from
668 !! the parallel decomposition
669
670 #ifdef IS_MPI
671 !! in MPI, we have to look up the unique IDs for each atom
672 unique_id_1 = tagRow(atom1)
673 #else
674 !! in the normal loop, the atom numbers are unique
675 unique_id_1 = atom1
676 #endif
677
678 !! We were called with only one atom, so just check the global exclude
679 !! list for this atom
680 if (.not. present(atom2)) then
681 do i = 1, nExcludes_global
682 if (excludesGlobal(i) == unique_id_1) then
683 skip_it = .true.
684 return
685 end if
686 end do
687 return
688 end if
689
690 #ifdef IS_MPI
691 unique_id_2 = tagColumn(atom2)
692 #else
693 unique_id_2 = atom2
694 #endif
695
696 #ifdef IS_MPI
697 !! this situation should only arise in MPI simulations
698 if (unique_id_1 == unique_id_2) then
699 skip_it = .true.
700 return
701 end if
702
703 !! this prevents us from doing the pair on multiple processors
704 if (unique_id_1 < unique_id_2) then
705 if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
706 return
707 else
708 if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
709 return
710 endif
711 #endif
712
713 !! the rest of these situations can happen in all simulations:
714 do i = 1, nExcludes_global
715 if ((excludesGlobal(i) == unique_id_1) .or. &
716 (excludesGlobal(i) == unique_id_2)) then
717 skip_it = .true.
718 return
719 endif
720 enddo
721
722 do i = 1, nExcludes_local
723 if (excludesLocal(1,i) == unique_id_1) then
724 if (excludesLocal(2,i) == unique_id_2) then
725 skip_it = .true.
726 return
727 endif
728 else
729 if (excludesLocal(1,i) == unique_id_2) then
730 if (excludesLocal(2,i) == unique_id_1) then
731 skip_it = .true.
732 return
733 endif
734 endif
735 endif
736 end do
737
738 return
739 end function skipThisPair
740
741 function FF_UsesDirectionalAtoms() result(doesit)
742 logical :: doesit
743 doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
744 FF_uses_GB .or. FF_uses_RF
745 end function FF_UsesDirectionalAtoms
746
747 function FF_RequiresPrepairCalc() result(doesit)
748 logical :: doesit
749 doesit = FF_uses_EAM
750 end function FF_RequiresPrepairCalc
751
752 function FF_RequiresPostpairCalc() result(doesit)
753 logical :: doesit
754 doesit = FF_uses_RF
755 end function FF_RequiresPostpairCalc
756
757 end module do_Forces