ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 358
Committed: Mon Mar 17 21:07:50 2003 UTC (21 years, 5 months ago) by gezelter
File size: 19732 byte(s)
Log Message:
bug fixes for new fForceField.h

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