ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 312
Committed: Tue Mar 11 17:46:18 2003 UTC (21 years, 6 months ago) by gezelter
File size: 15403 byte(s)
Log Message:
Bunch o' stuff, particularly the vector_class.F90 module

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.10 2003-03-11 17:46:18 gezelter Exp $, $Date: 2003-03-11 17:46:18 $, $Name: not supported by cvs2svn $, $Revision: 1.10 $
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
59 type(atype), pointer :: Atype_i
60 type(atype), pointer :: Atype_j
61
62
63
64
65
66
67 #ifdef IS_MPI
68 real( kind = DP ) :: pot_local
69
70 !! Local arrays needed for MPI
71
72 #endif
73
74
75
76 real( kind = DP ) :: pe
77 logical :: update_nlist
78
79
80 integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
81 integer :: nlist
82 integer :: j_start
83
84 real( kind = DP ) :: r_ij, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
85
86 real( kind = DP ) :: rx_ij, ry_ij, rz_ij, rijsq
87 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
88
89
90
91 ! a rig that need to be fixed.
92 #ifdef IS_MPI
93 real( kind = dp ) :: pe_local
94 integer :: nlocal
95 #endif
96 integer :: nrow
97 integer :: ncol
98 integer :: natoms
99 integer :: neighborListSize
100 integer :: listerror
101 !! should we calculate the stress tensor
102 logical :: do_stress = .false.
103
104
105 FFerror = 0
106
107 ! Make sure we are properly initialized.
108 if (.not. isFFInit) then
109 write(default_error,*) "ERROR: lj_FF has not been properly initialized"
110 FFerror = -1
111 return
112 endif
113 #ifdef IS_MPI
114 if (.not. isMPISimSet()) then
115 write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
116 FFerror = -1
117 return
118 endif
119 #endif
120
121 !! initialize local variables
122 natoms = getNlocal()
123 call getRcut(rcut,rcut2=rcutsq)
124 call getRlist(rlist,rlistsq)
125
126 !! Find ensemble
127 if (isEnsemble("NPT")) do_stress = .true.
128 !! set to wrap
129 if (isPBC()) wrap = .true.
130
131
132
133
134 !! See if we need to update neighbor lists
135 call check(q,update_nlist)
136
137 !--------------WARNING...........................
138 ! Zero variables, NOTE:::: Forces are zeroed in C
139 ! Zeroing them here could delete previously computed
140 ! Forces.
141 !------------------------------------------------
142 call zero_module_variables()
143
144
145 ! communicate MPI positions
146 #ifdef IS_MPI
147 call gather(q,q_Row,plan_row3d)
148 call gather(q,q_Col,plan_col3d)
149
150 call gather(u_l,u_l_Row,plan_row3d)
151 call gather(u_l,u_l_Col,plan_col3d)
152
153 call gather(A,A_Row,plan_row_rotation)
154 call gather(A,A_Col,plan_col_rotation)
155 #endif
156
157
158 #ifdef IS_MPI
159
160 if (update_nlist) then
161
162 ! save current configuration, contruct neighbor list,
163 ! and calculate forces
164 call save_neighborList(q)
165
166 neighborListSize = getNeighborListSize()
167 nlist = 0
168
169 nrow = getNrow(plan_row)
170 ncol = getNcol(plan_col)
171 nlocal = getNlocal()
172
173 do i = 1, nrow
174 point(i) = nlist + 1
175 Atype_i => identPtrListRow(i)%this
176
177 inner: do j = 1, ncol
178 Atype_j => identPtrListColumn(j)%this
179
180 call get_interatomic_vector(i,j,q_Row(:,i),q_Col(:,j),&
181 rxij,ryij,rzij,rijsq,r)
182
183 ! skip the loop if the atoms are identical
184 if (mpi_cycle_jLoop(i,j)) cycle inner:
185
186 if (rijsq < rlistsq) then
187
188 nlist = nlist + 1
189
190 if (nlist > neighborListSize) then
191 call expandList(listerror)
192 if (listerror /= 0) then
193 FFerror = -1
194 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
195 return
196 end if
197 endif
198
199 list(nlist) = j
200
201
202 if (rijsq < rcutsq) then
203 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
204 endif
205 endif
206 enddo inner
207 enddo
208
209 point(nrow + 1) = nlist + 1
210
211 else !! (update)
212
213 ! use the list to find the neighbors
214 do i = 1, nrow
215 JBEG = POINT(i)
216 JEND = POINT(i+1) - 1
217 ! check thiat molecule i has neighbors
218 if (jbeg .le. jend) then
219
220 Atype_i => identPtrListRow(i)%this
221 do jnab = jbeg, jend
222 j = list(jnab)
223 Atype_j = identPtrListColumn(j)%this
224 call get_interatomic_vector(i,j,q_Row(:,i),q_Col(:,j),&
225 rxij,ryij,rzij,rijsq,r)
226
227 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
228 enddo
229 endif
230 enddo
231 endif
232
233 #else
234
235 if (update_nlist) then
236
237 ! save current configuration, contruct neighbor list,
238 ! and calculate forces
239 call save_neighborList(q)
240
241 neighborListSize = getNeighborListSize()
242 nlist = 0
243
244
245 do i = 1, natoms-1
246 point(i) = nlist + 1
247 Atype_i => identPtrList(i)%this
248
249 inner: do j = i+1, natoms
250 Atype_j => identPtrList(j)%this
251 call get_interatomic_vector(i,j,q(:,i),q(:,j),&
252 rxij,ryij,rzij,rijsq,r)
253
254 if (rijsq < rlistsq) then
255
256 nlist = nlist + 1
257
258 if (nlist > neighborListSize) then
259 call expandList(listerror)
260 if (listerror /= 0) then
261 FFerror = -1
262 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
263 return
264 end if
265 endif
266
267 list(nlist) = j
268
269
270 if (rijsq < rcutsq) then
271 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
272 endif
273 endif
274 enddo inner
275 enddo
276
277 point(natoms) = nlist + 1
278
279 else !! (update)
280
281 ! use the list to find the neighbors
282 do i = 1, nrow
283 JBEG = POINT(i)
284 JEND = POINT(i+1) - 1
285 ! check thiat molecule i has neighbors
286 if (jbeg .le. jend) then
287
288 Atype_i => identPtrList(i)%this
289 do jnab = jbeg, jend
290 j = list(jnab)
291 Atype_j = identPtrList(j)%this
292 call get_interatomic_vector(i,j,q(:,i),q(:,j),&
293 rxij,ryij,rzij,rijsq,r)
294 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
295 enddo
296 endif
297 enddo
298 endif
299
300 #endif
301
302
303 #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 #ifdef IS_MPI
335 !!distribute forces
336
337 call scatter(f_Row,f,plan_row3d)
338 call scatter(f_Col,f_temp,plan_col3d)
339 do i = 1,nlocal
340 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
341 end do
342
343 if (doTorque()) then
344 call scatter(t_Row,t,plan_row3d)
345 call scatter(t_Col,t_temp,plan_col3d)
346
347 do i = 1,nlocal
348 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
349 end do
350 endif
351
352 if (do_pot) then
353 ! scatter/gather pot_row into the members of my column
354 call scatter(pot_Row, pot_Temp, plan_row)
355
356 ! scatter/gather pot_local into all other procs
357 ! add resultant to get total pot
358 do i = 1, nlocal
359 pot_local = pot_local + pot_Temp(i)
360 enddo
361
362 pot_Temp = 0.0_DP
363
364 call scatter(pot_Col, pot_Temp, plan_col)
365 do i = 1, nlocal
366 pot_local = pot_local + pot_Temp(i)
367 enddo
368
369 pot = pot_local
370 endif
371
372 if (doStress()) then
373 mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
374 mpi_comm_world,mpi_err)
375 mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
376 mpi_comm_world,mpi_err)
377 endif
378
379 #endif
380
381 if (doStress()) then
382 tau = tau_Temp
383 virial = virial_Temp
384 endif
385
386 end subroutine do_force_loop
387
388
389 !! Calculate any pre-force loop components and update nlist if necessary.
390 subroutine do_preForce(updateNlist)
391 logical, intent(inout) :: updateNlist
392
393
394
395 end subroutine do_preForce
396
397
398
399
400
401
402
403
404
405
406
407
408
409 !! Calculate any post force loop components, i.e. reaction field, etc.
410 subroutine do_postForce()
411
412
413
414 end subroutine do_postForce
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432 subroutine do_pair(atype_i,atype_j,i,j,r_ij,rx_ij,ry_ij,rz_ij)
433 type (atype ), pointer, intent(inout) :: atype_i
434 type (atype ), pointer, intent(inout) :: atype_j
435 integer :: i
436 integer :: j
437 real ( kind = dp ), intent(inout) :: rx_ij
438 real ( kind = dp ), intent(inout) :: ry_ij
439 real ( kind = dp ), intent(inout) :: rz_ij
440
441
442 real( kind = dp ) :: fx = 0.0_dp
443 real( kind = dp ) :: fy = 0.0_dp
444 real( kind = dp ) :: fz = 0.0_dp
445
446 real( kind = dp ) :: drdx = 0.0_dp
447 real( kind = dp ) :: drdy = 0.0_dp
448 real( kind = dp ) :: drdz = 0.0_dp
449
450
451 #ifdef IS_MPI
452
453 if (Atype_i%is_LJ .and. Atype_j%is_LJ) then
454 call do_lj_pair(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
455 pot, f)
456 endif
457
458 if (Atype_i%is_dp .and. Atype_j%is_dp) then
459
460 call do_dipole_pair(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
461 pot, u_l, f, t)
462
463 if (do_reaction_field) then
464 call accumulate_rf(i, j, r_ij, rt, rrf)
465 endif
466
467 endif
468
469 if (Atype_i%is_sticky .and. Atype_j%is_sticky) then
470 call getstickyforce(r, pot, dudr, Atype_i, Atype_j)
471 endif
472
473 #else
474
475 if (Atype_i%is_LJ .and. Atype_j%is_LJ) then
476 call do_lj_pair(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
477 pot, f)
478 endif
479
480 if (Atype_i%is_dp .and. Atype_j%is_dp) then
481 call do_dipole_pair(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
482 pot, u_l, f, t)
483
484 if (do_reaction_field) then
485 call accumulate_rf(i, j, r_ij, rt, rrf)
486 endif
487
488 endif
489
490 if (Atype_i%is_sticky .and. Atype_j%is_sticky) then
491 call getstickyforce(r,pot,dudr, Atype_i, Atype_j)
492 endif
493
494 #endif
495
496
497 end subroutine do_pair
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514 subroutine get_interatomic_vector(q_i,q_j,rx_ij,ry_ij,rz_ij,r_sq,r_ij)
515 !---------------- Arguments-------------------------------
516 !! index i
517
518 !! Position array
519 real (kind = dp), dimension(3) :: q_i
520 real (kind = dp), dimension(3) :: q_j
521 !! x component of vector between i and j
522 real ( kind = dp ), intent(out) :: rx_ij
523 !! y component of vector between i and j
524 real ( kind = dp ), intent(out) :: ry_ij
525 !! z component of vector between i and j
526 real ( kind = dp ), intent(out) :: rz_ij
527 !! magnitude of r squared
528 real ( kind = dp ), intent(out) :: r_sq
529 !! magnitude of vector r between atoms i and j.
530 real ( kind = dp ), intent(out) :: r_ij
531 !! wrap into periodic box.
532 logical, intent(in) :: wrap
533
534 !--------------- Local Variables---------------------------
535 !! Distance between i and j
536 real( kind = dp ) :: d(3)
537 !---------------- END DECLARATIONS-------------------------
538
539
540 ! Find distance between i and j
541 d(1:3) = q_i(1:3) - q_j(1:3)
542
543 ! Wrap back into periodic box if necessary
544 if ( wrap ) then
545 d(1:3) = d(1:3) - thisSim%box(1:3) * sign(1.0_dp,thisSim%box(1:3)) * &
546 int(abs(d(1:3)/thisSim%box(1:3) + 0.5_dp)
547 end if
548
549 ! Find Magnitude of the vector
550 r_sq = dot_product(d,d)
551 r_ij = sqrt(r_sq)
552
553 ! Set each component for force calculation
554 rx_ij = d(1)
555 ry_ij = d(2)
556 rz_ij = d(3)
557
558
559 end subroutine get_interatomic_vector
560
561 subroutine zero_module_variables()
562
563 #ifndef IS_MPI
564
565 pe = 0.0E0_DP
566 tauTemp = 0.0_dp
567 fTemp = 0.0_dp
568 tTemp = 0.0_dp
569 #else
570 qRow = 0.0_dp
571 qCol = 0.0_dp
572
573 muRow = 0.0_dp
574 muCol = 0.0_dp
575
576 u_lRow = 0.0_dp
577 u_lCol = 0.0_dp
578
579 ARow = 0.0_dp
580 ACol = 0.0_dp
581
582 fRow = 0.0_dp
583 fCol = 0.0_dp
584
585
586 tRow = 0.0_dp
587 tCol = 0.0_dp
588
589
590
591 eRow = 0.0_dp
592 eCol = 0.0_dp
593 eTemp = 0.0_dp
594 #endif
595
596 end subroutine zero_module_variables
597
598
599 !! Function to properly build neighbor lists in MPI using newtons 3rd law.
600 !! We don't want 2 processors doing the same i j pair twice.
601 !! Also checks to see if i and j are the same particle.
602 function checkExcludes(atom1,atom2) result(do_cycle)
603 !--------------- Arguments--------------------------
604 ! Index i
605 integer,intent(in) :: atom1
606 ! Index j
607 integer,intent(in), optional :: atom2
608 ! Result do_cycle
609 logical :: do_cycle
610 !--------------- Local variables--------------------
611 integer :: tag_i
612 integer :: tag_j
613 integer :: i
614 !--------------- END DECLARATIONS------------------
615 do_cycle = .false.
616
617 #ifdef IS_MPI
618 tag_i = tagRow(atom1)
619 #else
620 tag_i = tag(atom1)
621 #endif
622
623 !! Check global excludes first
624 if (.not. present(atom2)) then
625 do i = 1,nGlobalExcludes
626 if (excludeGlobal(i) == tag_i) then
627 do_cycle = .true.
628 return
629 end if
630 end do
631 return !! return after checking globals
632 end if
633
634 !! we return if j not present here.
635 tag_j = tagColumn(j)
636
637
638
639 if (tag_i == tag_j) then
640 do_cycle = .true.
641 return
642 end if
643
644 if (tag_i < tag_j) then
645 if (mod(tag_i + tag_j,2) == 0) do_cycle = .true.
646 return
647 else
648 if (mod(tag_i + tag_j,2) == 1) do_cycle = .true.
649 endif
650
651
652
653 do i = 1, nLocalExcludes
654 if (tag_i = excludes(1,i) .and. excludes(2,i) < 0) then
655 do_cycle = .true.
656 return
657 end if
658 end do
659
660
661 end function checkExcludes
662
663
664 end module do_Forces