ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/UseTheForce/doForces.F90
Revision: 1877
Committed: Thu Dec 9 21:15:19 2004 UTC (19 years, 8 months ago) by tim
File size: 36809 byte(s)
Log Message:
sticky module get compiled

File Contents

# Content
1 !! doForces.F90
2 !! module doForces
3 !! Calculates Long Range forces.
4
5 !! @author Charles F. Vardeman II
6 !! @author Matthew Meineke
7 !! @version $Id: doForces.F90,v 1.5.2.3 2004-12-09 21:15:19 tim Exp $, $Date: 2004-12-09 21:15:19 $, $Name: not supported by cvs2svn $, $Revision: 1.5.2.3 $
8
9 module doForces
10 use force_globals
11 use simulation
12 use definitions
13 use atype_module
14 use switcheroo
15 use neighborLists
16 use lj
17 use sticky
18 use dipole_dipole
19 use charge_charge
20 use reaction_field
21 use gb_pair
22 use shapes
23 use vector_class
24 use eam
25 use status
26 #ifdef IS_MPI
27 use mpiSimulation
28 #endif
29
30 implicit none
31 PRIVATE
32
33 #define __FORTRAN90
34 #include "UseTheForce/fSwitchingFunction.h"
35
36 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
37 INTEGER, PARAMETER:: PAIR_LOOP = 2
38
39 logical, save :: haveRlist = .false.
40 logical, save :: haveNeighborList = .false.
41 logical, save :: haveSIMvariables = .false.
42 logical, save :: havePropertyMap = .false.
43 logical, save :: haveSaneForceField = .false.
44
45 logical, save :: FF_uses_DirectionalAtoms
46 logical, save :: FF_uses_LennardJones
47 logical, save :: FF_uses_Electrostatic
48 logical, save :: FF_uses_charges
49 logical, save :: FF_uses_dipoles
50 logical, save :: FF_uses_sticky
51 logical, save :: FF_uses_GayBerne
52 logical, save :: FF_uses_EAM
53 logical, save :: FF_uses_Shapes
54 logical, save :: FF_uses_FLARB
55 logical, save :: FF_uses_RF
56
57 logical, save :: SIM_uses_DirectionalAtoms
58 logical, save :: SIM_uses_LennardJones
59 logical, save :: SIM_uses_Electrostatics
60 logical, save :: SIM_uses_Charges
61 logical, save :: SIM_uses_Dipoles
62 logical, save :: SIM_uses_Sticky
63 logical, save :: SIM_uses_GayBerne
64 logical, save :: SIM_uses_EAM
65 logical, save :: SIM_uses_Shapes
66 logical, save :: SIM_uses_FLARB
67 logical, save :: SIM_uses_RF
68 logical, save :: SIM_requires_postpair_calc
69 logical, save :: SIM_requires_prepair_calc
70 logical, save :: SIM_uses_PBC
71 logical, save :: SIM_uses_molecular_cutoffs
72
73 real(kind=dp), save :: rlist, rlistsq
74
75 public :: init_FF
76 public :: do_force_loop
77 public :: setRlistDF
78
79 #ifdef PROFILE
80 public :: getforcetime
81 real, save :: forceTime = 0
82 real :: forceTimeInitial, forceTimeFinal
83 integer :: nLoops
84 #endif
85
86 type :: Properties
87 logical :: is_Directional = .false.
88 logical :: is_LennardJones = .false.
89 logical :: is_Electrostatic = .false.
90 logical :: is_Charge = .false.
91 logical :: is_Dipole = .false.
92 logical :: is_Sticky = .false.
93 logical :: is_GayBerne = .false.
94 logical :: is_EAM = .false.
95 logical :: is_Shape = .false.
96 logical :: is_FLARB = .false.
97 end type Properties
98
99 type(Properties), dimension(:),allocatable :: PropertyMap
100
101 contains
102
103 subroutine setRlistDF( this_rlist )
104
105 real(kind=dp) :: this_rlist
106
107 rlist = this_rlist
108 rlistsq = rlist * rlist
109
110 haveRlist = .true.
111
112 end subroutine setRlistDF
113
114 subroutine createPropertyMap(status)
115 integer :: nAtypes
116 integer :: status
117 integer :: i
118 logical :: thisProperty
119 real (kind=DP) :: thisDPproperty
120
121 status = 0
122
123 nAtypes = getSize(atypes)
124
125 if (nAtypes == 0) then
126 status = -1
127 return
128 end if
129
130 if (.not. allocated(PropertyMap)) then
131 allocate(PropertyMap(nAtypes))
132 endif
133
134 do i = 1, nAtypes
135 call getElementProperty(atypes, i, "is_Directional", thisProperty)
136 PropertyMap(i)%is_Directional = thisProperty
137
138 call getElementProperty(atypes, i, "is_LennardJones", thisProperty)
139 PropertyMap(i)%is_LennardJones = thisProperty
140
141 call getElementProperty(atypes, i, "is_Electrostatic", thisProperty)
142 PropertyMap(i)%is_Electrostatic = thisProperty
143
144 call getElementProperty(atypes, i, "is_Charge", thisProperty)
145 PropertyMap(i)%is_Charge = thisProperty
146
147 call getElementProperty(atypes, i, "is_Dipole", thisProperty)
148 PropertyMap(i)%is_Dipole = thisProperty
149
150 call getElementProperty(atypes, i, "is_Sticky", thisProperty)
151 PropertyMap(i)%is_Sticky = thisProperty
152
153 call getElementProperty(atypes, i, "is_GayBerne", thisProperty)
154 PropertyMap(i)%is_GayBerne = thisProperty
155
156 call getElementProperty(atypes, i, "is_EAM", thisProperty)
157 PropertyMap(i)%is_EAM = thisProperty
158
159 call getElementProperty(atypes, i, "is_Shape", thisProperty)
160 PropertyMap(i)%is_Shape = thisProperty
161
162 call getElementProperty(atypes, i, "is_FLARB", thisProperty)
163 PropertyMap(i)%is_FLARB = thisProperty
164 end do
165
166 havePropertyMap = .true.
167
168 end subroutine createPropertyMap
169
170 subroutine setSimVariables()
171 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
172 SIM_uses_LennardJones = SimUsesLennardJones()
173 SIM_uses_Electrostatics = SimUsesElectrostatics()
174 SIM_uses_Charges = SimUsesCharges()
175 SIM_uses_Dipoles = SimUsesDipoles()
176 SIM_uses_Sticky = SimUsesSticky()
177 SIM_uses_GayBerne = SimUsesGayBerne()
178 SIM_uses_EAM = SimUsesEAM()
179 SIM_uses_Shapes = SimUsesShapes()
180 SIM_uses_FLARB = SimUsesFLARB()
181 SIM_uses_RF = SimUsesRF()
182 SIM_requires_postpair_calc = SimRequiresPostpairCalc()
183 SIM_requires_prepair_calc = SimRequiresPrepairCalc()
184 SIM_uses_PBC = SimUsesPBC()
185
186 haveSIMvariables = .true.
187
188 return
189 end subroutine setSimVariables
190
191 subroutine doReadyCheck(error)
192 integer, intent(out) :: error
193
194 integer :: myStatus
195
196 error = 0
197
198 if (.not. havePropertyMap) then
199
200 myStatus = 0
201
202 call createPropertyMap(myStatus)
203
204 if (myStatus .ne. 0) then
205 write(default_error, *) 'createPropertyMap failed in doForces!'
206 error = -1
207 return
208 endif
209 endif
210
211 if (.not. haveSIMvariables) then
212 call setSimVariables()
213 endif
214
215 if (.not. haveRlist) then
216 write(default_error, *) 'rList has not been set in doForces!'
217 error = -1
218 return
219 endif
220
221 if (.not. haveNeighborList) then
222 write(default_error, *) 'neighbor list has not been initialized in doForces!'
223 error = -1
224 return
225 end if
226
227 if (.not. haveSaneForceField) then
228 write(default_error, *) 'Force Field is not sane in doForces!'
229 error = -1
230 return
231 end if
232
233 #ifdef IS_MPI
234 if (.not. isMPISimSet()) then
235 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
236 error = -1
237 return
238 endif
239 #endif
240 return
241 end subroutine doReadyCheck
242
243
244 subroutine init_FF(use_RF_c, thisStat)
245
246 logical, intent(in) :: use_RF_c
247
248 integer, intent(out) :: thisStat
249 integer :: my_status, nMatches
250 integer, pointer :: MatchList(:) => null()
251 real(kind=dp) :: rcut, rrf, rt, dielect
252
253 !! assume things are copacetic, unless they aren't
254 thisStat = 0
255
256 !! Fortran's version of a cast:
257 FF_uses_RF = use_RF_c
258
259 !! init_FF is called *after* all of the atom types have been
260 !! defined in atype_module using the new_atype subroutine.
261 !!
262 !! this will scan through the known atypes and figure out what
263 !! interactions are used by the force field.
264
265 FF_uses_DirectionalAtoms = .false.
266 FF_uses_LennardJones = .false.
267 FF_uses_Electrostatic = .false.
268 FF_uses_Charges = .false.
269 FF_uses_Dipoles = .false.
270 FF_uses_Sticky = .false.
271 FF_uses_GayBerne = .false.
272 FF_uses_EAM = .false.
273 FF_uses_Shapes = .false.
274 FF_uses_FLARB = .false.
275
276 call getMatchingElementList(atypes, "is_Directional", .true., &
277 nMatches, MatchList)
278 if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true.
279
280 call getMatchingElementList(atypes, "is_LennardJones", .true., &
281 nMatches, MatchList)
282 if (nMatches .gt. 0) FF_uses_LennardJones = .true.
283
284 call getMatchingElementList(atypes, "is_Electrostatic", .true., &
285 nMatches, MatchList)
286 if (nMatches .gt. 0) then
287 FF_uses_Electrostatic = .true.
288 endif
289
290 call getMatchingElementList(atypes, "is_Charge", .true., &
291 nMatches, MatchList)
292 if (nMatches .gt. 0) then
293 FF_uses_charges = .true.
294 FF_uses_electrostatic = .true.
295 endif
296
297 call getMatchingElementList(atypes, "is_Dipole", .true., &
298 nMatches, MatchList)
299 if (nMatches .gt. 0) then
300 FF_uses_dipoles = .true.
301 FF_uses_electrostatic = .true.
302 FF_uses_DirectionalAtoms = .true.
303 endif
304
305 call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
306 MatchList)
307 if (nMatches .gt. 0) then
308 FF_uses_Sticky = .true.
309 FF_uses_DirectionalAtoms = .true.
310 endif
311
312 call getMatchingElementList(atypes, "is_GayBerne", .true., &
313 nMatches, MatchList)
314 if (nMatches .gt. 0) then
315 FF_uses_GayBerne = .true.
316 FF_uses_DirectionalAtoms = .true.
317 endif
318
319 call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
320 if (nMatches .gt. 0) FF_uses_EAM = .true.
321
322 call getMatchingElementList(atypes, "is_Shape", .true., &
323 nMatches, MatchList)
324 if (nMatches .gt. 0) then
325 FF_uses_Shapes = .true.
326 FF_uses_DirectionalAtoms = .true.
327 endif
328
329 call getMatchingElementList(atypes, "is_FLARB", .true., &
330 nMatches, MatchList)
331 if (nMatches .gt. 0) FF_uses_FLARB = .true.
332
333 !! Assume sanity (for the sake of argument)
334 haveSaneForceField = .true.
335
336 !! check to make sure the FF_uses_RF setting makes sense
337
338 if (FF_uses_dipoles) then
339 if (FF_uses_RF) then
340 dielect = getDielect()
341 call initialize_rf(dielect)
342 endif
343 else
344 if (FF_uses_RF) then
345 write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
346 thisStat = -1
347 haveSaneForceField = .false.
348 return
349 endif
350 endif
351
352 !sticky module does not contain check_sticky_FF anymore
353 !if (FF_uses_sticky) then
354 ! call check_sticky_FF(my_status)
355 ! if (my_status /= 0) then
356 ! thisStat = -1
357 ! haveSaneForceField = .false.
358 ! return
359 ! end if
360 !endif
361
362 if (FF_uses_EAM) then
363 call init_EAM_FF(my_status)
364 if (my_status /= 0) then
365 write(default_error, *) "init_EAM_FF returned a bad status"
366 thisStat = -1
367 haveSaneForceField = .false.
368 return
369 end if
370 endif
371
372 if (FF_uses_GayBerne) then
373 call check_gb_pair_FF(my_status)
374 if (my_status .ne. 0) then
375 thisStat = -1
376 haveSaneForceField = .false.
377 return
378 endif
379 endif
380
381 if (FF_uses_GayBerne .and. FF_uses_LennardJones) then
382 endif
383
384 if (.not. haveNeighborList) then
385 !! Create neighbor lists
386 call expandNeighborList(nLocal, my_status)
387 if (my_Status /= 0) then
388 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
389 thisStat = -1
390 return
391 endif
392 haveNeighborList = .true.
393 endif
394
395 end subroutine init_FF
396
397
398 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
399 !------------------------------------------------------------->
400 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
401 do_pot_c, do_stress_c, error)
402 !! Position array provided by C, dimensioned by getNlocal
403 real ( kind = dp ), dimension(3, nLocal) :: q
404 !! molecular center-of-mass position array
405 real ( kind = dp ), dimension(3, nGroups) :: q_group
406 !! Rotation Matrix for each long range particle in simulation.
407 real( kind = dp), dimension(9, nLocal) :: A
408 !! Unit vectors for dipoles (lab frame)
409 real( kind = dp ), dimension(9,nLocal) :: eFrame
410 !! Force array provided by C, dimensioned by getNlocal
411 real ( kind = dp ), dimension(3,nLocal) :: f
412 !! Torsion array provided by C, dimensioned by getNlocal
413 real( kind = dp ), dimension(3,nLocal) :: t
414
415 !! Stress Tensor
416 real( kind = dp), dimension(9) :: tau
417 real ( kind = dp ) :: pot
418 logical ( kind = 2) :: do_pot_c, do_stress_c
419 logical :: do_pot
420 logical :: do_stress
421 logical :: in_switching_region
422 #ifdef IS_MPI
423 real( kind = DP ) :: pot_local
424 integer :: nAtomsInRow
425 integer :: nAtomsInCol
426 integer :: nprocs
427 integer :: nGroupsInRow
428 integer :: nGroupsInCol
429 #endif
430 integer :: natoms
431 logical :: update_nlist
432 integer :: i, j, jstart, jend, jnab
433 integer :: istart, iend
434 integer :: ia, jb, atom1, atom2
435 integer :: nlist
436 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
437 real( kind = DP ) :: sw, dswdr, swderiv, mf
438 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij
439 real(kind=dp) :: rfpot, mu_i, virial
440 integer :: me_i, me_j, n_in_i, n_in_j
441 logical :: is_dp_i
442 integer :: neighborListSize
443 integer :: listerror, error
444 integer :: localError
445 integer :: propPack_i, propPack_j
446 integer :: loopStart, loopEnd, loop
447
448 real(kind=dp) :: listSkin = 1.0
449
450 !! initialize local variables
451
452 #ifdef IS_MPI
453 pot_local = 0.0_dp
454 nAtomsInRow = getNatomsInRow(plan_atom_row)
455 nAtomsInCol = getNatomsInCol(plan_atom_col)
456 nGroupsInRow = getNgroupsInRow(plan_group_row)
457 nGroupsInCol = getNgroupsInCol(plan_group_col)
458 #else
459 natoms = nlocal
460 #endif
461
462 call doReadyCheck(localError)
463 if ( localError .ne. 0 ) then
464 call handleError("do_force_loop", "Not Initialized")
465 error = -1
466 return
467 end if
468 call zero_work_arrays()
469
470 do_pot = do_pot_c
471 do_stress = do_stress_c
472
473 ! Gather all information needed by all force loops:
474
475 #ifdef IS_MPI
476
477 call gather(q, q_Row, plan_atom_row_3d)
478 call gather(q, q_Col, plan_atom_col_3d)
479
480 call gather(q_group, q_group_Row, plan_group_row_3d)
481 call gather(q_group, q_group_Col, plan_group_col_3d)
482
483 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
484 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
485 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
486
487 call gather(A, A_Row, plan_atom_row_rotation)
488 call gather(A, A_Col, plan_atom_col_rotation)
489 endif
490
491 #endif
492
493 !! Begin force loop timing:
494 #ifdef PROFILE
495 call cpu_time(forceTimeInitial)
496 nloops = nloops + 1
497 #endif
498
499 loopEnd = PAIR_LOOP
500 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
501 loopStart = PREPAIR_LOOP
502 else
503 loopStart = PAIR_LOOP
504 endif
505
506 do loop = loopStart, loopEnd
507
508 ! See if we need to update neighbor lists
509 ! (but only on the first time through):
510 if (loop .eq. loopStart) then
511 #ifdef IS_MPI
512 call checkNeighborList(nGroupsInRow, q_group_row, listSkin, &
513 update_nlist)
514 #else
515 call checkNeighborList(nGroups, q_group, listSkin, &
516 update_nlist)
517 #endif
518 endif
519
520 if (update_nlist) then
521 !! save current configuration and construct neighbor list
522 #ifdef IS_MPI
523 call saveNeighborList(nGroupsInRow, q_group_row)
524 #else
525 call saveNeighborList(nGroups, q_group)
526 #endif
527 neighborListSize = size(list)
528 nlist = 0
529 endif
530
531 istart = 1
532 #ifdef IS_MPI
533 iend = nGroupsInRow
534 #else
535 iend = nGroups - 1
536 #endif
537 outer: do i = istart, iend
538
539 if (update_nlist) point(i) = nlist + 1
540
541 n_in_i = groupStartRow(i+1) - groupStartRow(i)
542
543 if (update_nlist) then
544 #ifdef IS_MPI
545 jstart = 1
546 jend = nGroupsInCol
547 #else
548 jstart = i+1
549 jend = nGroups
550 #endif
551 else
552 jstart = point(i)
553 jend = point(i+1) - 1
554 ! make sure group i has neighbors
555 if (jstart .gt. jend) cycle outer
556 endif
557
558 do jnab = jstart, jend
559 if (update_nlist) then
560 j = jnab
561 else
562 j = list(jnab)
563 endif
564
565 #ifdef IS_MPI
566 call get_interatomic_vector(q_group_Row(:,i), &
567 q_group_Col(:,j), d_grp, rgrpsq)
568 #else
569 call get_interatomic_vector(q_group(:,i), &
570 q_group(:,j), d_grp, rgrpsq)
571 #endif
572
573 if (rgrpsq < rlistsq) then
574 if (update_nlist) then
575 nlist = nlist + 1
576
577 if (nlist > neighborListSize) then
578 #ifdef IS_MPI
579 call expandNeighborList(nGroupsInRow, listerror)
580 #else
581 call expandNeighborList(nGroups, listerror)
582 #endif
583 if (listerror /= 0) then
584 error = -1
585 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
586 return
587 end if
588 neighborListSize = size(list)
589 endif
590
591 list(nlist) = j
592 endif
593
594 if (loop .eq. PAIR_LOOP) then
595 vij = 0.0d0
596 fij(1:3) = 0.0d0
597 endif
598
599 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
600 in_switching_region)
601
602 n_in_j = groupStartCol(j+1) - groupStartCol(j)
603
604 do ia = groupStartRow(i), groupStartRow(i+1)-1
605
606 atom1 = groupListRow(ia)
607
608 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
609
610 atom2 = groupListCol(jb)
611
612 if (skipThisPair(atom1, atom2)) cycle inner
613
614 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
615 d_atm(1:3) = d_grp(1:3)
616 ratmsq = rgrpsq
617 else
618 #ifdef IS_MPI
619 call get_interatomic_vector(q_Row(:,atom1), &
620 q_Col(:,atom2), d_atm, ratmsq)
621 #else
622 call get_interatomic_vector(q(:,atom1), &
623 q(:,atom2), d_atm, ratmsq)
624 #endif
625 endif
626
627 if (loop .eq. PREPAIR_LOOP) then
628 #ifdef IS_MPI
629 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
630 rgrpsq, d_grp, do_pot, do_stress, &
631 eFrame, A, f, t, pot_local)
632 #else
633 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
634 rgrpsq, d_grp, do_pot, do_stress, &
635 eFrame, A, f, t, pot)
636 #endif
637 else
638 #ifdef IS_MPI
639 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
640 do_pot, &
641 eFrame, A, f, t, pot_local, vpair, fpair)
642 #else
643 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
644 do_pot, &
645 eFrame, A, f, t, pot, vpair, fpair)
646 #endif
647
648 vij = vij + vpair
649 fij(1:3) = fij(1:3) + fpair(1:3)
650 endif
651 enddo inner
652 enddo
653
654 if (loop .eq. PAIR_LOOP) then
655 if (in_switching_region) then
656 swderiv = vij*dswdr/rgrp
657 fij(1) = fij(1) + swderiv*d_grp(1)
658 fij(2) = fij(2) + swderiv*d_grp(2)
659 fij(3) = fij(3) + swderiv*d_grp(3)
660
661 do ia=groupStartRow(i), groupStartRow(i+1)-1
662 atom1=groupListRow(ia)
663 mf = mfactRow(atom1)
664 #ifdef IS_MPI
665 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
666 f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
667 f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
668 #else
669 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
670 f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
671 f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
672 #endif
673 enddo
674
675 do jb=groupStartCol(j), groupStartCol(j+1)-1
676 atom2=groupListCol(jb)
677 mf = mfactCol(atom2)
678 #ifdef IS_MPI
679 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
680 f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
681 f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
682 #else
683 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
684 f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
685 f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
686 #endif
687 enddo
688 endif
689
690 if (do_stress) call add_stress_tensor(d_grp, fij)
691 endif
692 end if
693 enddo
694 enddo outer
695
696 if (update_nlist) then
697 #ifdef IS_MPI
698 point(nGroupsInRow + 1) = nlist + 1
699 #else
700 point(nGroups) = nlist + 1
701 #endif
702 if (loop .eq. PREPAIR_LOOP) then
703 ! we just did the neighbor list update on the first
704 ! pass, so we don't need to do it
705 ! again on the second pass
706 update_nlist = .false.
707 endif
708 endif
709
710 if (loop .eq. PREPAIR_LOOP) then
711 call do_preforce(nlocal, pot)
712 endif
713
714 enddo
715
716 !! Do timing
717 #ifdef PROFILE
718 call cpu_time(forceTimeFinal)
719 forceTime = forceTime + forceTimeFinal - forceTimeInitial
720 #endif
721
722 #ifdef IS_MPI
723 !!distribute forces
724
725 f_temp = 0.0_dp
726 call scatter(f_Row,f_temp,plan_atom_row_3d)
727 do i = 1,nlocal
728 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
729 end do
730
731 f_temp = 0.0_dp
732 call scatter(f_Col,f_temp,plan_atom_col_3d)
733 do i = 1,nlocal
734 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
735 end do
736
737 if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then
738 t_temp = 0.0_dp
739 call scatter(t_Row,t_temp,plan_atom_row_3d)
740 do i = 1,nlocal
741 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
742 end do
743 t_temp = 0.0_dp
744 call scatter(t_Col,t_temp,plan_atom_col_3d)
745
746 do i = 1,nlocal
747 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
748 end do
749 endif
750
751 if (do_pot) then
752 ! scatter/gather pot_row into the members of my column
753 call scatter(pot_Row, pot_Temp, plan_atom_row)
754
755 ! scatter/gather pot_local into all other procs
756 ! add resultant to get total pot
757 do i = 1, nlocal
758 pot_local = pot_local + pot_Temp(i)
759 enddo
760
761 pot_Temp = 0.0_DP
762
763 call scatter(pot_Col, pot_Temp, plan_atom_col)
764 do i = 1, nlocal
765 pot_local = pot_local + pot_Temp(i)
766 enddo
767
768 endif
769 #endif
770
771 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
772
773 if (FF_uses_RF .and. SIM_uses_RF) then
774
775 #ifdef IS_MPI
776 call scatter(rf_Row,rf,plan_atom_row_3d)
777 call scatter(rf_Col,rf_Temp,plan_atom_col_3d)
778 do i = 1,nlocal
779 rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
780 end do
781 #endif
782
783 do i = 1, nLocal
784
785 rfpot = 0.0_DP
786 #ifdef IS_MPI
787 me_i = atid_row(i)
788 #else
789 me_i = atid(i)
790 #endif
791
792 if (PropertyMap(me_i)%is_Dipole) then
793
794 mu_i = getDipoleMoment(me_i)
795
796 !! The reaction field needs to include a self contribution
797 !! to the field:
798 call accumulate_self_rf(i, mu_i, eFrame)
799 !! Get the reaction field contribution to the
800 !! potential and torques:
801 call reaction_field_final(i, mu_i, eFrame, rfpot, t, do_pot)
802 #ifdef IS_MPI
803 pot_local = pot_local + rfpot
804 #else
805 pot = pot + rfpot
806
807 #endif
808 endif
809 enddo
810 endif
811 endif
812
813
814 #ifdef IS_MPI
815
816 if (do_pot) then
817 pot = pot + pot_local
818 !! we assume the c code will do the allreduce to get the total potential
819 !! we could do it right here if we needed to...
820 endif
821
822 if (do_stress) then
823 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
824 mpi_comm_world,mpi_err)
825 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
826 mpi_comm_world,mpi_err)
827 endif
828
829 #else
830
831 if (do_stress) then
832 tau = tau_Temp
833 virial = virial_Temp
834 endif
835
836 #endif
837
838 end subroutine do_force_loop
839
840 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
841 eFrame, A, f, t, pot, vpair, fpair)
842
843 real( kind = dp ) :: pot, vpair, sw
844 real( kind = dp ), dimension(3) :: fpair
845 real( kind = dp ), dimension(nLocal) :: mfact
846 real( kind = dp ), dimension(9,nLocal) :: eFrame
847 real( kind = dp ), dimension(9,nLocal) :: A
848 real( kind = dp ), dimension(3,nLocal) :: f
849 real( kind = dp ), dimension(3,nLocal) :: t
850
851 logical, intent(inout) :: do_pot
852 integer, intent(in) :: i, j
853 real ( kind = dp ), intent(inout) :: rijsq
854 real ( kind = dp ) :: r
855 real ( kind = dp ), intent(inout) :: d(3)
856 integer :: me_i, me_j
857
858 r = sqrt(rijsq)
859 vpair = 0.0d0
860 fpair(1:3) = 0.0d0
861
862 #ifdef IS_MPI
863 me_i = atid_row(i)
864 me_j = atid_col(j)
865 #else
866 me_i = atid(i)
867 me_j = atid(j)
868 #endif
869
870 if (FF_uses_LennardJones .and. SIM_uses_LennardJones) then
871
872 if ( PropertyMap(me_i)%is_LennardJones .and. &
873 PropertyMap(me_j)%is_LennardJones ) then
874 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
875 endif
876
877 endif
878
879 if (FF_uses_charges .and. SIM_uses_charges) then
880
881 if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
882 call do_charge_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
883 pot, f, do_pot)
884 endif
885
886 endif
887
888 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
889
890 if ( PropertyMap(me_i)%is_Dipole .and. PropertyMap(me_j)%is_Dipole) then
891 call do_dipole_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
892 pot, eFrame, f, t, do_pot)
893 if (FF_uses_RF .and. SIM_uses_RF) then
894 call accumulate_rf(i, j, r, eFrame, sw)
895 call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair)
896 endif
897 endif
898
899 endif
900
901 if (FF_uses_Sticky .and. SIM_uses_sticky) then
902
903 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
904 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
905 pot, A, f, t, do_pot)
906 endif
907
908 endif
909
910
911 if (FF_uses_GayBerne .and. SIM_uses_GayBerne) then
912
913 if ( PropertyMap(me_i)%is_GayBerne .and. &
914 PropertyMap(me_j)%is_GayBerne) then
915 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
916 pot, A, f, t, do_pot)
917 endif
918
919 endif
920
921 if (FF_uses_EAM .and. SIM_uses_EAM) then
922
923 if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
924 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
925 do_pot)
926 endif
927
928 endif
929
930 if (FF_uses_Shapes .and. SIM_uses_Shapes) then
931
932 if ( PropertyMap(me_i)%is_Shape .and. &
933 PropertyMap(me_j)%is_Shape ) then
934 call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, &
935 pot, A, f, t, do_pot)
936 endif
937
938 endif
939
940 end subroutine do_pair
941
942 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
943 do_pot, do_stress, eFrame, A, f, t, pot)
944
945 real( kind = dp ) :: pot, sw
946 real( kind = dp ), dimension(9,nLocal) :: eFrame
947 real (kind=dp), dimension(9,nLocal) :: A
948 real (kind=dp), dimension(3,nLocal) :: f
949 real (kind=dp), dimension(3,nLocal) :: t
950
951 logical, intent(inout) :: do_pot, do_stress
952 integer, intent(in) :: i, j
953 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
954 real ( kind = dp ) :: r, rc
955 real ( kind = dp ), intent(inout) :: d(3), dc(3)
956
957 logical :: is_EAM_i, is_EAM_j
958
959 integer :: me_i, me_j
960
961
962 r = sqrt(rijsq)
963 if (SIM_uses_molecular_cutoffs) then
964 rc = sqrt(rcijsq)
965 else
966 rc = r
967 endif
968
969
970 #ifdef IS_MPI
971 me_i = atid_row(i)
972 me_j = atid_col(j)
973 #else
974 me_i = atid(i)
975 me_j = atid(j)
976 #endif
977
978 if (FF_uses_EAM .and. SIM_uses_EAM) then
979
980 if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
981 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
982
983 endif
984
985 end subroutine do_prepair
986
987
988 subroutine do_preforce(nlocal,pot)
989 integer :: nlocal
990 real( kind = dp ) :: pot
991
992 if (FF_uses_EAM .and. SIM_uses_EAM) then
993 call calc_EAM_preforce_Frho(nlocal,pot)
994 endif
995
996
997 end subroutine do_preforce
998
999
1000 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1001
1002 real (kind = dp), dimension(3) :: q_i
1003 real (kind = dp), dimension(3) :: q_j
1004 real ( kind = dp ), intent(out) :: r_sq
1005 real( kind = dp ) :: d(3), scaled(3)
1006 integer i
1007
1008 d(1:3) = q_j(1:3) - q_i(1:3)
1009
1010 ! Wrap back into periodic box if necessary
1011 if ( SIM_uses_PBC ) then
1012
1013 if( .not.boxIsOrthorhombic ) then
1014 ! calc the scaled coordinates.
1015
1016 scaled = matmul(HmatInv, d)
1017
1018 ! wrap the scaled coordinates
1019
1020 scaled = scaled - anint(scaled)
1021
1022
1023 ! calc the wrapped real coordinates from the wrapped scaled
1024 ! coordinates
1025
1026 d = matmul(Hmat,scaled)
1027
1028 else
1029 ! calc the scaled coordinates.
1030
1031 do i = 1, 3
1032 scaled(i) = d(i) * HmatInv(i,i)
1033
1034 ! wrap the scaled coordinates
1035
1036 scaled(i) = scaled(i) - anint(scaled(i))
1037
1038 ! calc the wrapped real coordinates from the wrapped scaled
1039 ! coordinates
1040
1041 d(i) = scaled(i)*Hmat(i,i)
1042 enddo
1043 endif
1044
1045 endif
1046
1047 r_sq = dot_product(d,d)
1048
1049 end subroutine get_interatomic_vector
1050
1051 subroutine zero_work_arrays()
1052
1053 #ifdef IS_MPI
1054
1055 q_Row = 0.0_dp
1056 q_Col = 0.0_dp
1057
1058 q_group_Row = 0.0_dp
1059 q_group_Col = 0.0_dp
1060
1061 eFrame_Row = 0.0_dp
1062 eFrame_Col = 0.0_dp
1063
1064 A_Row = 0.0_dp
1065 A_Col = 0.0_dp
1066
1067 f_Row = 0.0_dp
1068 f_Col = 0.0_dp
1069 f_Temp = 0.0_dp
1070
1071 t_Row = 0.0_dp
1072 t_Col = 0.0_dp
1073 t_Temp = 0.0_dp
1074
1075 pot_Row = 0.0_dp
1076 pot_Col = 0.0_dp
1077 pot_Temp = 0.0_dp
1078
1079 rf_Row = 0.0_dp
1080 rf_Col = 0.0_dp
1081 rf_Temp = 0.0_dp
1082
1083 #endif
1084
1085 if (FF_uses_EAM .and. SIM_uses_EAM) then
1086 call clean_EAM()
1087 endif
1088
1089 rf = 0.0_dp
1090 tau_Temp = 0.0_dp
1091 virial_Temp = 0.0_dp
1092 end subroutine zero_work_arrays
1093
1094 function skipThisPair(atom1, atom2) result(skip_it)
1095 integer, intent(in) :: atom1
1096 integer, intent(in), optional :: atom2
1097 logical :: skip_it
1098 integer :: unique_id_1, unique_id_2
1099 integer :: me_i,me_j
1100 integer :: i
1101
1102 skip_it = .false.
1103
1104 !! there are a number of reasons to skip a pair or a particle
1105 !! mostly we do this to exclude atoms who are involved in short
1106 !! range interactions (bonds, bends, torsions), but we also need
1107 !! to exclude some overcounted interactions that result from
1108 !! the parallel decomposition
1109
1110 #ifdef IS_MPI
1111 !! in MPI, we have to look up the unique IDs for each atom
1112 unique_id_1 = AtomRowToGlobal(atom1)
1113 #else
1114 !! in the normal loop, the atom numbers are unique
1115 unique_id_1 = atom1
1116 #endif
1117
1118 !! We were called with only one atom, so just check the global exclude
1119 !! list for this atom
1120 if (.not. present(atom2)) then
1121 do i = 1, nExcludes_global
1122 if (excludesGlobal(i) == unique_id_1) then
1123 skip_it = .true.
1124 return
1125 end if
1126 end do
1127 return
1128 end if
1129
1130 #ifdef IS_MPI
1131 unique_id_2 = AtomColToGlobal(atom2)
1132 #else
1133 unique_id_2 = atom2
1134 #endif
1135
1136 #ifdef IS_MPI
1137 !! this situation should only arise in MPI simulations
1138 if (unique_id_1 == unique_id_2) then
1139 skip_it = .true.
1140 return
1141 end if
1142
1143 !! this prevents us from doing the pair on multiple processors
1144 if (unique_id_1 < unique_id_2) then
1145 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1146 skip_it = .true.
1147 return
1148 endif
1149 else
1150 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1151 skip_it = .true.
1152 return
1153 endif
1154 endif
1155 #endif
1156
1157 !! the rest of these situations can happen in all simulations:
1158 do i = 1, nExcludes_global
1159 if ((excludesGlobal(i) == unique_id_1) .or. &
1160 (excludesGlobal(i) == unique_id_2)) then
1161 skip_it = .true.
1162 return
1163 endif
1164 enddo
1165
1166 do i = 1, nSkipsForAtom(atom1)
1167 if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1168 skip_it = .true.
1169 return
1170 endif
1171 end do
1172
1173 return
1174 end function skipThisPair
1175
1176 function FF_UsesDirectionalAtoms() result(doesit)
1177 logical :: doesit
1178 doesit = FF_uses_DirectionalAtoms .or. FF_uses_Dipoles .or. &
1179 FF_uses_Sticky .or. FF_uses_GayBerne .or. FF_uses_Shapes
1180 end function FF_UsesDirectionalAtoms
1181
1182 function FF_RequiresPrepairCalc() result(doesit)
1183 logical :: doesit
1184 doesit = FF_uses_EAM
1185 end function FF_RequiresPrepairCalc
1186
1187 function FF_RequiresPostpairCalc() result(doesit)
1188 logical :: doesit
1189 doesit = FF_uses_RF
1190 end function FF_RequiresPostpairCalc
1191
1192 #ifdef PROFILE
1193 function getforcetime() result(totalforcetime)
1194 real(kind=dp) :: totalforcetime
1195 totalforcetime = forcetime
1196 end function getforcetime
1197 #endif
1198
1199 !! This cleans componets of force arrays belonging only to fortran
1200
1201 subroutine add_stress_tensor(dpair, fpair)
1202
1203 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1204
1205 ! because the d vector is the rj - ri vector, and
1206 ! because fx, fy, fz are the force on atom i, we need a
1207 ! negative sign here:
1208
1209 tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1210 tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1211 tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1212 tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1213 tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1214 tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1215 tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1216 tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1217 tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1218
1219 virial_Temp = virial_Temp + &
1220 (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1221
1222 end subroutine add_stress_tensor
1223
1224 end module doForces
1225
1226 !! Interfaces for C programs to module....
1227
1228 subroutine initFortranFF(use_RF_c, thisStat)
1229 use doForces, ONLY: init_FF
1230 logical, intent(in) :: use_RF_c
1231
1232 integer, intent(out) :: thisStat
1233 call init_FF(use_RF_c, thisStat)
1234
1235 end subroutine initFortranFF
1236
1237 subroutine doForceloop(q, q_group, A, eFrame, f, t, tau, pot, &
1238 do_pot_c, do_stress_c, error)
1239
1240 use definitions, ONLY: dp
1241 use simulation
1242 use doForces, ONLY: do_force_loop
1243 !! Position array provided by C, dimensioned by getNlocal
1244 real ( kind = dp ), dimension(3, nLocal) :: q
1245 !! molecular center-of-mass position array
1246 real ( kind = dp ), dimension(3, nGroups) :: q_group
1247 !! Rotation Matrix for each long range particle in simulation.
1248 real( kind = dp), dimension(9, nLocal) :: A
1249 !! Unit vectors for dipoles (lab frame)
1250 real( kind = dp ), dimension(9,nLocal) :: eFrame
1251 !! Force array provided by C, dimensioned by getNlocal
1252 real ( kind = dp ), dimension(3,nLocal) :: f
1253 !! Torsion array provided by C, dimensioned by getNlocal
1254 real( kind = dp ), dimension(3,nLocal) :: t
1255
1256 !! Stress Tensor
1257 real( kind = dp), dimension(9) :: tau
1258 real ( kind = dp ) :: pot
1259 logical ( kind = 2) :: do_pot_c, do_stress_c
1260 integer :: error
1261
1262 call do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, &
1263 do_pot_c, do_stress_c, error)
1264
1265 end subroutine doForceloop