ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 323
Committed: Wed Mar 12 15:39:01 2003 UTC (21 years, 3 months ago) by gezelter
File size: 13367 byte(s)
Log Message:
bunch of fixes

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