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