ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 388
Committed: Fri Mar 21 22:11:50 2003 UTC (21 years, 3 months ago) by chuckv
File size: 20532 byte(s)
Log Message:
various write statements for debugging

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 chuckv 388 !! @version $Id: do_Forces.F90,v 1.2 2003-03-21 22:11:50 chuckv Exp $, $Date: 2003-03-21 22:11:50 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $
8 mmeineke 377
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 chuckv 388 if (skipThisPair(i,j)) cycle inner
327    
328 mmeineke 377 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
329    
330 chuckv 388
331 mmeineke 377 if (rijsq < rlistsq) then
332    
333     nlist = nlist + 1
334    
335     if (nlist > neighborListSize) then
336     call expandNeighborList(natoms, listerror)
337     if (listerror /= 0) then
338     error = -1
339     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
340     return
341     end if
342     neighborListSize = size(list)
343     endif
344    
345     list(nlist) = j
346    
347     if (rijsq < rcutsq) then
348     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
349     u_l, A, f, t,pot)
350     endif
351     endif
352     enddo inner
353     enddo
354    
355     point(natoms) = nlist + 1
356    
357     else !! (update)
358    
359     ! use the list to find the neighbors
360     do i = 1, natoms-1
361     JBEG = POINT(i)
362     JEND = POINT(i+1) - 1
363     ! check thiat molecule i has neighbors
364     if (jbeg .le. jend) then
365    
366     do jnab = jbeg, jend
367     j = list(jnab)
368    
369     call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
370     call do_pair(i, j, rijsq, d, do_pot, do_stress, &
371     u_l, A, f, t,pot)
372    
373     enddo
374     endif
375     enddo
376     endif
377    
378     #endif
379    
380     ! phew, done with main loop.
381    
382     #ifdef IS_MPI
383     !!distribute forces
384    
385     call scatter(f_Row,f,plan_row3d)
386     call scatter(f_Col,f_temp,plan_col3d)
387     do i = 1,nlocal
388     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
389     end do
390    
391     if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
392     call scatter(t_Row,t,plan_row3d)
393     call scatter(t_Col,t_temp,plan_col3d)
394    
395     do i = 1,nlocal
396     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
397     end do
398     endif
399    
400     if (do_pot) then
401     ! scatter/gather pot_row into the members of my column
402     call scatter(pot_Row, pot_Temp, plan_row)
403    
404     ! scatter/gather pot_local into all other procs
405     ! add resultant to get total pot
406     do i = 1, nlocal
407     pot_local = pot_local + pot_Temp(i)
408     enddo
409    
410     pot_Temp = 0.0_DP
411    
412     call scatter(pot_Col, pot_Temp, plan_col)
413     do i = 1, nlocal
414     pot_local = pot_local + pot_Temp(i)
415     enddo
416    
417     endif
418     #endif
419    
420     if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
421    
422     if (FF_uses_RF .and. SimUsesRF()) then
423    
424     #ifdef IS_MPI
425     call scatter(rf_Row,rf,plan_row3d)
426     call scatter(rf_Col,rf_Temp,plan_col3d)
427     do i = 1,nlocal
428     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
429     end do
430     #endif
431    
432     do i = 1, getNlocal()
433    
434     rfpot = 0.0_DP
435     #ifdef IS_MPI
436     me_i = atid_row(i)
437     #else
438     me_i = atid(i)
439     #endif
440     call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
441     if ( is_DP_i ) then
442     call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
443     !! The reaction field needs to include a self contribution
444     !! to the field:
445     call accumulate_self_rf(i, mu_i, u_l)
446     !! Get the reaction field contribution to the
447     !! potential and torques:
448     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
449     #ifdef IS_MPI
450     pot_local = pot_local + rfpot
451     #else
452     pot = pot + rfpot
453    
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,pot)
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     if (FF_uses_RF .and. SimUsesRF()) then
537     call accumulate_rf(i, j, r, u_l)
538     call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
539     endif
540    
541     endif
542     endif
543    
544     if (FF_uses_Sticky .and. SimUsesSticky()) then
545    
546     call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
547     call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
548 chuckv 388
549 mmeineke 377 if ( is_Sticky_i .and. is_Sticky_j ) then
550     call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
551     do_pot, do_stress)
552     endif
553     endif
554    
555    
556     if (FF_uses_GB .and. SimUsesGB()) then
557    
558     call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
559     call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
560    
561     if ( is_GB_i .and. is_GB_j ) then
562     call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
563     do_pot, do_stress)
564     endif
565     endif
566    
567     end subroutine do_pair
568    
569    
570     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
571    
572     real (kind = dp), dimension(3) :: q_i
573     real (kind = dp), dimension(3) :: q_j
574     real ( kind = dp ), intent(out) :: r_sq
575     real( kind = dp ) :: d(3)
576     real( kind = dp ) :: d_old(3)
577     d(1:3) = q_i(1:3) - q_j(1:3)
578     d_old = d
579     ! Wrap back into periodic box if necessary
580     if ( SimUsesPBC() ) then
581    
582     d(1:3) = d(1:3) - box(1:3) * sign(1.0_dp,d(1:3)) * &
583     int(abs(d(1:3)/box(1:3)) + 0.5_dp)
584    
585     endif
586     r_sq = dot_product(d,d)
587    
588     end subroutine get_interatomic_vector
589    
590     subroutine check_initialization(error)
591     integer, intent(out) :: error
592    
593     error = 0
594     ! Make sure we are properly initialized.
595     if (.not. do_forces_initialized) then
596     error = -1
597     return
598     endif
599    
600     #ifdef IS_MPI
601     if (.not. isMPISimSet()) then
602     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
603     error = -1
604     return
605     endif
606     #endif
607    
608     return
609     end subroutine check_initialization
610    
611    
612     subroutine zero_work_arrays()
613    
614     #ifdef IS_MPI
615    
616     q_Row = 0.0_dp
617     q_Col = 0.0_dp
618    
619     u_l_Row = 0.0_dp
620     u_l_Col = 0.0_dp
621    
622     A_Row = 0.0_dp
623     A_Col = 0.0_dp
624    
625     f_Row = 0.0_dp
626     f_Col = 0.0_dp
627     f_Temp = 0.0_dp
628    
629     t_Row = 0.0_dp
630     t_Col = 0.0_dp
631     t_Temp = 0.0_dp
632    
633     pot_Row = 0.0_dp
634     pot_Col = 0.0_dp
635     pot_Temp = 0.0_dp
636    
637     rf_Row = 0.0_dp
638     rf_Col = 0.0_dp
639     rf_Temp = 0.0_dp
640    
641     #endif
642    
643     rf = 0.0_dp
644     tau_Temp = 0.0_dp
645     virial_Temp = 0.0_dp
646    
647     end subroutine zero_work_arrays
648    
649     function skipThisPair(atom1, atom2) result(skip_it)
650     integer, intent(in) :: atom1
651     integer, intent(in), optional :: atom2
652     logical :: skip_it
653     integer :: unique_id_1, unique_id_2
654 chuckv 388 integer :: me_i,me_j
655 mmeineke 377 integer :: i
656    
657     skip_it = .false.
658    
659     !! there are a number of reasons to skip a pair or a particle
660     !! mostly we do this to exclude atoms who are involved in short
661     !! range interactions (bonds, bends, torsions), but we also need
662     !! to exclude some overcounted interactions that result from
663     !! the parallel decomposition
664    
665     #ifdef IS_MPI
666     !! in MPI, we have to look up the unique IDs for each atom
667     unique_id_1 = tagRow(atom1)
668     #else
669     !! in the normal loop, the atom numbers are unique
670     unique_id_1 = atom1
671     #endif
672 chuckv 388
673 mmeineke 377 !! We were called with only one atom, so just check the global exclude
674     !! list for this atom
675     if (.not. present(atom2)) then
676     do i = 1, nExcludes_global
677     if (excludesGlobal(i) == unique_id_1) then
678     skip_it = .true.
679     return
680     end if
681     end do
682     return
683     end if
684    
685     #ifdef IS_MPI
686     unique_id_2 = tagColumn(atom2)
687     #else
688     unique_id_2 = atom2
689     #endif
690    
691     #ifdef IS_MPI
692     !! this situation should only arise in MPI simulations
693     if (unique_id_1 == unique_id_2) then
694     skip_it = .true.
695     return
696     end if
697    
698     !! this prevents us from doing the pair on multiple processors
699     if (unique_id_1 < unique_id_2) then
700     if (mod(unique_id_1 + unique_id_2,2) == 0) skip_it = .true.
701     return
702     else
703     if (mod(unique_id_1 + unique_id_2,2) == 1) skip_it = .true.
704     return
705     endif
706     #endif
707    
708     !! the rest of these situations can happen in all simulations:
709     do i = 1, nExcludes_global
710     if ((excludesGlobal(i) == unique_id_1) .or. &
711     (excludesGlobal(i) == unique_id_2)) then
712     skip_it = .true.
713     return
714     endif
715     enddo
716 chuckv 388
717 mmeineke 377 do i = 1, nExcludes_local
718     if (excludesLocal(1,i) == unique_id_1) then
719     if (excludesLocal(2,i) == unique_id_2) then
720     skip_it = .true.
721     return
722     endif
723     else
724     if (excludesLocal(1,i) == unique_id_2) then
725     if (excludesLocal(2,i) == unique_id_1) then
726     skip_it = .true.
727     return
728     endif
729     endif
730     endif
731     end do
732    
733     return
734     end function skipThisPair
735    
736     function FF_UsesDirectionalAtoms() result(doesit)
737     logical :: doesit
738     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
739     FF_uses_GB .or. FF_uses_RF
740     end function FF_UsesDirectionalAtoms
741    
742     function FF_RequiresPrepairCalc() result(doesit)
743     logical :: doesit
744     doesit = FF_uses_EAM
745     end function FF_RequiresPrepairCalc
746    
747     function FF_RequiresPostpairCalc() result(doesit)
748     logical :: doesit
749     doesit = FF_uses_RF
750     end function FF_RequiresPostpairCalc
751    
752     end module do_Forces