ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 897
Committed: Mon Jan 5 22:18:52 2004 UTC (20 years, 6 months ago) by chuckv
File size: 31007 byte(s)
Log Message:
mangling forces even further

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