ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 673
Committed: Fri Aug 8 21:22:37 2003 UTC (20 years, 11 months ago) by chuckv
File size: 28630 byte(s)
Log Message:
EAM works...... Neighbor list also works.....

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.28 2003-08-08 21:22:37 chuckv Exp $, $Date: 2003-08-08 21:22:37 $, $Name: not supported by cvs2svn $, $Revision: 1.28 $
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 use vector_class
21 use eam
22 use status
23 #ifdef IS_MPI
24 use mpiSimulation
25 #endif
26
27 implicit none
28 PRIVATE
29
30 #define __FORTRAN90
31 #include "fForceField.h"
32
33 logical, save :: do_forces_initialized = .false., haveRlist = .false.
34 logical, save :: havePolicies = .false.
35 logical, save :: FF_uses_LJ
36 logical, save :: FF_uses_sticky
37 logical, save :: FF_uses_dipoles
38 logical, save :: FF_uses_RF
39 logical, save :: FF_uses_GB
40 logical, save :: FF_uses_EAM
41
42 real(kind=dp), save :: rlist, rlistsq
43
44 public :: init_FF
45 public :: do_force_loop
46 public :: setRlistDF
47
48 contains
49
50 subroutine setRlistDF( this_rlist )
51
52 real(kind=dp) :: this_rlist
53
54 rlist = this_rlist
55 rlistsq = rlist * rlist
56
57 haveRlist = .true.
58 if( havePolicies ) do_forces_initialized = .true.
59
60 end subroutine setRlistDF
61
62 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
63
64 integer, intent(in) :: LJMIXPOLICY
65 logical, intent(in) :: use_RF_c
66
67 integer, intent(out) :: thisStat
68 integer :: my_status, nMatches
69 integer, pointer :: MatchList(:) => null()
70 real(kind=dp) :: rcut, rrf, rt, dielect
71
72 !! assume things are copacetic, unless they aren't
73 thisStat = 0
74
75 !! Fortran's version of a cast:
76 FF_uses_RF = use_RF_c
77
78 !! init_FF is called *after* all of the atom types have been
79 !! defined in atype_module using the new_atype subroutine.
80 !!
81 !! this will scan through the known atypes and figure out what
82 !! interactions are used by the force field.
83
84 FF_uses_LJ = .false.
85 FF_uses_sticky = .false.
86 FF_uses_dipoles = .false.
87 FF_uses_GB = .false.
88 FF_uses_EAM = .false.
89
90 call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
91 if (nMatches .gt. 0) FF_uses_LJ = .true.
92
93 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
94 if (nMatches .gt. 0) FF_uses_dipoles = .true.
95
96 call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
97 MatchList)
98 if (nMatches .gt. 0) FF_uses_Sticky = .true.
99
100 call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
101 if (nMatches .gt. 0) FF_uses_GB = .true.
102
103 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
104 if (nMatches .gt. 0) FF_uses_EAM = .true.
105
106 !! check to make sure the FF_uses_RF setting makes sense
107
108 if (FF_uses_dipoles) then
109 if (FF_uses_RF) then
110 dielect = getDielect()
111 call initialize_rf(dielect)
112 endif
113 else
114 if (FF_uses_RF) then
115 write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
116 thisStat = -1
117 return
118 endif
119 endif
120
121 if (FF_uses_LJ) then
122
123 select case (LJMIXPOLICY)
124 case (LB_MIXING_RULE)
125 call init_lj_FF(LB_MIXING_RULE, my_status)
126 case (EXPLICIT_MIXING_RULE)
127 call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
128 case default
129 write(default_error,*) 'unknown LJ Mixing Policy!'
130 thisStat = -1
131 return
132 end select
133 if (my_status /= 0) then
134 thisStat = -1
135 return
136 end if
137 endif
138
139 if (FF_uses_sticky) then
140 call check_sticky_FF(my_status)
141 if (my_status /= 0) then
142 thisStat = -1
143 return
144 end if
145 endif
146
147
148 if (FF_uses_EAM) then
149 call init_EAM_FF(my_status)
150 if (my_status /= 0) then
151 thisStat = -1
152 return
153 end if
154 endif
155
156
157
158 if (FF_uses_GB) then
159 call check_gb_pair_FF(my_status)
160 if (my_status .ne. 0) then
161 thisStat = -1
162 return
163 endif
164 endif
165
166 if (FF_uses_GB .and. FF_uses_LJ) then
167 endif
168 if (.not. do_forces_initialized) then
169 !! Create neighbor lists
170 call expandNeighborList(getNlocal(), my_status)
171 if (my_Status /= 0) then
172 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
173 thisStat = -1
174 return
175 endif
176 endif
177
178
179 havePolicies = .true.
180 if( haveRlist ) do_forces_initialized = .true.
181
182 end subroutine init_FF
183
184
185 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
186 !------------------------------------------------------------->
187 subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, &
188 error)
189 !! Position array provided by C, dimensioned by getNlocal
190 real ( kind = dp ), dimension(3,getNlocal()) :: q
191 !! Rotation Matrix for each long range particle in simulation.
192 real( kind = dp), dimension(9,getNlocal()) :: A
193 !! Unit vectors for dipoles (lab frame)
194 real( kind = dp ), dimension(3,getNlocal()) :: u_l
195 !! Force array provided by C, dimensioned by getNlocal
196 real ( kind = dp ), dimension(3,getNlocal()) :: f
197 !! Torsion array provided by C, dimensioned by getNlocal
198 real( kind = dp ), dimension(3,getNlocal()) :: t
199 !! Stress Tensor
200 real( kind = dp), dimension(9) :: tau
201 real ( kind = dp ) :: pot
202 logical ( kind = 2) :: do_pot_c, do_stress_c
203 logical :: do_pot
204 logical :: do_stress
205 #ifdef IS_MPI
206 real( kind = DP ) :: pot_local
207 integer :: nrow
208 integer :: ncol
209 #endif
210 integer :: nlocal
211 integer :: natoms
212 logical :: update_nlist
213 integer :: i, j, jbeg, jend, jnab
214 integer :: nlist
215 real( kind = DP ) :: rijsq
216 real(kind=dp),dimension(3) :: d
217 real(kind=dp) :: rfpot, mu_i, virial
218 integer :: me_i
219 logical :: is_dp_i
220 integer :: neighborListSize
221 integer :: listerror, error
222 integer :: localError
223
224 real(kind=dp) :: listSkin = 1.0
225
226
227 !! initialize local variables
228
229 #ifdef IS_MPI
230 pot_local = 0.0_dp
231 nlocal = getNlocal()
232 nrow = getNrow(plan_row)
233 ncol = getNcol(plan_col)
234 #else
235 nlocal = getNlocal()
236 natoms = nlocal
237 #endif
238
239 call check_initialization(localError)
240 if ( localError .ne. 0 ) then
241 call handleError("do_force_loop","Not Initialized")
242 error = -1
243 return
244 end if
245 call zero_work_arrays()
246
247 do_pot = do_pot_c
248 do_stress = do_stress_c
249
250
251 ! Gather all information needed by all force loops:
252
253 #ifdef IS_MPI
254
255 call gather(q,q_Row,plan_row3d)
256 call gather(q,q_Col,plan_col3d)
257
258 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
259 call gather(u_l,u_l_Row,plan_row3d)
260 call gather(u_l,u_l_Col,plan_col3d)
261
262 call gather(A,A_Row,plan_row_rotation)
263 call gather(A,A_Col,plan_col_rotation)
264 endif
265
266 #endif
267
268 if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
269 !! See if we need to update neighbor lists
270 call checkNeighborList(nlocal, q, listSkin, update_nlist)
271 !! if_mpi_gather_stuff_for_prepair
272 !! do_prepair_loop_if_needed
273 !! if_mpi_scatter_stuff_from_prepair
274 !! if_mpi_gather_stuff_from_prepair_to_main_loop
275
276 !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
277 #ifdef IS_MPI
278
279 if (update_nlist) then
280
281 !! save current configuration, construct neighbor list,
282 !! and calculate forces
283 call saveNeighborList(nlocal, q)
284
285 neighborListSize = size(list)
286 nlist = 0
287
288 do i = 1, nrow
289 point(i) = nlist + 1
290
291 prepair_inner: do j = 1, ncol
292
293 if (skipThisPair(i,j)) cycle prepair_inner
294
295 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
296
297 if (rijsq < rlistsq) then
298
299 nlist = nlist + 1
300
301 if (nlist > neighborListSize) then
302 call expandNeighborList(nlocal, listerror)
303 if (listerror /= 0) then
304 error = -1
305 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
306 return
307 end if
308 neighborListSize = size(list)
309 endif
310
311 list(nlist) = j
312 call do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot_local)
313 endif
314 enddo prepair_inner
315 enddo
316
317 point(nrow + 1) = nlist + 1
318
319 else !! (of update_check)
320
321 ! use the list to find the neighbors
322 do i = 1, nrow
323 JBEG = POINT(i)
324 JEND = POINT(i+1) - 1
325 ! check thiat molecule i has neighbors
326 if (jbeg .le. jend) then
327
328 do jnab = jbeg, jend
329 j = list(jnab)
330
331 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
332 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
333 u_l, A, f, t, pot_local)
334
335 enddo
336 endif
337 enddo
338 endif
339
340 #else
341
342 if (update_nlist) then
343
344 ! save current configuration, contruct neighbor list,
345 ! and calculate forces
346 call saveNeighborList(natoms, q)
347
348 neighborListSize = size(list)
349
350 nlist = 0
351
352 do i = 1, natoms-1
353 point(i) = nlist + 1
354
355 prepair_inner: do j = i+1, natoms
356
357 if (skipThisPair(i,j)) cycle prepair_inner
358
359 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
360
361
362 if (rijsq < rlistsq) then
363
364
365 nlist = nlist + 1
366
367 if (nlist > neighborListSize) then
368 call expandNeighborList(natoms, listerror)
369 if (listerror /= 0) then
370 error = -1
371 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
372 return
373 end if
374 neighborListSize = size(list)
375 endif
376
377 list(nlist) = j
378
379 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
380 u_l, A, f, t, pot)
381
382 endif
383 enddo prepair_inner
384 enddo
385
386 point(natoms) = nlist + 1
387
388 else !! (update)
389
390 ! use the list to find the neighbors
391 do i = 1, natoms-1
392 JBEG = POINT(i)
393 JEND = POINT(i+1) - 1
394 ! check thiat molecule i has neighbors
395 if (jbeg .le. jend) then
396
397 do jnab = jbeg, jend
398 j = list(jnab)
399
400 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
401 call do_prepair(i, j, rijsq, d, do_pot, do_stress, &
402 u_l, A, f, t, pot)
403
404 enddo
405 endif
406 enddo
407 endif
408 #endif
409 !! Do rest of preforce calculations
410 !! do necessary preforce calculations
411 call do_preforce(nlocal,pot)
412 ! we have already updated the neighbor list set it to false...
413 update_nlist = .false.
414 else
415 !! See if we need to update neighbor lists for non pre-pair
416 call checkNeighborList(nlocal, q, listSkin, update_nlist)
417 endif
418
419
420
421
422
423 !---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>>>>>>>>>>>>>>>>
424
425
426
427
428
429 #ifdef IS_MPI
430
431 if (update_nlist) then
432
433 !! save current configuration, construct neighbor list,
434 !! and calculate forces
435 call saveNeighborList(nlocal, q)
436
437 neighborListSize = size(list)
438 nlist = 0
439
440 do i = 1, nrow
441 point(i) = nlist + 1
442
443 inner: do j = 1, ncol
444
445 if (skipThisPair(i,j)) cycle inner
446
447 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
448
449 if (rijsq < rlistsq) then
450
451 nlist = nlist + 1
452
453 if (nlist > neighborListSize) then
454 call expandNeighborList(nlocal, listerror)
455 if (listerror /= 0) then
456 error = -1
457 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
458 return
459 end if
460 neighborListSize = size(list)
461 endif
462
463 list(nlist) = j
464
465 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
466 u_l, A, f, t, pot_local)
467
468 endif
469 enddo inner
470 enddo
471
472 point(nrow + 1) = nlist + 1
473
474 else !! (of update_check)
475
476 ! use the list to find the neighbors
477 do i = 1, nrow
478 JBEG = POINT(i)
479 JEND = POINT(i+1) - 1
480 ! check thiat molecule i has neighbors
481 if (jbeg .le. jend) then
482
483 do jnab = jbeg, jend
484 j = list(jnab)
485
486 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
487 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
488 u_l, A, f, t, pot_local)
489
490 enddo
491 endif
492 enddo
493 endif
494
495 #else
496
497 if (update_nlist) then
498
499 ! save current configuration, contruct neighbor list,
500 ! and calculate forces
501 call saveNeighborList(natoms, q)
502
503 neighborListSize = size(list)
504
505 nlist = 0
506
507 do i = 1, natoms-1
508 point(i) = nlist + 1
509
510 inner: do j = i+1, natoms
511
512 if (skipThisPair(i,j)) cycle inner
513
514 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
515
516
517 if (rijsq < rlistsq) then
518
519 nlist = nlist + 1
520
521 if (nlist > neighborListSize) then
522 call expandNeighborList(natoms, listerror)
523 if (listerror /= 0) then
524 error = -1
525 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
526 return
527 end if
528 neighborListSize = size(list)
529 endif
530
531 list(nlist) = j
532
533 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
534 u_l, A, f, t, pot)
535
536 endif
537 enddo inner
538 enddo
539
540 point(natoms) = nlist + 1
541
542 else !! (update)
543
544 ! use the list to find the neighbors
545 do i = 1, natoms-1
546 JBEG = POINT(i)
547 JEND = POINT(i+1) - 1
548 ! check thiat molecule i has neighbors
549 if (jbeg .le. jend) then
550
551 do jnab = jbeg, jend
552 j = list(jnab)
553
554 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
555 call do_pair(i, j, rijsq, d, do_pot, do_stress, &
556 u_l, A, f, t, pot)
557
558 enddo
559 endif
560 enddo
561 endif
562
563 #endif
564
565 ! phew, done with main loop.
566
567 #ifdef IS_MPI
568 !!distribute forces
569
570 f_temp = 0.0_dp
571 call scatter(f_Row,f_temp,plan_row3d)
572 do i = 1,nlocal
573 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
574 end do
575
576 f_temp = 0.0_dp
577 call scatter(f_Col,f_temp,plan_col3d)
578 do i = 1,nlocal
579 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
580 end do
581
582 if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then
583 t_temp = 0.0_dp
584 call scatter(t_Row,t_temp,plan_row3d)
585 do i = 1,nlocal
586 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
587 end do
588 t_temp = 0.0_dp
589 call scatter(t_Col,t_temp,plan_col3d)
590
591 do i = 1,nlocal
592 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
593 end do
594 endif
595
596 if (do_pot) then
597 ! scatter/gather pot_row into the members of my column
598 call scatter(pot_Row, pot_Temp, plan_row)
599
600 ! scatter/gather pot_local into all other procs
601 ! add resultant to get total pot
602 do i = 1, nlocal
603 pot_local = pot_local + pot_Temp(i)
604 enddo
605
606 pot_Temp = 0.0_DP
607
608 call scatter(pot_Col, pot_Temp, plan_col)
609 do i = 1, nlocal
610 pot_local = pot_local + pot_Temp(i)
611 enddo
612
613 endif
614 #endif
615
616 if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then
617
618 if (FF_uses_RF .and. SimUsesRF()) then
619
620 #ifdef IS_MPI
621 call scatter(rf_Row,rf,plan_row3d)
622 call scatter(rf_Col,rf_Temp,plan_col3d)
623 do i = 1,nlocal
624 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
625 end do
626 #endif
627
628 do i = 1, getNlocal()
629
630 rfpot = 0.0_DP
631 #ifdef IS_MPI
632 me_i = atid_row(i)
633 #else
634 me_i = atid(i)
635 #endif
636 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
637 if ( is_DP_i ) then
638 call getElementProperty(atypes, me_i, "dipole_moment", mu_i)
639 !! The reaction field needs to include a self contribution
640 !! to the field:
641 call accumulate_self_rf(i, mu_i, u_l)
642 !! Get the reaction field contribution to the
643 !! potential and torques:
644 call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
645 #ifdef IS_MPI
646 pot_local = pot_local + rfpot
647 #else
648 pot = pot + rfpot
649
650 #endif
651 endif
652 enddo
653 endif
654 endif
655
656
657 #ifdef IS_MPI
658
659 if (do_pot) then
660 pot = pot + pot_local
661 !! we assume the c code will do the allreduce to get the total potential
662 !! we could do it right here if we needed to...
663 endif
664
665 if (do_stress) then
666 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
667 mpi_comm_world,mpi_err)
668 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
669 mpi_comm_world,mpi_err)
670 endif
671
672 #else
673
674 if (do_stress) then
675 tau = tau_Temp
676 virial = virial_Temp
677 endif
678
679 #endif
680
681 end subroutine do_force_loop
682
683 subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
684
685 real( kind = dp ) :: pot
686 real( kind = dp ), dimension(3,getNlocal()) :: u_l
687 real (kind=dp), dimension(9,getNlocal()) :: A
688 real (kind=dp), dimension(3,getNlocal()) :: f
689 real (kind=dp), dimension(3,getNlocal()) :: t
690
691 logical, intent(inout) :: do_pot, do_stress
692 integer, intent(in) :: i, j
693 real ( kind = dp ), intent(inout) :: rijsq
694 real ( kind = dp ) :: r
695 real ( kind = dp ), intent(inout) :: d(3)
696 logical :: is_LJ_i, is_LJ_j
697 logical :: is_DP_i, is_DP_j
698 logical :: is_GB_i, is_GB_j
699 logical :: is_EAM_i,is_EAM_j
700 logical :: is_Sticky_i, is_Sticky_j
701 integer :: me_i, me_j
702
703 r = sqrt(rijsq)
704
705 #ifdef IS_MPI
706 if (tagRow(i) .eq. tagColumn(j)) then
707 write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
708 endif
709
710 me_i = atid_row(i)
711 me_j = atid_col(j)
712
713 #else
714
715 me_i = atid(i)
716 me_j = atid(j)
717
718 #endif
719
720 if (FF_uses_LJ .and. SimUsesLJ()) then
721 call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
722 call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
723
724 if ( is_LJ_i .and. is_LJ_j ) &
725 call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
726 endif
727
728 if (FF_uses_dipoles .and. SimUsesDipoles()) then
729 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
730 call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
731
732 if ( is_DP_i .and. is_DP_j ) then
733
734 call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
735 do_pot, do_stress)
736 if (FF_uses_RF .and. SimUsesRF()) then
737 call accumulate_rf(i, j, r, u_l)
738 call rf_correct_forces(i, j, d, r, u_l, f, do_stress)
739 endif
740
741 endif
742 endif
743
744 if (FF_uses_Sticky .and. SimUsesSticky()) then
745
746 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
747 call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
748
749 if ( is_Sticky_i .and. is_Sticky_j ) then
750 call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, &
751 do_pot, do_stress)
752 endif
753 endif
754
755
756 if (FF_uses_GB .and. SimUsesGB()) then
757
758 call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
759 call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
760
761 if ( is_GB_i .and. is_GB_j ) then
762 call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, &
763 do_pot, do_stress)
764 endif
765 endif
766
767
768
769 if (FF_uses_EAM .and. SimUsesEAM()) then
770 call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
771 call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
772
773 if ( is_EAM_i .and. is_EAM_j ) &
774 call do_eam_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
775 endif
776
777
778
779
780 end subroutine do_pair
781
782
783
784 subroutine do_prepair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
785 real( kind = dp ) :: pot
786 real( kind = dp ), dimension(3,getNlocal()) :: u_l
787 real (kind=dp), dimension(9,getNlocal()) :: A
788 real (kind=dp), dimension(3,getNlocal()) :: f
789 real (kind=dp), dimension(3,getNlocal()) :: t
790
791 logical, intent(inout) :: do_pot, do_stress
792 integer, intent(in) :: i, j
793 real ( kind = dp ), intent(inout) :: rijsq
794 real ( kind = dp ) :: r
795 real ( kind = dp ), intent(inout) :: d(3)
796
797 logical :: is_EAM_i, is_EAM_j
798
799 integer :: me_i, me_j
800
801 r = sqrt(rijsq)
802
803
804 #ifdef IS_MPI
805 if (tagRow(i) .eq. tagColumn(j)) then
806 write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
807 endif
808
809 me_i = atid_row(i)
810 me_j = atid_col(j)
811
812 #else
813
814 me_i = atid(i)
815 me_j = atid(j)
816
817 #endif
818
819 if (FF_uses_EAM .and. SimUsesEAM()) then
820 call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
821 call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
822
823 if ( is_EAM_i .and. is_EAM_j ) &
824 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
825 endif
826
827 end subroutine do_prepair
828
829
830
831
832 subroutine do_preforce(nlocal,pot)
833 integer :: nlocal
834 real( kind = dp ) :: pot
835
836 if (FF_uses_EAM .and. SimUsesEAM()) then
837 call calc_EAM_preforce_Frho(nlocal,pot)
838 endif
839
840
841 end subroutine do_preforce
842
843
844 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
845
846 real (kind = dp), dimension(3) :: q_i
847 real (kind = dp), dimension(3) :: q_j
848 real ( kind = dp ), intent(out) :: r_sq
849 real( kind = dp ) :: d(3), scaled(3)
850 integer i
851
852 d(1:3) = q_j(1:3) - q_i(1:3)
853
854 ! Wrap back into periodic box if necessary
855 if ( SimUsesPBC() ) then
856
857 if( .not.boxIsOrthorhombic ) then
858 ! calc the scaled coordinates.
859
860 scaled = matmul(HmatInv, d)
861
862 ! wrap the scaled coordinates
863
864 scaled = scaled - anint(scaled)
865
866
867 ! calc the wrapped real coordinates from the wrapped scaled
868 ! coordinates
869
870 d = matmul(Hmat,scaled)
871
872 else
873 ! calc the scaled coordinates.
874
875 do i = 1, 3
876 scaled(i) = d(i) * HmatInv(i,i)
877
878 ! wrap the scaled coordinates
879
880 scaled(i) = scaled(i) - anint(scaled(i))
881
882 ! calc the wrapped real coordinates from the wrapped scaled
883 ! coordinates
884
885 d(i) = scaled(i)*Hmat(i,i)
886 enddo
887 endif
888
889 endif
890
891 r_sq = dot_product(d,d)
892
893 end subroutine get_interatomic_vector
894
895 subroutine check_initialization(error)
896 integer, intent(out) :: error
897
898 error = 0
899 ! Make sure we are properly initialized.
900 if (.not. do_forces_initialized) then
901 write(*,*) "Forces not initialized"
902 error = -1
903 return
904 endif
905
906 #ifdef IS_MPI
907 if (.not. isMPISimSet()) then
908 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
909 error = -1
910 return
911 endif
912 #endif
913
914 return
915 end subroutine check_initialization
916
917
918 subroutine zero_work_arrays()
919
920 #ifdef IS_MPI
921
922 q_Row = 0.0_dp
923 q_Col = 0.0_dp
924
925 u_l_Row = 0.0_dp
926 u_l_Col = 0.0_dp
927
928 A_Row = 0.0_dp
929 A_Col = 0.0_dp
930
931 f_Row = 0.0_dp
932 f_Col = 0.0_dp
933 f_Temp = 0.0_dp
934
935 t_Row = 0.0_dp
936 t_Col = 0.0_dp
937 t_Temp = 0.0_dp
938
939 pot_Row = 0.0_dp
940 pot_Col = 0.0_dp
941 pot_Temp = 0.0_dp
942
943 rf_Row = 0.0_dp
944 rf_Col = 0.0_dp
945 rf_Temp = 0.0_dp
946
947 #endif
948
949
950 if (FF_uses_EAM .and. SimUsesEAM()) then
951 call clean_EAM()
952 endif
953
954
955
956
957
958 rf = 0.0_dp
959 tau_Temp = 0.0_dp
960 virial_Temp = 0.0_dp
961 end subroutine zero_work_arrays
962
963 function skipThisPair(atom1, atom2) result(skip_it)
964 integer, intent(in) :: atom1
965 integer, intent(in), optional :: atom2
966 logical :: skip_it
967 integer :: unique_id_1, unique_id_2
968 integer :: me_i,me_j
969 integer :: i
970
971 skip_it = .false.
972
973 !! there are a number of reasons to skip a pair or a particle
974 !! mostly we do this to exclude atoms who are involved in short
975 !! range interactions (bonds, bends, torsions), but we also need
976 !! to exclude some overcounted interactions that result from
977 !! the parallel decomposition
978
979 #ifdef IS_MPI
980 !! in MPI, we have to look up the unique IDs for each atom
981 unique_id_1 = tagRow(atom1)
982 #else
983 !! in the normal loop, the atom numbers are unique
984 unique_id_1 = atom1
985 #endif
986
987 !! We were called with only one atom, so just check the global exclude
988 !! list for this atom
989 if (.not. present(atom2)) then
990 do i = 1, nExcludes_global
991 if (excludesGlobal(i) == unique_id_1) then
992 skip_it = .true.
993 return
994 end if
995 end do
996 return
997 end if
998
999 #ifdef IS_MPI
1000 unique_id_2 = tagColumn(atom2)
1001 #else
1002 unique_id_2 = atom2
1003 #endif
1004
1005 #ifdef IS_MPI
1006 !! this situation should only arise in MPI simulations
1007 if (unique_id_1 == unique_id_2) then
1008 skip_it = .true.
1009 return
1010 end if
1011
1012 !! this prevents us from doing the pair on multiple processors
1013 if (unique_id_1 < unique_id_2) then
1014 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1015 skip_it = .true.
1016 return
1017 endif
1018 else
1019 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1020 skip_it = .true.
1021 return
1022 endif
1023 endif
1024 #endif
1025
1026 !! the rest of these situations can happen in all simulations:
1027 do i = 1, nExcludes_global
1028 if ((excludesGlobal(i) == unique_id_1) .or. &
1029 (excludesGlobal(i) == unique_id_2)) then
1030 skip_it = .true.
1031 return
1032 endif
1033 enddo
1034
1035 do i = 1, nExcludes_local
1036 if (excludesLocal(1,i) == unique_id_1) then
1037 if (excludesLocal(2,i) == unique_id_2) then
1038 skip_it = .true.
1039 return
1040 endif
1041 else
1042 if (excludesLocal(1,i) == unique_id_2) then
1043 if (excludesLocal(2,i) == unique_id_1) then
1044 skip_it = .true.
1045 return
1046 endif
1047 endif
1048 endif
1049 end do
1050
1051 return
1052 end function skipThisPair
1053
1054 function FF_UsesDirectionalAtoms() result(doesit)
1055 logical :: doesit
1056 doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1057 FF_uses_GB .or. FF_uses_RF
1058 end function FF_UsesDirectionalAtoms
1059
1060 function FF_RequiresPrepairCalc() result(doesit)
1061 logical :: doesit
1062 doesit = FF_uses_EAM
1063 end function FF_RequiresPrepairCalc
1064
1065 function FF_RequiresPostpairCalc() result(doesit)
1066 logical :: doesit
1067 doesit = FF_uses_RF
1068 end function FF_RequiresPostpairCalc
1069
1070 !! This cleans componets of force arrays belonging only to fortran
1071
1072 end module do_Forces