ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 356
Committed: Mon Mar 17 20:42:57 2003 UTC (21 years, 5 months ago) by gezelter
File size: 19711 byte(s)
Log Message:
fortran-side fixes to play nice with C

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