ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 330
Committed: Wed Mar 12 23:15:46 2003 UTC (21 years, 4 months ago) by gezelter
File size: 16620 byte(s)
Log Message:
forceGlobals is gone (part3)

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