ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/do_Forces.F90
Revision: 297
Committed: Thu Mar 6 22:08:29 2003 UTC (21 years, 6 months ago) by gezelter
File size: 23620 byte(s)
Log Message:
Dan Goes Crazy

File Contents

# Content
1 !! Calculates Long Range forces.
2 !! @author Charles F. Vardeman II
3 !! @author Matthew Meineke
4 !! @version $Id: do_Forces.F90,v 1.4 2003-03-06 22:08:29 gezelter Exp $, $Date: 2003-03-06 22:08:29 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
5
6
7
8 module do_Forces
9 use simulation
10 use definitions
11 use generic_atypes
12 use neighborLists
13
14 use lj_FF
15 use sticky_FF
16 use dp_FF
17 use gb_FF
18
19 #ifdef IS_MPI
20 use mpiSimulation
21 #endif
22 implicit none
23 PRIVATE
24
25 !! Number of lj_atypes in lj_atype_list
26 integer, save :: n_atypes = 0
27
28 !! Global list of lj atypes in simulation
29 type (atype), pointer :: ListHead => null()
30 type (atype), pointer :: ListTail => null()
31
32
33
34
35 logical, save :: firstTime = .True.
36
37 !! Atype identity pointer lists
38 #ifdef IS_MPI
39 !! Row lj_atype pointer list
40 type (identPtrList), dimension(:), pointer :: identPtrListRow => null()
41 !! Column lj_atype pointer list
42 type (identPtrList), dimension(:), pointer :: identPtrListColumn => null()
43 #else
44 type(identPtrList ), dimension(:), pointer :: identPtrList => null()
45 #endif
46
47
48 !! Logical has lj force field module been initialized?
49 logical, save :: isFFinit = .false.
50
51 !! Use periodic boundry conditions
52 logical :: wrap = .false.
53
54 !! Potential energy global module variables
55 #ifdef IS_MPI
56 real(kind = dp), dimension(3,getNrow(plan_row)) :: qRow = 0.0_dp
57 real(kind = dp), dimension(3,getNcol(plan_col)) :: qCol = 0.0_dp
58
59 real(kind = dp), dimension(3,getNrow(plan_row)) :: muRow = 0.0_dp
60 real(kind = dp), dimension(3,getNcol(plan_col)) :: muCol = 0.0_dp
61
62 real(kind = dp), dimension(3,getNrow(plan_row)) :: u_lRow = 0.0_dp
63 real(kind = dp), dimension(3,getNcol(plan_col)) :: u_lCol = 0.0_dp
64
65 real(kind = dp), dimension(9,getNrow(plan_row)) :: ARow = 0.0_dp
66 real(kind = dp), dimension(9,getNcol(plan_col)) :: ACol = 0.0_dp
67
68 real(kind = dp), dimension(3,getNrow(plan_row)) :: fRow = 0.0_dp
69 real(kind = dp), dimension(3,getNcol(plan_col)) :: fCol = 0.0_dp
70 real(kind = dp), dimension(3,getNlocal()) :: fTemp1 = 0.0_dp
71 real(kind = dp), dimension(3,getNlocal()) :: tTemp1 = 0.0_dp
72 real(kind = dp), dimension(3,getNlocal()) :: fTemp2 = 0.0_dp
73 real(kind = dp), dimension(3,getNlocal()) :: tTemp2 = 0.0_dp
74 real(kind = dp), dimension(3,getNlocal()) :: fTemp = 0.0_dp
75 real(kind = dp), dimension(3,getNlocal()) :: tTemp = 0.0_dp
76
77 real(kind = dp), dimension(3,getNrow(plan_row)) :: tRow = 0.0_dp
78 real(kind = dp), dimension(3,getNcol(plan_col)) :: tCol = 0.0_dp
79
80 real(kind = dp), dimension(3,getNrow(plan_row)) :: rflRow = 0.0_dp
81 real(kind = dp), dimension(3,getNcol(plan_col)) :: rflCol = 0.0_dp
82 real(kind = dp), dimension(3,getNlocal()) :: rflTemp = 0.0_dp
83
84 real(kind = dp), dimension(getNrow(plan_row)) :: eRow = 0.0_dp
85 real(kind = dp), dimension(getNcol(plan_col)) :: eCol = 0.0_dp
86
87 real(kind = dp), dimension(getNlocal()) :: eTemp = 0.0_dp
88 #endif
89 real(kind = dp) :: pe = 0.0_dp
90 real(kind = dp), dimension(3,natoms) :: fTemp = 0.0_dp
91 real(kind = dp), dimension(3,natoms) :: tTemp = 0.0_dp
92 real(kind = dp), dimension(3,natoms) :: rflTemp = 0.0_dp
93 real(kind = dp), dimension(9) :: tauTemp = 0.0_dp
94
95 logical :: do_preForce = .false.
96 logical :: do_postForce = .false.
97
98
99
100 !! Public methods and data
101 public :: new_atype
102 public :: do_forceLoop
103 public :: init_FF
104
105
106
107
108 contains
109
110 !! Adds a new lj_atype to the list.
111 subroutine new_atype(ident,mass,epsilon,sigma, &
112 is_LJ,is_Sticky,is_DP,is_GB,w0,v0,dipoleMoment,status)
113 real( kind = dp ), intent(in) :: mass
114 real( kind = dp ), intent(in) :: epsilon
115 real( kind = dp ), intent(in) :: sigma
116 real( kind = dp ), intent(in) :: w0
117 real( kind = dp ), intent(in) :: v0
118 real( kind = dp ), intent(in) :: dipoleMoment
119
120 integer, intent(in) :: ident
121 integer, intent(out) :: status
122 integer, intent(in) :: is_Sticky
123 integer, intent(in) :: is_DP
124 integer, intent(in) :: is_GB
125 integer, intent(in) :: is_LJ
126
127
128 type (atype), pointer :: the_new_atype
129 integer :: alloc_error
130 integer :: atype_counter = 0
131 integer :: alloc_size
132 integer :: err_stat
133 status = 0
134
135
136
137 ! allocate a new atype
138 allocate(the_new_atype,stat=alloc_error)
139 if (alloc_error /= 0 ) then
140 status = -1
141 return
142 end if
143
144 ! assign our new atype information
145 the_new_atype%mass = mass
146 the_new_atype%epsilon = epsilon
147 the_new_atype%sigma = sigma
148 the_new_atype%sigma2 = sigma * sigma
149 the_new_atype%sigma6 = the_new_atype%sigma2 * the_new_atype%sigma2 &
150 * the_new_atype%sigma2
151 the_new_atype%w0 = w0
152 the_new_atype%v0 = v0
153 the_new_atype%dipoleMoment = dipoleMoment
154
155
156 ! assume that this atype will be successfully added
157 the_new_atype%atype_ident = ident
158 the_new_atype%atype_number = n_lj_atypes + 1
159
160 if ( is_Sticky /= 0 ) the_new_atype%is_Sticky = .true.
161 if ( is_GB /= 0 ) the_new_atype%is_GB = .true.
162 if ( is_LJ /= 0 ) the_new_atype%is_LJ = .true.
163 if ( is_DP /= 0 ) the_new_atype%is_DP = .true.
164
165 call add_atype(the_new_atype,ListHead,ListTail,err_stat)
166 if (err_stat /= 0 ) then
167 status = -1
168 return
169 endif
170
171 n_atypes = n_atypes + 1
172
173
174 end subroutine new_atype
175
176
177 subroutine init_FF(nComponents,ident, status)
178 !! Number of components in ident array
179 integer, intent(inout) :: nComponents
180 !! Array of identities nComponents long corresponding to
181 !! ljatype ident.
182 integer, dimension(nComponents),intent(inout) :: ident
183 !! Result status, success = 0, error = -1
184 integer, intent(out) :: Status
185
186 integer :: alloc_stat
187
188 integer :: thisStat
189 integer :: i
190
191 integer :: myNode
192 #ifdef IS_MPI
193 integer, allocatable, dimension(:) :: identRow
194 integer, allocatable, dimension(:) :: identCol
195 integer :: nrow
196 integer :: ncol
197 #endif
198 status = 0
199
200
201
202
203 !! if were're not in MPI, we just update ljatypePtrList
204 #ifndef IS_MPI
205 call create_IdentPtrlst(ident,ListHead,identPtrList,thisStat)
206 if ( thisStat /= 0 ) then
207 status = -1
208 return
209 endif
210
211
212 ! if were're in MPI, we also have to worry about row and col lists
213 #else
214
215 ! We can only set up forces if mpiSimulation has been setup.
216 if (.not. isMPISimSet()) then
217 write(default_error,*) "MPI is not set"
218 status = -1
219 return
220 endif
221 nrow = getNrow(plan_row)
222 ncol = getNcol(plan_col)
223 mynode = getMyNode()
224 !! Allocate temperary arrays to hold gather information
225 allocate(identRow(nrow),stat=alloc_stat)
226 if (alloc_stat /= 0 ) then
227 status = -1
228 return
229 endif
230
231 allocate(identCol(ncol),stat=alloc_stat)
232 if (alloc_stat /= 0 ) then
233 status = -1
234 return
235 endif
236
237 !! Gather idents into row and column idents
238
239 call gather(ident,identRow,plan_row)
240 call gather(ident,identCol,plan_col)
241
242
243 !! Create row and col pointer lists
244
245 call create_IdentPtrlst(identRow,ListHead,identPtrListRow,thisStat)
246 if (thisStat /= 0 ) then
247 status = -1
248 return
249 endif
250
251 call create_IdentPtrlst(identCol,ListHead,identPtrListColumn,thisStat)
252 if (thisStat /= 0 ) then
253 status = -1
254 return
255 endif
256
257 !! free temporary ident arrays
258 if (allocated(identCol)) then
259 deallocate(identCol)
260 end if
261 if (allocated(identCol)) then
262 deallocate(identRow)
263 endif
264
265 #endif
266
267 call initForce_Modules(thisStat)
268 if (thisStat /= 0) then
269 status = -1
270 return
271 endif
272
273 !! Create neighbor lists
274 call expandList(thisStat)
275 if (thisStat /= 0) then
276 status = -1
277 return
278 endif
279
280 isFFinit = .true.
281
282
283 end subroutine init_FF
284
285
286
287
288 subroutine initForce_Modules(thisStat)
289 integer, intent(out) :: thisStat
290 integer :: my_status
291
292 thisStat = 0
293 call init_lj_FF(ListHead,my_status)
294 if (my_status /= 0) then
295 thisStat = -1
296 return
297 end if
298
299 end subroutine initForce_Modules
300
301
302
303
304 !! FORCE routine Calculates Lennard Jones forces.
305 !------------------------------------------------------------->
306 subroutine do_force_loop(q,A,mu,u_l,f,t,tau,potE,do_pot,FFerror)
307 !! Position array provided by C, dimensioned by getNlocal
308 real ( kind = dp ), dimension(3,getNlocal()) :: q
309 !! Rotation Matrix for each long range particle in simulation.
310 real( kind = dp), dimension(9,getNlocal()) :: A
311
312 !! Magnitude dipole moment
313 real( kind = dp ), dimension(3,getNlocal()) :: mu
314 !! Unit vectors for dipoles (lab frame)
315 real( kind = dp ), dimension(3,getNlocal()) :: u_l
316 !! Force array provided by C, dimensioned by getNlocal
317 real ( kind = dp ), dimension(3,getNlocal()) :: f
318 !! Torsion array provided by C, dimensioned by getNlocal
319 real( kind = dp ), dimension(3,getNlocal()) :: t
320
321 !! Stress Tensor
322 real( kind = dp), dimension(9) :: tau
323
324 real ( kind = dp ) :: potE
325 logical ( kind = 2) :: do_pot
326 integer :: FFerror
327
328
329 type(atype), pointer :: Atype_i
330 type(atype), pointer :: Atype_j
331
332
333
334
335
336
337 #ifdef IS_MPI
338 real( kind = DP ) :: pot_local
339
340 !! Local arrays needed for MPI
341
342 #endif
343
344
345
346 real( kind = DP ) :: pe
347 logical :: update_nlist
348
349
350 integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
351 integer :: nlist
352 integer :: j_start
353
354 real( kind = DP ) :: r_ij, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
355
356 real( kind = DP ) :: rx_ij, ry_ij, rz_ij, rijsq
357 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
358
359 real( kind = DP ) :: dielectric = 0.0_dp
360
361 ! a rig that need to be fixed.
362 #ifdef IS_MPI
363 real( kind = dp ) :: pe_local
364 integer :: nlocal
365 #endif
366 integer :: nrow
367 integer :: ncol
368 integer :: natoms
369 integer :: neighborListSize
370 integer :: listerror
371 !! should we calculate the stress tensor
372 logical :: do_stress = .false.
373
374
375 FFerror = 0
376
377 ! Make sure we are properly initialized.
378 if (.not. isFFInit) then
379 write(default_error,*) "ERROR: lj_FF has not been properly initialized"
380 FFerror = -1
381 return
382 endif
383 #ifdef IS_MPI
384 if (.not. isMPISimSet()) then
385 write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
386 FFerror = -1
387 return
388 endif
389 #endif
390
391 !! initialize local variables
392 natoms = getNlocal()
393 call getRcut(rcut,rcut2=rcutsq)
394 call getRlist(rlist,rlistsq)
395
396 !! Find ensemble
397 if (isEnsemble("NPT")) do_stress = .true.
398 !! set to wrap
399 if (isPBC()) wrap = .true.
400
401
402
403
404 !! See if we need to update neighbor lists
405 call check(q,update_nlist)
406
407 !--------------WARNING...........................
408 ! Zero variables, NOTE:::: Forces are zeroed in C
409 ! Zeroing them here could delete previously computed
410 ! Forces.
411 !------------------------------------------------
412 call zero_module_variables()
413
414
415 ! communicate MPI positions
416 #ifdef IS_MPI
417 call gather(q,qRow,plan_row3d)
418 call gather(q,qCol,plan_col3d)
419
420 call gather(mu,muRow,plan_row3d)
421 call gather(mu,muCol,plan_col3d)
422
423 call gather(u_l,u_lRow,plan_row3d)
424 call gather(u_l,u_lCol,plan_col3d)
425
426 call gather(A,ARow,plan_row_rotation)
427 call gather(A,ACol,plan_col_rotation)
428 #endif
429
430
431 #ifdef IS_MPI
432
433 if (update_nlist) then
434
435 ! save current configuration, contruct neighbor list,
436 ! and calculate forces
437 call save_neighborList(q)
438
439 neighborListSize = getNeighborListSize()
440 nlist = 0
441
442 nrow = getNrow(plan_row)
443 ncol = getNcol(plan_col)
444 nlocal = getNlocal()
445
446 do i = 1, nrow
447 point(i) = nlist + 1
448 Atype_i => identPtrListRow(i)%this
449
450 inner: do j = 1, ncol
451 Atype_j => identPtrListColumn(j)%this
452
453 call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),&
454 rxij,ryij,rzij,rijsq,r)
455
456 ! skip the loop if the atoms are identical
457 if (mpi_cycle_jLoop(i,j)) cycle inner:
458
459 if (rijsq < rlistsq) then
460
461 nlist = nlist + 1
462
463 if (nlist > neighborListSize) then
464 call expandList(listerror)
465 if (listerror /= 0) then
466 FFerror = -1
467 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
468 return
469 end if
470 endif
471
472 list(nlist) = j
473
474
475 if (rijsq < rcutsq) then
476 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
477 endif
478 endif
479 enddo inner
480 enddo
481
482 point(nrow + 1) = nlist + 1
483
484 else !! (update)
485
486 ! use the list to find the neighbors
487 do i = 1, nrow
488 JBEG = POINT(i)
489 JEND = POINT(i+1) - 1
490 ! check thiat molecule i has neighbors
491 if (jbeg .le. jend) then
492
493 Atype_i => identPtrListRow(i)%this
494 do jnab = jbeg, jend
495 j = list(jnab)
496 Atype_j = identPtrListColumn(j)%this
497 call get_interatomic_vector(i,j,qRow(:,i),qCol(:,j),&
498 rxij,ryij,rzij,rijsq,r)
499
500 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
501 enddo
502 endif
503 enddo
504 endif
505
506 #else
507
508 if (update_nlist) then
509
510 ! save current configuration, contruct neighbor list,
511 ! and calculate forces
512 call save_neighborList(q)
513
514 neighborListSize = getNeighborListSize()
515 nlist = 0
516
517
518 do i = 1, natoms-1
519 point(i) = nlist + 1
520 Atype_i => identPtrList(i)%this
521
522 inner: do j = i+1, natoms
523 Atype_j => identPtrList(j)%this
524 call get_interatomic_vector(i,j,q(:,i),q(:,j),&
525 rxij,ryij,rzij,rijsq,r)
526
527 if (rijsq < rlistsq) then
528
529 nlist = nlist + 1
530
531 if (nlist > neighborListSize) then
532 call expandList(listerror)
533 if (listerror /= 0) then
534 FFerror = -1
535 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
536 return
537 end if
538 endif
539
540 list(nlist) = j
541
542
543 if (rijsq < rcutsq) then
544 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
545 endif
546 endif
547 enddo inner
548 enddo
549
550 point(natoms) = nlist + 1
551
552 else !! (update)
553
554 ! use the list to find the neighbors
555 do i = 1, nrow
556 JBEG = POINT(i)
557 JEND = POINT(i+1) - 1
558 ! check thiat molecule i has neighbors
559 if (jbeg .le. jend) then
560
561 Atype_i => identPtrList(i)%this
562 do jnab = jbeg, jend
563 j = list(jnab)
564 Atype_j = identPtrList(j)%this
565 call get_interatomic_vector(i,j,q(:,i),q(:,j),&
566 rxij,ryij,rzij,rijsq,r)
567 call do_pair(Atype_i,Atype_j,i,j,r,rxij,ryij,rzij)
568 enddo
569 endif
570 enddo
571 endif
572
573 #endif
574
575
576 #ifdef IS_MPI
577 !! distribute all reaction field stuff (or anything for post-pair):
578 call scatter(rflRow,rflTemp1,plan_row3d)
579 call scatter(rflCol,rflTemp2,plan_col3d)
580 do i = 1,nlocal
581 rflTemp(1:3,i) = rflTemp1(1:3,i) + rflTemp2(1:3,i)
582 end do
583 #endif
584
585 ! This is the post-pair loop:
586 #ifdef IS_MPI
587
588 if (system_has_postpair_atoms) then
589 do i = 1, nlocal
590 Atype_i => identPtrListRow(i)%this
591 call do_postpair(i, Atype_i)
592 enddo
593 endif
594
595 #else
596
597 if (system_has_postpair_atoms) then
598 do i = 1, natoms
599 Atype_i => identPtr(i)%this
600 call do_postpair(i, Atype_i)
601 enddo
602 endif
603
604 #endif
605
606
607
608
609 #ifdef IS_MPI
610 !!distribute forces
611
612 call scatter(fRow,fTemp1,plan_row3d)
613 call scatter(fCol,fTemp2,plan_col3d)
614
615
616 do i = 1,nlocal
617 fTemp(1:3,i) = fTemp1(1:3,i) + fTemp2(1:3,i)
618 end do
619
620 if (do_torque) then
621 call scatter(tRow,tTemp1,plan_row3d)
622 call scatter(tCol,tTemp2,plan_col3d)
623
624 do i = 1,nlocal
625 tTemp(1:3,i) = tTemp1(1:3,i) + tTemp2(1:3,i)
626 end do
627 endif
628
629 if (do_pot) then
630 ! scatter/gather pot_row into the members of my column
631 call scatter(eRow,eTemp,plan_row)
632
633 ! scatter/gather pot_local into all other procs
634 ! add resultant to get total pot
635 do i = 1, nlocal
636 pe_local = pe_local + eTemp(i)
637 enddo
638
639 eTemp = 0.0E0_DP
640 call scatter(eCol,eTemp,plan_col)
641 do i = 1, nlocal
642 pe_local = pe_local + eTemp(i)
643 enddo
644
645 pe = pe_local
646 endif
647 #else
648 ! Copy local array into return array for c
649 f = f+fTemp
650 t = t+tTemp
651 #endif
652
653 potE = pe
654
655
656 if (do_stress) then
657 #ifdef IS_MPI
658 mpi_allreduce = (tau,tauTemp,9,mpi_double_precision,mpi_sum, &
659 mpi_comm_world,mpi_err)
660 #else
661 tau = tauTemp
662 #endif
663 endif
664
665 end subroutine do_force_loop
666
667
668
669
670
671
672
673
674
675
676 !! Calculate any pre-force loop components and update nlist if necessary.
677 subroutine do_preForce(updateNlist)
678 logical, intent(inout) :: updateNlist
679
680
681
682 end subroutine do_preForce
683
684
685
686
687
688
689
690
691
692
693
694
695
696 !! Calculate any post force loop components, i.e. reaction field, etc.
697 subroutine do_postForce()
698
699
700
701 end subroutine do_postForce
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719 subroutine do_pair(atype_i,atype_j,i,j,r_ij,rx_ij,ry_ij,rz_ij)
720 type (atype ), pointer, intent(inout) :: atype_i
721 type (atype ), pointer, intent(inout) :: atype_j
722 integer :: i
723 integer :: j
724 real ( kind = dp ), intent(inout) :: rx_ij
725 real ( kind = dp ), intent(inout) :: ry_ij
726 real ( kind = dp ), intent(inout) :: rz_ij
727
728
729 real( kind = dp ) :: fx = 0.0_dp
730 real( kind = dp ) :: fy = 0.0_dp
731 real( kind = dp ) :: fz = 0.0_dp
732
733 real( kind = dp ) :: drdx = 0.0_dp
734 real( kind = dp ) :: drdy = 0.0_dp
735 real( kind = dp ) :: drdz = 0.0_dp
736
737
738 if (Atype_i%is_LJ .and. Atype_j%is_LJ) then
739 call getLJForce(r,pot,dudr,ljAtype_i,ljAtype_j,fx,fy,fz)
740 endif
741
742 if (Atype_i%is_dp .and. Atype_j%is_dp) then
743
744 #ifdef IS_MPI
745 call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
746 ulRow(:,i), ulCol(:,j), rt, rrf, pot)
747 #else
748 call dipole_dipole(i, j, atype_i, atype_j, rx_ij, ry_ij, rz_ij, r_ij, &
749 ul(:,i), ul(:,j), rt, rrf, pot)
750 #endif
751
752 if (do_reaction_field) then
753 #ifdef IS_MPI
754 call accumulate_rf(i, j, r_ij, rflRow(:,i), rflCol(:j), &
755 ulRow(:i), ulCol(:,j), rt, rrf)
756 #else
757 call accumulate_rf(i, j, r_ij, rfl(:,i), rfl(:j), &
758 ul(:,i), ul(:,j), rt, rrf)
759 #endif
760 endif
761
762
763 endif
764
765 if (Atype_i%is_sticky .and. Atype_j%is_sticky) then
766 call getstickyforce(r,pot,dudr,ljAtype_i,ljAtype_j)
767 endif
768
769
770 #ifdef IS_MPI
771 eRow(i) = eRow(i) + pot*0.5
772 eCol(i) = eCol(i) + pot*0.5
773 #else
774 pe = pe + pot
775 #endif
776
777 drdx = -rxij / r
778 drdy = -ryij / r
779 drdz = -rzij / r
780
781 fx = dudr * drdx
782 fy = dudr * drdy
783 fz = dudr * drdz
784
785
786
787
788
789
790
791 #ifdef IS_MPI
792 fCol(1,j) = fCol(1,j) - fx
793 fCol(2,j) = fCol(2,j) - fy
794 fCol(3,j) = fCol(3,j) - fz
795
796 fRow(1,j) = fRow(1,j) + fx
797 fRow(2,j) = fRow(2,j) + fy
798 fRow(3,j) = fRow(3,j) + fz
799 #else
800 fTemp(1,j) = fTemp(1,j) - fx
801 fTemp(2,j) = fTemp(2,j) - fy
802 fTemp(3,j) = fTemp(3,j) - fz
803 fTemp(1,i) = fTemp(1,i) + fx
804 fTemp(2,i) = fTemp(2,i) + fy
805 fTemp(3,i) = fTemp(3,i) + fz
806 #endif
807
808 if (do_stress) then
809 tauTemp(1) = tauTemp(1) + fx * rxij
810 tauTemp(2) = tauTemp(2) + fx * ryij
811 tauTemp(3) = tauTemp(3) + fx * rzij
812 tauTemp(4) = tauTemp(4) + fy * rxij
813 tauTemp(5) = tauTemp(5) + fy * ryij
814 tauTemp(6) = tauTemp(6) + fy * rzij
815 tauTemp(7) = tauTemp(7) + fz * rxij
816 tauTemp(8) = tauTemp(8) + fz * ryij
817 tauTemp(9) = tauTemp(9) + fz * rzij
818 endif
819
820
821
822 end subroutine do_pair
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839 subroutine get_interatomic_vector(q_i,q_j,rx_ij,ry_ij,rz_ij,r_sq,r_ij)
840 !---------------- Arguments-------------------------------
841 !! index i
842
843 !! Position array
844 real (kind = dp), dimension(3) :: q_i
845 real (kind = dp), dimension(3) :: q_j
846 !! x component of vector between i and j
847 real ( kind = dp ), intent(out) :: rx_ij
848 !! y component of vector between i and j
849 real ( kind = dp ), intent(out) :: ry_ij
850 !! z component of vector between i and j
851 real ( kind = dp ), intent(out) :: rz_ij
852 !! magnitude of r squared
853 real ( kind = dp ), intent(out) :: r_sq
854 !! magnitude of vector r between atoms i and j.
855 real ( kind = dp ), intent(out) :: r_ij
856 !! wrap into periodic box.
857 logical, intent(in) :: wrap
858
859 !--------------- Local Variables---------------------------
860 !! Distance between i and j
861 real( kind = dp ) :: d(3)
862 !---------------- END DECLARATIONS-------------------------
863
864
865 ! Find distance between i and j
866 d(1:3) = q_i(1:3) - q_j(1:3)
867
868 ! Wrap back into periodic box if necessary
869 if ( wrap ) then
870 d(1:3) = d(1:3) - thisSim%box(1:3) * sign(1.0_dp,thisSim%box(1:3)) * &
871 int(abs(d(1:3)/thisSim%box(1:3) + 0.5_dp)
872 end if
873
874 ! Find Magnitude of the vector
875 r_sq = dot_product(d,d)
876 r_ij = sqrt(r_sq)
877
878 ! Set each component for force calculation
879 rx_ij = d(1)
880 ry_ij = d(2)
881 rz_ij = d(3)
882
883
884 end subroutine get_interatomic_vector
885
886 subroutine zero_module_variables()
887
888 #ifndef IS_MPI
889
890 pe = 0.0E0_DP
891 tauTemp = 0.0_dp
892 fTemp = 0.0_dp
893 tTemp = 0.0_dp
894 #else
895 qRow = 0.0_dp
896 qCol = 0.0_dp
897
898 muRow = 0.0_dp
899 muCol = 0.0_dp
900
901 u_lRow = 0.0_dp
902 u_lCol = 0.0_dp
903
904 ARow = 0.0_dp
905 ACol = 0.0_dp
906
907 fRow = 0.0_dp
908 fCol = 0.0_dp
909
910
911 tRow = 0.0_dp
912 tCol = 0.0_dp
913
914
915
916 eRow = 0.0_dp
917 eCol = 0.0_dp
918 eTemp = 0.0_dp
919 #endif
920
921 end subroutine zero_module_variables
922
923 #ifdef IS_MPI
924 !! Function to properly build neighbor lists in MPI using newtons 3rd law.
925 !! We don't want 2 processors doing the same i j pair twice.
926 !! Also checks to see if i and j are the same particle.
927 function mpi_cycle_jLoop(i,j) result(do_cycle)
928 !--------------- Arguments--------------------------
929 ! Index i
930 integer,intent(in) :: i
931 ! Index j
932 integer,intent(in) :: j
933 ! Result do_cycle
934 logical :: do_cycle
935 !--------------- Local variables--------------------
936 integer :: tag_i
937 integer :: tag_j
938 !--------------- END DECLARATIONS------------------
939 tag_i = tagRow(i)
940 tag_j = tagColumn(j)
941
942 do_cycle = .false.
943
944 if (tag_i == tag_j) then
945 do_cycle = .true.
946 return
947 end if
948
949 if (tag_i < tag_j) then
950 if (mod(tag_i + tag_j,2) == 0) do_cycle = .true.
951 return
952 else
953 if (mod(tag_i + tag_j,2) == 1) do_cycle = .true.
954 endif
955 end function mpi_cycle_jLoop
956 #endif
957
958 end module do_Forces