ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 328
Committed: Wed Mar 12 20:00:58 2003 UTC (21 years, 3 months ago) by gezelter
File size: 13128 byte(s)
Log Message:
bye bye atypeGlobals (part2)

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.14 2003-03-12 20:00:58 gezelter Exp $, $Date: 2003-03-12 20:00:58 $, $Name: not supported by cvs2svn $, $Revision: 1.14 $
8
9
10
11 module do_Forces
12 use simulation
13 use definitions
14 use forceGlobals
15 use atype_module
16 use neighborLists
17 use lj_FF
18 use sticky_FF
19 use dipole_dipole
20 use gb_FF
21
22 #ifdef IS_MPI
23 use mpiSimulation
24 #endif
25 implicit none
26 PRIVATE
27
28 public :: do_force_loop
29
30 logical :: do_pot
31 logical :: do_stress
32
33 contains
34
35 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
36 !------------------------------------------------------------->
37 subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, &
38 FFerror)
39 !! Position array provided by C, dimensioned by getNlocal
40 real ( kind = dp ), dimension(3,getNlocal()) :: q
41 !! Rotation Matrix for each long range particle in simulation.
42 real( kind = dp), dimension(9,getNlocal()) :: A
43 !! Unit vectors for dipoles (lab frame)
44 real( kind = dp ), dimension(3,getNlocal()) :: u_l
45 !! Force array provided by C, dimensioned by getNlocal
46 real ( kind = dp ), dimension(3,getNlocal()) :: f
47 !! Torsion array provided by C, dimensioned by getNlocal
48 real( kind = dp ), dimension(3,getNlocal()) :: t
49 !! Stress Tensor
50 real( kind = dp), dimension(9) :: tau
51 real ( kind = dp ) :: pot
52 logical ( kind = 2) :: do_pot_c, do_stress_c
53 integer :: FFerror
54
55 #ifdef IS_MPI
56 real( kind = DP ) :: pot_local
57 #endif
58
59 logical :: update_nlist
60 integer :: i, j, jbeg, jend, jnab
61 integer :: nlist
62 real( kind = DP ) :: rijsq, rlistsq, rcutsq, rlist, rcut
63
64 #ifdef IS_MPI
65 integer :: nlocal
66 #endif
67 integer :: nrow
68 integer :: ncol
69 integer :: natoms
70 integer :: neighborListSize
71 integer :: listerror
72 FFerror = 0
73
74 do_pot = do_pot_c
75 do_stress = do_stress_c
76
77 ! Make sure we are properly initialized.
78 if (.not. isFFInit) then
79 write(default_error,*) "ERROR: lj_FF has not been properly initialized"
80 FFerror = -1
81 return
82 endif
83 #ifdef IS_MPI
84 if (.not. isMPISimSet()) then
85 write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
86 FFerror = -1
87 return
88 endif
89 #endif
90
91 !! initialize local variables
92 natoms = getNlocal()
93 call getRcut(rcut,rcut2=rcutsq)
94 call getRlist(rlist,rlistsq)
95
96 !! See if we need to update neighbor lists
97 call checkNeighborList(natoms, q, rcut, rlist, update_nlist)
98
99 !--------------WARNING...........................
100 ! Zero variables, NOTE:::: Forces are zeroed in C
101 ! Zeroing them here could delete previously computed
102 ! Forces.
103 !------------------------------------------------
104 call zero_module_variables()
105
106 ! communicate MPI positions
107 #ifdef IS_MPI
108 call gather(q,q_Row,plan_row3d)
109 call gather(q,q_Col,plan_col3d)
110
111 call gather(u_l,u_l_Row,plan_row3d)
112 call gather(u_l,u_l_Col,plan_col3d)
113
114 call gather(A,A_Row,plan_row_rotation)
115 call gather(A,A_Col,plan_col_rotation)
116 #endif
117
118
119 #ifdef IS_MPI
120
121 if (update_nlist) then
122
123 ! save current configuration, contruct neighbor list,
124 ! and calculate forces
125 call save_neighborList(q)
126
127 neighborListSize = getNeighborListSize()
128 nlist = 0
129
130 nrow = getNrow(plan_row)
131 ncol = getNcol(plan_col)
132 nlocal = getNlocal()
133
134 do i = 1, nrow
135 point(i) = nlist + 1
136
137 inner: do j = 1, ncol
138
139 if (check_exclude(i,j)) cycle inner:
140
141 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
142
143 if (rijsq < rlistsq) then
144
145 nlist = nlist + 1
146
147 if (nlist > neighborListSize) then
148 call expandNeighborList(nlocal, listerror)
149 if (listerror /= 0) then
150 FFerror = -1
151 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
152 return
153 end if
154 endif
155
156 list(nlist) = j
157
158 if (rijsq < rcutsq) then
159 call do_pair(i, j, rijsq, d)
160 endif
161 endif
162 enddo inner
163 enddo
164
165 point(nrow + 1) = nlist + 1
166
167 else !! (update)
168
169 ! use the list to find the neighbors
170 do i = 1, nrow
171 JBEG = POINT(i)
172 JEND = POINT(i+1) - 1
173 ! check thiat molecule i has neighbors
174 if (jbeg .le. jend) then
175
176 do jnab = jbeg, jend
177 j = list(jnab)
178
179 call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
180 call do_pair(i, j, rijsq, d)
181
182 enddo
183 endif
184 enddo
185 endif
186
187 #else
188
189 if (update_nlist) then
190
191 ! save current configuration, contruct neighbor list,
192 ! and calculate forces
193 call save_neighborList(q)
194
195 neighborListSize = getNeighborListSize()
196 nlist = 0
197
198 do i = 1, natoms-1
199 point(i) = nlist + 1
200
201 inner: do j = i+1, natoms
202
203 if (check_exclude(i,j)) cycle inner:
204
205 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
206
207 if (rijsq < rlistsq) then
208
209 nlist = nlist + 1
210
211 if (nlist > neighborListSize) then
212 call expandList(natoms, listerror)
213 if (listerror /= 0) then
214 FFerror = -1
215 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
216 return
217 end if
218 endif
219
220 list(nlist) = j
221
222 if (rijsq < rcutsq) then
223 call do_pair(i, j, rijsq, d)
224 endif
225 endif
226 enddo inner
227 enddo
228
229 point(natoms) = nlist + 1
230
231 else !! (update)
232
233 ! use the list to find the neighbors
234 do i = 1, nrow
235 JBEG = POINT(i)
236 JEND = POINT(i+1) - 1
237 ! check thiat molecule i has neighbors
238 if (jbeg .le. jend) then
239
240 do jnab = jbeg, jend
241 j = list(jnab)
242
243 call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
244 call do_pair(i, j, rijsq, d)
245
246 enddo
247 endif
248 enddo
249 endif
250
251 #endif
252
253
254 #ifdef IS_MPI
255 !! distribute all reaction field stuff (or anything for post-pair):
256 call scatter(rflRow,rflTemp1,plan_row3d)
257 call scatter(rflCol,rflTemp2,plan_col3d)
258 do i = 1,nlocal
259 rflTemp(1:3,i) = rflTemp1(1:3,i) + rflTemp2(1:3,i)
260 end do
261 #endif
262
263 ! This is the post-pair loop:
264 #ifdef IS_MPI
265
266 if (system_has_postpair_atoms) then
267 do i = 1, nlocal
268 Atype_i => identPtrListRow(i)%this
269 call do_postpair(i, Atype_i)
270 enddo
271 endif
272
273 #else
274
275 if (system_has_postpair_atoms) then
276 do i = 1, natoms
277 Atype_i => identPtr(i)%this
278 call do_postpair(i, Atype_i)
279 enddo
280 endif
281
282 #endif
283
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 (doTorque()) 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 pot = pot_local
321 endif
322
323 if (doStress()) then
324 mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
325 mpi_comm_world,mpi_err)
326 mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
327 mpi_comm_world,mpi_err)
328 endif
329
330 #endif
331
332 if (doStress()) then
333 tau = tau_Temp
334 virial = virial_Temp
335 endif
336
337 end subroutine do_force_loop
338
339
340 !! Calculate any pre-force loop components and update nlist if necessary.
341 subroutine do_preForce(updateNlist)
342 logical, intent(inout) :: updateNlist
343
344
345
346 end subroutine do_preForce
347
348 !! Calculate any post force loop components, i.e. reaction field, etc.
349 subroutine do_postForce()
350
351
352
353 end subroutine do_postForce
354
355 subroutine do_pair(i, j, rijsq, d)
356
357 integer, intent(in) :: i, j
358 real ( kind = dp ), intent(in) :: rijsq
359 real ( kind = dp ) :: r
360 real ( kind = dp ), intent(inout) :: d(3)
361
362 r = sqrt(rijsq)
363
364 logical :: is_LJ_i, is_LJ_j
365 logical :: is_DP_i, is_DP_j
366 logical :: is_Sticky_i, is_Sticky_j
367 integer :: me_i, me_j
368
369 #ifdef IS_MPI
370
371 me_i = atid_row(i)
372 me_j = atid_col(j)
373
374 #else
375
376 me_i = atid(i)
377 me_j = atid(j)
378
379 #endif
380
381 call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i)
382 call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j)
383
384 if ( is_LJ_i .and. is_LJ_j ) &
385 call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress)
386
387 call getElementProperty(atypes, me_i, "is_DP", is_DP_i)
388 call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
389
390 if ( is_DP_i .and. is_DP_j ) then
391
392 call do_dipole_pair(i, j, d, r, pot, u_l, f, t, do_pot, do_stress)
393
394 if (do_reaction_field) then
395 call accumulate_rf(i, j, r)
396 endif
397
398 endif
399
400 call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i)
401 call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j)
402
403 if ( is_Sticky_i .and. is_Sticky_j ) then
404 call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, do_pot, do_stress)
405 endif
406
407
408 end subroutine do_pair
409
410
411 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
412
413 real (kind = dp), dimension(3) :: q_i
414 real (kind = dp), dimension(3) :: q_j
415 real ( kind = dp ), intent(out) :: r_sq
416 real( kind = dp ) :: d(3)
417
418 d(1:3) = q_i(1:3) - q_j(1:3)
419
420 ! Wrap back into periodic box if necessary
421 if ( isPBC() ) then
422 d(1:3) = d(1:3) - thisSim%box(1:3) * sign(1.0_dp,thisSim%box(1:3)) * &
423 int(abs(d(1:3)/thisSim%box(1:3) + 0.5_dp)
424 endif
425
426 r_sq = dot_product(d,d)
427
428 end subroutine get_interatomic_vector
429
430 subroutine zero_work_arrays()
431
432 #ifdef IS_MPI
433
434 q_Row = 0.0_dp
435 q_Col = 0.0_dp
436
437 u_l_Row = 0.0_dp
438 u_l_Col = 0.0_dp
439
440 A_Row = 0.0_dp
441 A_Col = 0.0_dp
442
443 f_Row = 0.0_dp
444 f_Col = 0.0_dp
445 f_Temp = 0.0_dp
446
447 t_Row = 0.0_dp
448 t_Col = 0.0_dp
449 t_Temp = 0.0_dp
450
451 pot_Row = 0.0_dp
452 pot_Col = 0.0_dp
453 pot_Temp = 0.0_dp
454
455 #endif
456
457 tau_Temp = 0.0_dp
458 virial_Temp = 0.0_dp
459
460 end subroutine zero_work_arrays
461
462
463 !! Function to properly build neighbor lists in MPI using newtons 3rd law.
464 !! We don't want 2 processors doing the same i j pair twice.
465 !! Also checks to see if i and j are the same particle.
466 function checkExcludes(atom1,atom2) result(do_cycle)
467 !--------------- Arguments--------------------------
468 ! Index i
469 integer,intent(in) :: atom1
470 ! Index j
471 integer,intent(in), optional :: atom2
472 ! Result do_cycle
473 logical :: do_cycle
474 !--------------- Local variables--------------------
475 integer :: tag_i
476 integer :: tag_j
477 integer :: i
478 !--------------- END DECLARATIONS------------------
479 do_cycle = .false.
480
481 #ifdef IS_MPI
482 tag_i = tagRow(atom1)
483 #else
484 tag_i = tag(atom1)
485 #endif
486
487 !! Check global excludes first
488 if (.not. present(atom2)) then
489 do i = 1,nGlobalExcludes
490 if (excludeGlobal(i) == tag_i) then
491 do_cycle = .true.
492 return
493 end if
494 end do
495 return !! return after checking globals
496 end if
497
498 !! we return if j not present here.
499 tag_j = tagColumn(j)
500
501
502
503 if (tag_i == tag_j) then
504 do_cycle = .true.
505 return
506 end if
507
508 if (tag_i < tag_j) then
509 if (mod(tag_i + tag_j,2) == 0) do_cycle = .true.
510 return
511 else
512 if (mod(tag_i + tag_j,2) == 1) do_cycle = .true.
513 endif
514
515
516
517 do i = 1, nLocalExcludes
518 if (tag_i = excludes(1,i) .and. excludes(2,i) < 0) then
519 do_cycle = .true.
520 return
521 end if
522 end do
523
524
525 end function checkExcludes
526
527
528 end module do_Forces