ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 377
Committed: Fri Mar 21 17:42:12 2003 UTC (21 years, 3 months ago) by mmeineke
Original Path: branches/mmeineke/OOPSE/libmdtools/do_Forces.F90
File size: 20558 byte(s)
Log Message:
New OOPSE Tree

File Contents

# User Rev Content
1 mmeineke 377 !! 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.1.1.1 2003-03-21 17:42:12 mmeineke Exp $, $Date: 2003-03-21 17:42:12 $, $Name: not supported by cvs2svn $, $Revision: 1.1.1.1 $
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 = size(list)
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     neighborListSize = size(list)
274     endif
275    
276     list(nlist) = j
277    
278     if (rijsq < rcutsq) then
279     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
280     u_l, A, f, t,pot)
281     endif
282     endif
283     enddo inner
284     enddo
285    
286     point(nrow + 1) = nlist + 1
287    
288     else !! (of update_check)
289    
290     ! use the list to find the neighbors
291     do i = 1, nrow
292     JBEG = POINT(i)
293     JEND = POINT(i+1) - 1
294     ! check thiat molecule i has neighbors
295     if (jbeg .le. jend) then
296    
297     do jnab = jbeg, jend
298     j = list(jnab)
299    
300     call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
301     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
302     u_l, A, f, t,pot)
303    
304     enddo
305     endif
306     enddo
307     endif
308    
309     #else
310    
311     if (update_nlist) then
312    
313     ! save current configuration, contruct neighbor list,
314     ! and calculate forces
315     call saveNeighborList(q)
316    
317     neighborListSize = size(list)
318    
319     nlist = 0
320    
321     do i = 1, natoms-1
322     point(i) = nlist + 1
323    
324     inner: do j = i+1, natoms
325    
326     if (skipThisPair(i,j)) cycle inner
327    
328     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
329    
330     if (rijsq < rlistsq) then
331    
332     nlist = nlist + 1
333    
334     if (nlist > neighborListSize) then
335     call expandNeighborList(natoms, listerror)
336     if (listerror /= 0) then
337     error = -1
338     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
339     return
340     end if
341     neighborListSize = size(list)
342     endif
343    
344     list(nlist) = j
345    
346     if (rijsq < rcutsq) then
347     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
348     u_l, A, f, t,pot)
349     endif
350     endif
351     enddo inner
352     enddo
353    
354     point(natoms) = nlist + 1
355    
356     else !! (update)
357    
358     ! use the list to find the neighbors
359     do i = 1, natoms-1
360     JBEG = POINT(i)
361     JEND = POINT(i+1) - 1
362     ! check thiat molecule i has neighbors
363     if (jbeg .le. jend) then
364    
365     do jnab = jbeg, jend
366     j = list(jnab)
367    
368     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
369     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
370     u_l, A, f, t,pot)
371    
372     enddo
373     endif
374     enddo
375     endif
376    
377     #endif
378    
379     ! phew, done with main loop.
380    
381     #ifdef IS_MPI
382     !!distribute forces
383    
384     call scatter(f_Row,f,plan_row3d)
385     call scatter(f_Col,f_temp,plan_col3d)
386     do i = 1,nlocal
387     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
388     end do
389    
390     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
391     call scatter(t_Row,t,plan_row3d)
392     call scatter(t_Col,t_temp,plan_col3d)
393    
394     do i = 1,nlocal
395     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
396     end do
397     endif
398    
399     if (do_pot) then
400     ! scatter/gather pot_row into the members of my column
401     call scatter(pot_Row, pot_Temp, plan_row)
402    
403     ! scatter/gather pot_local into all other procs
404     ! add resultant to get total pot
405     do i = 1, nlocal
406     pot_local = pot_local + pot_Temp(i)
407     enddo
408    
409     pot_Temp = 0.0_DP
410    
411     call scatter(pot_Col, pot_Temp, plan_col)
412     do i = 1, nlocal
413     pot_local = pot_local + pot_Temp(i)
414     enddo
415    
416     endif
417     #endif
418    
419     if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
420    
421     if (FF_uses_RF .and. SimUsesRF()) then
422    
423     #ifdef IS_MPI
424     call scatter(rf_Row,rf,plan_row3d)
425     call scatter(rf_Col,rf_Temp,plan_col3d)
426     do i = 1,nlocal
427     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
428     end do
429     #endif
430    
431     do i = 1, getNlocal()
432    
433     rfpot = 0.0_DP
434     #ifdef IS_MPI
435     me_i = atid_row(i)
436     #else
437     me_i = atid(i)
438     #endif
439     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
440     if ( is_DP_i ) then
441     call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
442     !! The reaction field needs to include a self contribution
443     !! to the field:
444     call accumulate_self_rf(i, mu_i, u_l)
445     !! Get the reaction field contribution to the
446     !! potential and torques:
447     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
448     #ifdef IS_MPI
449     pot_local = pot_local + rfpot
450     #else
451     pot = pot + rfpot
452    
453     #endif
454     endif
455     enddo
456     endif
457     endif
458    
459    
460     #ifdef IS_MPI
461    
462     if (do_pot) then
463     pot = pot_local
464     !! we assume the c code will do the allreduce to get the total potential
465     !! we could do it right here if we needed to...
466     endif
467    
468     if (do_stress) then
469     call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
470     mpi_comm_world,mpi_err)
471     call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
472     mpi_comm_world,mpi_err)
473     endif
474    
475     #else
476    
477     if (do_stress) then
478     tau = tau_Temp
479     virial = virial_Temp
480     endif
481    
482     #endif
483    
484     end subroutine do_force_loop
485    
486     subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t,pot)
487    
488     real( kind = dp ) :: pot
489     real( kind = dp ), dimension(:,:) :: u_l
490     real (kind=dp), dimension(:,:) :: A
491     real (kind=dp), dimension(:,:) :: f
492     real (kind=dp), dimension(:,:) :: t
493    
494     logical, intent(inout) :: do_pot, do_stress
495     integer, intent(in) :: i, j
496     real ( kind = dp ), intent(inout) :: rijsq
497     real ( kind = dp ) :: r
498     real ( kind = dp ), intent(inout) :: d(3)
499     logical :: is_LJ_i, is_LJ_j
500     logical :: is_DP_i, is_DP_j
501     logical :: is_GB_i, is_GB_j
502     logical :: is_Sticky_i, is_Sticky_j
503     integer :: me_i, me_j
504    
505     r = sqrt(rijsq)
506    
507     #ifdef IS_MPI
508    
509     me_i = atid_row(i)
510     me_j = atid_col(j)
511    
512     #else
513    
514     me_i = atid(i)
515     me_j = atid(j)
516    
517     #endif
518    
519     if (FF_uses_LJ .and. SimUsesLJ()) then
520     call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
521     call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
522    
523     if ( is_LJ_i .and. is_LJ_j ) &
524     call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
525     endif
526    
527     if (FF_uses_dipoles .and. SimUsesDipoles()) then
528     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
529     call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
530    
531     if ( is_DP_i .and. is_DP_j ) then
532    
533     call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
534     do_pot, do_stress)
535    
536     if (FF_uses_RF .and. SimUsesRF()) then
537    
538     call accumulate_rf(i, j, r, u_l)
539     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
540    
541     endif
542    
543     endif
544     endif
545    
546     if (FF_uses_Sticky .and. SimUsesSticky()) then
547    
548     call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
549     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
550    
551     if ( is_Sticky_i .and. is_Sticky_j ) then
552     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
553     do_pot, do_stress)
554     endif
555     endif
556    
557    
558     if (FF_uses_GB .and. SimUsesGB()) then
559    
560     call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
561     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
562    
563     if ( is_GB_i .and. is_GB_j ) then
564     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
565     do_pot, do_stress)
566     endif
567     endif
568    
569     end subroutine do_pair
570    
571    
572     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
573    
574     real (kind = dp), dimension(3) :: q_i
575     real (kind = dp), dimension(3) :: q_j
576     real ( kind = dp ), intent(out) :: r_sq
577     real( kind = dp ) :: d(3)
578     real( kind = dp ) :: d_old(3)
579     d(1:3) = q_i(1:3) - q_j(1:3)
580     d_old = d
581     ! Wrap back into periodic box if necessary
582     if ( SimUsesPBC() ) then
583    
584     d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
585     int(abs(d(1:3)/box(1:3)) + 0.5_dp)
586    
587     endif
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     error = -1
599     return
600     endif
601    
602     #ifdef IS_MPI
603     if (.not. isMPISimSet()) then
604     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
605     error = -1
606     return
607     endif
608     #endif
609    
610     return
611     end subroutine check_initialization
612    
613    
614     subroutine zero_work_arrays()
615    
616     #ifdef IS_MPI
617    
618     q_Row = 0.0_dp
619     q_Col = 0.0_dp
620    
621     u_l_Row = 0.0_dp
622     u_l_Col = 0.0_dp
623    
624     A_Row = 0.0_dp
625     A_Col = 0.0_dp
626    
627     f_Row = 0.0_dp
628     f_Col = 0.0_dp
629     f_Temp = 0.0_dp
630    
631     t_Row = 0.0_dp
632     t_Col = 0.0_dp
633     t_Temp = 0.0_dp
634    
635     pot_Row = 0.0_dp
636     pot_Col = 0.0_dp
637     pot_Temp = 0.0_dp
638    
639     rf_Row = 0.0_dp
640     rf_Col = 0.0_dp
641     rf_Temp = 0.0_dp
642    
643     #endif
644    
645     rf = 0.0_dp
646     tau_Temp = 0.0_dp
647     virial_Temp = 0.0_dp
648    
649     end subroutine zero_work_arrays
650    
651     function skipThisPair(atom1, atom2) result(skip_it)
652    
653     integer, intent(in) :: atom1
654     integer, intent(in), optional :: atom2
655     logical :: skip_it
656     integer :: unique_id_1, unique_id_2
657     integer :: i
658    
659     skip_it = .false.
660    
661     !! there are a number of reasons to skip a pair or a particle
662     !! mostly we do this to exclude atoms who are involved in short
663     !! range interactions (bonds, bends, torsions), but we also need
664     !! to exclude some overcounted interactions that result from
665     !! the parallel decomposition
666    
667     #ifdef IS_MPI
668     !! in MPI, we have to look up the unique IDs for each atom
669     unique_id_1 = tagRow(atom1)
670     #else
671     !! in the normal loop, the atom numbers are unique
672     unique_id_1 = atom1
673     #endif
674    
675     !! We were called with only one atom, so just check the global exclude
676     !! list for this atom
677     if (.not. present(atom2)) then
678     do i = 1, nExcludes_global
679     if (excludesGlobal(i) == unique_id_1) then
680     skip_it = .true.
681     return
682     end if
683     end do
684     return
685     end if
686    
687     #ifdef IS_MPI
688     unique_id_2 = tagColumn(atom2)
689     #else
690     unique_id_2 = atom2
691     #endif
692    
693     #ifdef IS_MPI
694     !! this situation should only arise in MPI simulations
695     if (unique_id_1 == unique_id_2) then
696     skip_it = .true.
697     return
698     end if
699    
700     !! this prevents us from doing the pair on multiple processors
701     if (unique_id_1 < unique_id_2) then
702     if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
703     return
704     else
705     if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
706     return
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