ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/doForces.F90
Revision: 1528
Committed: Fri Dec 17 20:11:05 2010 UTC (15 years, 1 month ago) by gezelter
File size: 40745 byte(s)
Log Message:
Doesn't build yet, but getting much closer...


File Contents

# Content
1 !!
2 !! Copyright (c) 2005, 2009 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41
42 !! doForces.F90
43 !! module doForces
44 !! Calculates Long Range forces.
45
46 !! @author Charles F. Vardeman II
47 !! @author Matthew Meineke
48 !! @version $Id$, $Date$, $Name: not supported by cvs2svn $, $Revision$
49
50
51 module doForces
52 use force_globals
53 use fForceOptions
54 use simulation
55 use definitions
56 use switcheroo
57 use neighborLists
58 use vector_class
59 use status
60 use ISO_C_BINDING
61
62 #ifdef IS_MPI
63 use mpiSimulation
64 #endif
65
66 implicit none
67 PRIVATE
68
69 real(kind=dp), external :: get_cutoff
70
71 #define __FORTRAN90
72 #include "UseTheForce/DarkSide/fInteractionMap.h"
73
74 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
75 INTEGER, PARAMETER:: PAIR_LOOP = 2
76
77 logical, save :: haveNeighborList = .false.
78 logical, save :: haveSIMvariables = .false.
79 logical, save :: haveCutoffs = .false.
80 logical, save :: haveSkinThickness = .false.
81
82 logical, save :: SIM_uses_DirectionalAtoms
83 logical, save :: SIM_uses_MetallicAtoms
84 logical, save :: SIM_requires_skip_correction
85 logical, save :: SIM_requires_self_correction
86 logical, save :: SIM_uses_PBC
87 logical, save :: SIM_uses_AtomicVirial
88
89 real(kind=dp), save :: rCut, rSwitch, rList, rListSq, rCutSq
90 real(kind=dp), save :: skinThickness
91 logical, save :: shifted_pot
92 logical, save :: shifted_force
93
94 public :: init_FF
95 public :: setCutoffs
96 public :: setSkinThickness
97 public :: do_force_loop
98
99 #ifdef PROFILE
100 public :: getforcetime
101 real, save :: forceTime = 0
102 real :: forceTimeInitial, forceTimeFinal
103 integer :: nLoops
104 #endif
105
106 contains
107
108 subroutine setCutoffs(rc, rsw, sp, sf)
109
110 real(kind=dp),intent(in) :: rc, rsw
111 integer, intent(in) :: sp, sf
112 character(len = statusMsgSize) :: errMsg
113 integer :: localError
114
115 rCut = rc
116 rSwitch = rsw
117 rCutsq = rCut * rCut
118
119 if (haveSkinThickness) then
120 rList = rCut + skinThickness
121 rListSq = rListSq
122 endif
123
124 if (sp .ne. 0) then
125 shifted_pot = .true.
126 else
127 shifted_pot = .false.
128 endif
129 if (sf .ne. 0) then
130 shifted_force = .true.
131 else
132 shifted_force = .false.
133 endif
134
135 if (abs(rCut - rSwitch) .lt. 0.0001) then
136 if (shifted_force) then
137 write(errMsg, *) &
138 'cutoffRadius and switchingRadius are set to the', newline &
139 // tab, 'same value. OpenMD will use shifted force', newline &
140 // tab, 'potentials instead of switching functions.'
141
142 call handleInfo("setCutoffs", errMsg)
143 else
144 write(errMsg, *) &
145 'cutoffRadius and switchingRadius are set to the', newline &
146 // tab, 'same value. OpenMD will use shifted', newline &
147 // tab, 'potentials instead of switching functions.'
148
149 call handleInfo("setCutoffs", errMsg)
150
151 shifted_pot = .true.
152 endif
153
154 endif
155
156 localError = 0
157 call set_switch(rSwitch, rCut)
158 call setHmatDangerousRcutValue(rCut)
159
160 haveCutoffs = .true.
161 end subroutine setCutoffs
162
163 subroutine setSkinThickness( thisSkin )
164
165 real(kind=dp), intent(in) :: thisSkin
166
167 skinThickness = thisSkin
168 haveSkinThickness = .true.
169
170 if (haveCutoffs) then
171 rList = rCut + skinThickness
172 rListSq = rList * rList
173 endif
174
175 end subroutine setSkinThickness
176
177 subroutine setSimVariables()
178 SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms()
179 SIM_uses_MetallicAtoms = SimUsesMetallicAtoms()
180 SIM_requires_skip_correction = SimRequiresSkipCorrection()
181 SIM_requires_self_correction = SimRequiresSelfCorrection()
182 SIM_uses_PBC = SimUsesPBC()
183 SIM_uses_AtomicVirial = SimUsesAtomicVirial()
184
185 haveSIMvariables = .true.
186 return
187 end subroutine setSimVariables
188
189 subroutine doReadyCheck(error)
190 integer, intent(out) :: error
191 integer :: myStatus
192
193 error = 0
194
195 if (.not. haveSIMvariables) then
196 call setSimVariables()
197 endif
198
199 if (.not. haveCutoffs) then
200 write(default_error, *) 'cutoffs have not been set in doForces!'
201 error = -1
202 return
203 endif
204
205 if (.not. haveSkinThickness) then
206 write(default_error, *) 'skin thickness has not been set in doForces!'
207 error = -1
208 return
209 endif
210
211 if (.not. haveNeighborList) then
212 write(default_error, *) 'neighbor list has not been initialized in doForces!'
213 error = -1
214 return
215 end if
216
217
218 #ifdef IS_MPI
219 if (.not. isMPISimSet()) then
220 write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
221 error = -1
222 return
223 endif
224 #endif
225 return
226 end subroutine doReadyCheck
227
228 subroutine init_FF(thisStat)
229
230 integer :: my_status
231 integer, intent(out) :: thisStat
232
233 !! assume things are copacetic, unless they aren't
234 thisStat = 0
235
236 if (.not. haveNeighborList) then
237 !! Create neighbor lists
238 call expandNeighborList(nLocal, my_status)
239 if (my_Status /= 0) then
240 write(default_error,*) "SimSetup: ExpandNeighborList returned error."
241 thisStat = -1
242 return
243 endif
244 haveNeighborList = .true.
245 endif
246
247 end subroutine init_FF
248
249
250 !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
251 !------------------------------------------------------------->
252 subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, &
253 error)
254 !! Position array provided by C, dimensioned by getNlocal
255 real ( kind = dp ), dimension(3, nLocal) :: q
256 !! molecular center-of-mass position array
257 real ( kind = dp ), dimension(3, nGroups) :: q_group
258 !! Rotation Matrix for each long range particle in simulation.
259 real( kind = dp), dimension(9, nLocal) :: A
260 !! Unit vectors for dipoles (lab frame)
261 real( kind = dp ), dimension(9,nLocal) :: eFrame
262 !! Force array provided by C, dimensioned by getNlocal
263 real ( kind = dp ), dimension(3,nLocal) :: f
264 !! Torsion array provided by C, dimensioned by getNlocal
265 real( kind = dp ), dimension(3,nLocal) :: t
266
267 !! Stress Tensor
268 real( kind = dp), dimension(9) :: tau
269 real ( kind = dp ),dimension(LR_POT_TYPES) :: pot
270 real( kind = dp ), dimension(nLocal) :: particle_pot
271 real( kind = dp ), dimension(nLocal) :: skipped_charge
272
273 logical :: in_switching_region
274 #ifdef IS_MPI
275 real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local
276 integer :: nAtomsInRow
277 integer :: nAtomsInCol
278 integer :: nprocs
279 integer :: nGroupsInRow
280 integer :: nGroupsInCol
281 #endif
282 integer :: natoms
283 logical :: update_nlist
284 integer :: i, j, jstart, jend, jnab
285 integer :: istart, iend
286 integer :: ia, jb, atom1, atom2
287 integer :: nlist
288 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij
289 real( kind = DP ) :: sw, dswdr, swderiv, mf
290 real( kind = DP ) :: rVal
291 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag
292 real(kind=dp) :: rfpot, mu_i
293 real(kind=dp):: rCut
294 integer :: n_in_i, n_in_j, iG, j1
295 logical :: is_dp_i
296 integer :: neighborListSize
297 integer :: listerror, error
298 integer :: localError
299 integer :: propPack_i, propPack_j
300 integer :: loopStart, loopEnd, loop
301 integer :: i1, topoDist
302
303 real(kind=dp) :: skch
304
305 !! initialize local variables
306
307 #ifdef IS_MPI
308 pot_local = 0.0_dp
309 nAtomsInRow = getNatomsInRow(plan_atom_row)
310 nAtomsInCol = getNatomsInCol(plan_atom_col)
311 nGroupsInRow = getNgroupsInRow(plan_group_row)
312 nGroupsInCol = getNgroupsInCol(plan_group_col)
313 #else
314 natoms = nlocal
315 #endif
316
317 call doReadyCheck(localError)
318 if ( localError .ne. 0 ) then
319 call handleError("do_force_loop", "Not Initialized")
320 error = -1
321 return
322 end if
323 call zero_work_arrays()
324
325 ! Gather all information needed by all force loops:
326
327 #ifdef IS_MPI
328
329 call gather(q, q_Row, plan_atom_row_3d)
330 call gather(q, q_Col, plan_atom_col_3d)
331
332 call gather(q_group, q_group_Row, plan_group_row_3d)
333 call gather(q_group, q_group_Col, plan_group_col_3d)
334
335 if (SIM_uses_DirectionalAtoms) then
336 call gather(eFrame, eFrame_Row, plan_atom_row_rotation)
337 call gather(eFrame, eFrame_Col, plan_atom_col_rotation)
338
339 call gather(A, A_Row, plan_atom_row_rotation)
340 call gather(A, A_Col, plan_atom_col_rotation)
341 endif
342
343 #endif
344
345 !! Begin force loop timing:
346 #ifdef PROFILE
347 call cpu_time(forceTimeInitial)
348 nloops = nloops + 1
349 #endif
350
351 loopEnd = PAIR_LOOP
352 if (SIM_uses_MetallicAtoms) then
353 loopStart = PREPAIR_LOOP
354 else
355 loopStart = PAIR_LOOP
356 endif
357
358 do loop = loopStart, loopEnd
359
360 ! See if we need to update neighbor lists
361 ! (but only on the first time through):
362 if (loop .eq. loopStart) then
363 #ifdef IS_MPI
364 call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, &
365 update_nlist)
366 #else
367 call checkNeighborList(nGroups, q_group, skinThickness, &
368 update_nlist)
369 #endif
370 endif
371
372 if (update_nlist) then
373 !! save current configuration and construct neighbor list
374 #ifdef IS_MPI
375 call saveNeighborList(nGroupsInRow, q_group_row)
376 #else
377 call saveNeighborList(nGroups, q_group)
378 #endif
379 neighborListSize = size(list)
380 nlist = 0
381 endif
382
383 istart = 1
384 #ifdef IS_MPI
385 iend = nGroupsInRow
386 #else
387 iend = nGroups - 1
388 #endif
389 outer: do i = istart, iend
390
391 if (update_nlist) point(i) = nlist + 1
392
393 n_in_i = groupStartRow(i+1) - groupStartRow(i)
394
395 if (update_nlist) then
396 #ifdef IS_MPI
397 jstart = 1
398 jend = nGroupsInCol
399 #else
400 jstart = i+1
401 jend = nGroups
402 #endif
403 else
404 jstart = point(i)
405 jend = point(i+1) - 1
406 ! make sure group i has neighbors
407 if (jstart .gt. jend) cycle outer
408 endif
409
410 do jnab = jstart, jend
411 if (update_nlist) then
412 j = jnab
413 else
414 j = list(jnab)
415 endif
416
417 #ifdef IS_MPI
418 call get_interatomic_vector(q_group_Row(:,i), &
419 q_group_Col(:,j), d_grp, rgrpsq)
420 #else
421 call get_interatomic_vector(q_group(:,i), &
422 q_group(:,j), d_grp, rgrpsq)
423 #endif
424
425 if (rgrpsq < rListsq) then
426 if (update_nlist) then
427 nlist = nlist + 1
428
429 if (nlist > neighborListSize) then
430 #ifdef IS_MPI
431 call expandNeighborList(nGroupsInRow, listerror)
432 #else
433 call expandNeighborList(nGroups, listerror)
434 #endif
435 if (listerror /= 0) then
436 error = -1
437 write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
438 return
439 end if
440 neighborListSize = size(list)
441 endif
442
443 list(nlist) = j
444 endif
445
446 if (rgrpsq < rCutsq) then
447
448 if (loop .eq. PAIR_LOOP) then
449 vij = 0.0_dp
450 fij(1) = 0.0_dp
451 fij(2) = 0.0_dp
452 fij(3) = 0.0_dp
453 endif
454
455 call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region)
456
457 n_in_j = groupStartCol(j+1) - groupStartCol(j)
458
459 do ia = groupStartRow(i), groupStartRow(i+1)-1
460
461 atom1 = groupListRow(ia)
462
463 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
464
465 atom2 = groupListCol(jb)
466
467 if (skipThisPair(atom1, atom2)) cycle inner
468
469 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
470 d_atm(1) = d_grp(1)
471 d_atm(2) = d_grp(2)
472 d_atm(3) = d_grp(3)
473 ratmsq = rgrpsq
474 else
475 #ifdef IS_MPI
476 call get_interatomic_vector(q_Row(:,atom1), &
477 q_Col(:,atom2), d_atm, ratmsq)
478 #else
479 call get_interatomic_vector(q(:,atom1), &
480 q(:,atom2), d_atm, ratmsq)
481 #endif
482 endif
483
484 topoDist = getTopoDistance(atom1, atom2)
485
486 if (loop .eq. PREPAIR_LOOP) then
487 #ifdef IS_MPI
488 call f_do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
489 rgrpsq, d_grp, rCut, &
490 eFrame, A, f, t, pot_local)
491 #else
492 call f_do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
493 rgrpsq, d_grp, rCut, &
494 eFrame, A, f, t, pot)
495 #endif
496 else
497 #ifdef IS_MPI
498 call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, &
499 eFrame, A, f, t, pot_local, particle_pot, vpair, &
500 fpair, d_grp, rgrp, rCut, topoDist)
501 ! particle_pot will be accumulated from row & column
502 ! arrays later
503 #else
504 call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, &
505 eFrame, A, f, t, pot, particle_pot, vpair, &
506 fpair, d_grp, rgrp, rCut, topoDist)
507 #endif
508 vij = vij + vpair
509 fij(1) = fij(1) + fpair(1)
510 fij(2) = fij(2) + fpair(2)
511 fij(3) = fij(3) + fpair(3)
512 call add_stress_tensor(d_atm, fpair, tau)
513 endif
514 enddo inner
515 enddo
516
517 if (loop .eq. PAIR_LOOP) then
518 if (in_switching_region) then
519 swderiv = vij*dswdr/rgrp
520 fg = swderiv*d_grp
521
522 fij(1) = fij(1) + fg(1)
523 fij(2) = fij(2) + fg(2)
524 fij(3) = fij(3) + fg(3)
525
526 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
527 call add_stress_tensor(d_atm, fg, tau)
528 endif
529
530 do ia=groupStartRow(i), groupStartRow(i+1)-1
531 atom1=groupListRow(ia)
532 mf = mfactRow(atom1)
533 ! fg is the force on atom ia due to cutoff group's
534 ! presence in switching region
535 fg = swderiv*d_grp*mf
536 #ifdef IS_MPI
537 f_Row(1,atom1) = f_Row(1,atom1) + fg(1)
538 f_Row(2,atom1) = f_Row(2,atom1) + fg(2)
539 f_Row(3,atom1) = f_Row(3,atom1) + fg(3)
540 #else
541 f(1,atom1) = f(1,atom1) + fg(1)
542 f(2,atom1) = f(2,atom1) + fg(2)
543 f(3,atom1) = f(3,atom1) + fg(3)
544 #endif
545 if (n_in_i .gt. 1) then
546 if (SIM_uses_AtomicVirial) then
547 ! find the distance between the atom
548 ! and the center of the cutoff group:
549 #ifdef IS_MPI
550 call get_interatomic_vector(q_Row(:,atom1), &
551 q_group_Row(:,i), dag, rag)
552 #else
553 call get_interatomic_vector(q(:,atom1), &
554 q_group(:,i), dag, rag)
555 #endif
556 call add_stress_tensor(dag,fg,tau)
557 endif
558 endif
559 enddo
560
561 do jb=groupStartCol(j), groupStartCol(j+1)-1
562 atom2=groupListCol(jb)
563 mf = mfactCol(atom2)
564 ! fg is the force on atom jb due to cutoff group's
565 ! presence in switching region
566 fg = -swderiv*d_grp*mf
567 #ifdef IS_MPI
568 f_Col(1,atom2) = f_Col(1,atom2) + fg(1)
569 f_Col(2,atom2) = f_Col(2,atom2) + fg(2)
570 f_Col(3,atom2) = f_Col(3,atom2) + fg(3)
571 #else
572 f(1,atom2) = f(1,atom2) + fg(1)
573 f(2,atom2) = f(2,atom2) + fg(2)
574 f(3,atom2) = f(3,atom2) + fg(3)
575 #endif
576 if (n_in_j .gt. 1) then
577 if (SIM_uses_AtomicVirial) then
578 ! find the distance between the atom
579 ! and the center of the cutoff group:
580 #ifdef IS_MPI
581 call get_interatomic_vector(q_Col(:,atom2), &
582 q_group_Col(:,j), dag, rag)
583 #else
584 call get_interatomic_vector(q(:,atom2), &
585 q_group(:,j), dag, rag)
586 #endif
587 call add_stress_tensor(dag,fg,tau)
588 endif
589 endif
590 enddo
591 endif
592 !if (.not.SIM_uses_AtomicVirial) then
593 ! call add_stress_tensor(d_grp, fij, tau)
594 !endif
595 endif
596 endif
597 endif
598 enddo
599
600 enddo outer
601
602 if (update_nlist) then
603 #ifdef IS_MPI
604 point(nGroupsInRow + 1) = nlist + 1
605 #else
606 point(nGroups) = nlist + 1
607 #endif
608 if (loop .eq. PREPAIR_LOOP) then
609 ! we just did the neighbor list update on the first
610 ! pass, so we don't need to do it
611 ! again on the second pass
612 update_nlist = .false.
613 endif
614 endif
615
616 if (loop .eq. PREPAIR_LOOP) then
617 #ifdef IS_MPI
618 call f_do_preforce(nlocal, pot_local, particle_pot)
619 #else
620 call f_do_preforce(nlocal, pot, particle_pot)
621 #endif
622 endif
623
624 enddo
625
626 !! Do timing
627 #ifdef PROFILE
628 call cpu_time(forceTimeFinal)
629 forceTime = forceTime + forceTimeFinal - forceTimeInitial
630 #endif
631
632 #ifdef IS_MPI
633 !!distribute forces
634
635 f_temp = 0.0_dp
636 call scatter(f_Row,f_temp,plan_atom_row_3d)
637 do i = 1,nlocal
638 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
639 end do
640
641 f_temp = 0.0_dp
642 call scatter(f_Col,f_temp,plan_atom_col_3d)
643 do i = 1,nlocal
644 f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
645 end do
646
647 if (SIM_uses_DirectionalAtoms) then
648 t_temp = 0.0_dp
649 call scatter(t_Row,t_temp,plan_atom_row_3d)
650 do i = 1,nlocal
651 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
652 end do
653 t_temp = 0.0_dp
654 call scatter(t_Col,t_temp,plan_atom_col_3d)
655
656 do i = 1,nlocal
657 t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
658 end do
659 endif
660
661 ! scatter/gather pot_row into the members of my column
662 do i = 1,LR_POT_TYPES
663 call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row)
664 end do
665 ! scatter/gather pot_local into all other procs
666 ! add resultant to get total pot
667 do i = 1, nlocal
668 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) &
669 + pot_Temp(1:LR_POT_TYPES,i)
670 enddo
671
672 do i = 1,LR_POT_TYPES
673 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
674 enddo
675
676 pot_Temp = 0.0_DP
677
678 do i = 1,LR_POT_TYPES
679 call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col)
680 end do
681
682 do i = 1, nlocal
683 pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)&
684 + pot_Temp(1:LR_POT_TYPES,i)
685 enddo
686
687 do i = 1,LR_POT_TYPES
688 particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal)
689 enddo
690
691 ppot_Temp = 0.0_DP
692
693 call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row)
694 do i = 1, nlocal
695 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
696 enddo
697
698 ppot_Temp = 0.0_DP
699
700 call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col)
701 do i = 1, nlocal
702 particle_pot(i) = particle_pot(i) + ppot_Temp(i)
703 enddo
704
705 #endif
706
707 if (SIM_requires_skip_correction) then
708 do i = 1, nlocal
709
710 do i1 = 1, nSkipsForLocalAtom(i)
711 j = skipsForLocalAtom(i, i1)
712
713 ! prevent overcounting the skips
714 if (i.lt.j) then
715
716 call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq)
717 rVal = sqrt(ratmsq)
718 call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region)
719 #ifdef IS_MPI
720 call do_skip_correction(c_idents_local(i), c_idents_local(j), &
721 d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, &
722 1.0_dp, pot_local(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j))
723 # else
724 call do_skip_correction(c_idents_local(i), c_idents_local(j), &
725 d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, &
726 1.0_dp, pot(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j))
727 #endif
728 endif
729 enddo
730 enddo
731 endif
732
733 if (SIM_requires_self_correction) then
734
735 do i = 1, nlocal
736
737 #ifdef IS_MPI
738 call do_self_correction(c_idents_local(i), eFrame(:,i), &
739 skipped_charge(i), pot_local(ELECTROSTATIC_POT), t(:,i))
740 #else
741 call do_self_correction(c_idents_local(i), eFrame(:,i), &
742 skipped_charge(i), pot(ELECTROSTATIC_POT), t(:,i))
743 #endif
744 enddo
745 endif
746
747 #ifdef IS_MPI
748 #ifdef SINGLE_PRECISION
749 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, &
750 mpi_comm_world,mpi_err)
751 #else
752 call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, &
753 mpi_sum, mpi_comm_world,mpi_err)
754 #endif
755 #endif
756
757 end subroutine do_force_loop
758
759 subroutine f_do_pair(i, j, rijsq, d, sw, &
760 eFrame, A, f, t, pot, particle_pot, vpair, &
761 fpair, d_grp, r_grp, rCut, topoDist)
762
763 real( kind = dp ) :: vpair, sw
764 real( kind = dp ), dimension(LR_POT_TYPES) :: pot, pairpot
765 real( kind = dp ), dimension(nLocal) :: particle_pot
766 real( kind = dp ), dimension(3) :: fpair
767 real( kind = dp ), dimension(nLocal) :: mfact
768 real( kind = dp ), dimension(9,nLocal) :: eFrame
769 real( kind = dp ), dimension(9,nLocal) :: A
770 real( kind = dp ), dimension(3,nLocal) :: f
771 real( kind = dp ), dimension(3,nLocal) :: t
772
773 integer, intent(in) :: i, j
774 real ( kind = dp ), intent(inout) :: rijsq
775 real ( kind = dp ), intent(inout) :: r_grp
776 real ( kind = dp ), intent(inout) :: d(3)
777 real ( kind = dp ), intent(inout) :: d_grp(3)
778 real ( kind = dp ), intent(inout) :: rCut
779 integer, intent(inout) :: topoDist
780 real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult
781 real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx
782
783 real( kind = dp), dimension(3) :: f1, t1, t2
784 real( kind = dp), dimension(9) :: A1, A2, eF1, eF2
785 real( kind = dp) :: dfrhodrho_i, dfrhodrho_j
786 real( kind = dp) :: rho_i, rho_j
787 real( kind = dp) :: fshift_i, fshift_j
788 real( kind = dp) :: p_vdw, p_elect, p_hb, p_met
789 integer :: id1, id2, idx
790 integer :: k
791 integer :: c_ident_i, c_ident_j
792
793 integer :: iHash
794
795 r = sqrt(rijsq)
796
797 vpair = 0.0_dp
798 fpair(1:3) = 0.0_dp
799
800 p_vdw = 0.0
801 p_elect = 0.0
802 p_hb = 0.0
803 p_met = 0.0
804
805 f1(1:3) = 0.0
806 t1(1:3) = 0.0
807 t2(1:3) = 0.0
808
809 #ifdef IS_MPI
810 c_ident_i = c_idents_row(i)
811 c_ident_j = c_idents_col(j)
812
813 A1(:) = A_Row(:,i)
814 A2(:) = A_Col(:,j)
815 eF1(:) = eFrame_Row(:,i)
816 eF2(:) = eFrame_Col(:,j)
817
818 dfrhodrho_i = dfrhodrho_row(i)
819 dfrhodrho_j = dfrhodrho_col(j)
820 rho_i = rho_row(i)
821 rho_j = rho_col(j)
822 #else
823 c_ident_i = c_idents_local(i)
824 c_ident_j = c_idents_local(j)
825
826 A1(:) = A(:,i)
827 A2(:) = A(:,j)
828 eF1(:) = eFrame(:,i)
829 eF2(:) = eFrame(:,j)
830
831 dfrhodrho_i = dfrhodrho(i)
832 dfrhodrho_j = dfrhodrho(j)
833 rho_i = rho(i)
834 rho_j = rho(j)
835 #endif
836
837 vdwMult = vdwScale(topoDist)
838 electroMult = electrostaticScale(topoDist)
839
840 call do_pair(c_ident_i, c_ident_j, d, r, rijsq, sw, vpair, &
841 vdwMult, electroMult, A1, A2, eF1, eF2, &
842 pairpot, f1, t1, t2, &
843 rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j)
844
845 #ifdef IS_MPI
846 id1 = AtomRowToGlobal(i)
847 id2 = AtomColToGlobal(j)
848
849 pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*pairpot(VDW_POT)
850 pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*pairpot(VDW_POT)
851 pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*pairpot(ELECTROSTATIC_POT)
852 pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*pairpot(ELECTROSTATIC_POT)
853 pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*pairpot(HB_POT)
854 pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*pairpot(HB_POT)
855 pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*pairpot(METALLIC_POT)
856 pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*pairpot(METALLIC_POT)
857
858 do idx = 1, 3
859 f_Row(idx,i) = f_Row(idx,i) + f1(idx)
860 f_Col(idx,j) = f_Col(idx,j) - f1(idx)
861
862 t_Row(idx,i) = t_Row(idx,i) + t1(idx)
863 t_Col(idx,j) = t_Col(idx,j) + t2(idx)
864 enddo
865 ! particle_pot is the difference between the full potential
866 ! and the full potential without the presence of a particular
867 ! particle (atom1).
868 !
869 ! This reduces the density at other particle locations, so
870 ! we need to recompute the density at atom2 assuming atom1
871 ! didn't contribute. This then requires recomputing the
872 ! density functional for atom2 as well.
873 !
874 ! Most of the particle_pot heavy lifting comes from the
875 ! pair interaction, and will be handled by vpair. Parallel version.
876
877 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
878 ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j
879 ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i
880 end if
881
882 #else
883 id1 = i
884 id2 = j
885
886 pot(VDW_POT) = pot(VDW_POT) + pairpot(VDW_POT)
887 pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + pairpot(ELECTROSTATIC_POT)
888 pot(HB_POT) = pot(HB_POT) + pairpot(HB_POT)
889 pot(METALLIC_POT) = pot(METALLIC_POT) + pairpot(METALLIC_POT)
890
891 do idx = 1, 3
892 f(idx,i) = f(idx,i) + f1(idx)
893 f(idx,j) = f(idx,j) - f1(idx)
894
895 t(idx,i) = t(idx,i) + t1(idx)
896 t(idx,j) = t(idx,j) + t2(idx)
897 enddo
898 ! particle_pot is the difference between the full potential
899 ! and the full potential without the presence of a particular
900 ! particle (atom1).
901 !
902 ! This reduces the density at other particle locations, so
903 ! we need to recompute the density at atom2 assuming atom1
904 ! didn't contribute. This then requires recomputing the
905 ! density functional for atom2 as well.
906 !
907 ! Most of the particle_pot heavy lifting comes from the
908 ! pair interaction, and will be handled by vpair. NonParallel version.
909
910 if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then
911 particle_pot(i) = particle_pot(i) - frho(j) + fshift_j
912 particle_pot(j) = particle_pot(j) - frho(i) + fshift_i
913 end if
914
915
916 #endif
917
918 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
919
920 fpair(1) = fpair(1) + f1(1)
921 fpair(2) = fpair(2) + f1(2)
922 fpair(3) = fpair(3) + f1(3)
923
924 endif
925 end subroutine f_do_pair
926
927 subroutine f_do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, &
928 eFrame, A, f, t, pot)
929
930 real( kind = dp ) :: sw
931 real( kind = dp ), dimension(LR_POT_TYPES) :: pot
932 real( kind = dp ), dimension(9,nLocal) :: eFrame
933 real (kind=dp), dimension(9,nLocal) :: A
934 real (kind=dp), dimension(3,nLocal) :: f
935 real (kind=dp), dimension(3,nLocal) :: t
936
937 integer, intent(in) :: i, j
938 real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut
939 real ( kind = dp ) :: r, rc
940 real ( kind = dp ), intent(inout) :: d(3), dc(3)
941 real ( kind = dp ) :: rho_i_at_j, rho_j_at_i
942 integer :: c_ident_i, c_ident_j
943
944 r = sqrt(rijsq)
945
946 #ifdef IS_MPI
947 c_ident_i = c_idents_row(i)
948 c_ident_j = c_idents_col(j)
949 #else
950 c_ident_i = c_idents_local(i)
951 c_ident_j = c_idents_local(j)
952 #endif
953 rho_i_at_j = 0.0_dp
954 rho_j_at_i = 0.0_dp
955
956 call do_prepair(c_ident_i, c_ident_j, r, &
957 rho_i_at_j, rho_j_at_i)
958
959 #ifdef IS_MPI
960 rho_col(j) = rho_col(j) + rho_i_at_j
961 rho_row(i) = rho_row(i) + rho_j_at_i
962 #else
963 rho(j) = rho(j) + rho_i_at_j
964 rho(i) = rho(i) + rho_j_at_i
965 #endif
966
967 end subroutine f_do_prepair
968
969
970 subroutine f_do_preforce(nlocal, pot, particle_pot)
971 integer :: nlocal
972 real( kind = dp ),dimension(LR_POT_TYPES) :: pot
973 real( kind = dp ),dimension(nlocal) :: particle_pot
974 integer :: sc_err = 0
975 integer :: atom, c_ident1
976
977 if (SIM_uses_MetallicAtoms) then
978
979 #ifdef IS_MPI
980 call scatter(rho_row,rho,plan_atom_row,sc_err)
981 if (sc_err /= 0 ) then
982 call handleError("do_preforce()", "Error scattering rho_row into rho")
983 endif
984 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
985 if (sc_err /= 0 ) then
986 call handleError("do_preforce()", "Error scattering rho_col into rho")
987 endif
988 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
989 #endif
990
991
992 do atom = 1, nlocal
993 c_ident1 = c_idents_local(atom)
994
995 call do_preforce(c_ident1, rho(atom), frho(atom), dfrhodrho(atom))
996 pot(METALLIC_POT) = pot(METALLIC_POT) + frho(atom)
997 particle_pot(atom) = particle_pot(atom) + frho(atom)
998 end do
999
1000 #ifdef IS_MPI
1001 !! communicate f(rho) and derivatives back into row and column arrays
1002 call gather(frho,frho_row,plan_atom_row, sc_err)
1003 if (sc_err /= 0) then
1004 call handleError("do_preforce()","MPI gather frho_row failure")
1005 endif
1006 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
1007 if (sc_err /= 0) then
1008 call handleError("do_preforce()","MPI gather dfrhodrho_row failure")
1009 endif
1010 call gather(frho,frho_col,plan_atom_col, sc_err)
1011 if (sc_err /= 0) then
1012 call handleError("do_preforce()","MPI gather frho_col failure")
1013 endif
1014 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
1015 if (sc_err /= 0) then
1016 call handleError("do_preforce()","MPI gather dfrhodrho_col failure")
1017 endif
1018 #endif
1019
1020 end if
1021 end subroutine f_do_preforce
1022
1023
1024 subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1025
1026 real(kind = dp), dimension(3) :: q_i
1027 real(kind = dp), dimension(3) :: q_j
1028 real(kind = dp), intent(out) :: r_sq
1029 real(kind = dp) :: d(3), scaled(3)
1030 real(kind = dp) :: t
1031 integer i
1032
1033 d(1) = q_j(1) - q_i(1)
1034 d(2) = q_j(2) - q_i(2)
1035 d(3) = q_j(3) - q_i(3)
1036
1037 ! Wrap back into periodic box if necessary
1038 if ( SIM_uses_PBC ) then
1039
1040 if( .not.boxIsOrthorhombic ) then
1041 ! calc the scaled coordinates.
1042 ! scaled = matmul(HmatInv, d)
1043 ! unwrap the matmul and do things explicitly:
1044
1045 scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3)
1046 scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3)
1047 scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3)
1048
1049 ! wrap the scaled coordinates (but don't use anint for speed)
1050
1051 t = scaled(1)
1052 if (t .ge. 0.0) then
1053 scaled(1) = t - floor(t + 0.5)
1054 else
1055 scaled(1) = t + ceiling(t - 0.5)
1056 endif
1057
1058 t = scaled(2)
1059 if (t .ge. 0.0) then
1060 scaled(2) = t - floor(t + 0.5)
1061 else
1062 scaled(2) = t + ceiling(t - 0.5)
1063 endif
1064
1065 t = scaled(3)
1066 if (t .ge. 0.0) then
1067 scaled(3) = t - floor(t + 0.5)
1068 else
1069 scaled(3) = t + ceiling(t - 0.5)
1070 endif
1071
1072 ! calc the wrapped real coordinates from the wrapped scaled
1073 ! coordinates
1074 ! d = matmul(Hmat,scaled)
1075 ! unwrap the matmul and do things explicitly:
1076
1077 d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3)
1078 d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3)
1079 d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3)
1080
1081 else
1082 ! calc the scaled coordinates.
1083
1084 scaled(1) = d(1) * HmatInv(1,1)
1085 scaled(2) = d(2) * HmatInv(2,2)
1086 scaled(3) = d(3) * HmatInv(3,3)
1087
1088 ! wrap the scaled coordinates (but don't use anint for speed)
1089
1090 t = scaled(1)
1091 if (t .ge. 0.0) then
1092 scaled(1) = t - floor(t + 0.5)
1093 else
1094 scaled(1) = t + ceiling(t - 0.5)
1095 endif
1096
1097 t = scaled(2)
1098 if (t .ge. 0.0) then
1099 scaled(2) = t - floor(t + 0.5)
1100 else
1101 scaled(2) = t + ceiling(t - 0.5)
1102 endif
1103
1104 t = scaled(3)
1105 if (t .ge. 0.0) then
1106 scaled(3) = t - floor(t + 0.5)
1107 else
1108 scaled(3) = t + ceiling(t - 0.5)
1109 endif
1110
1111 ! calc the wrapped real coordinates from the wrapped scaled
1112 ! coordinates
1113
1114 d(1) = scaled(1)*Hmat(1,1)
1115 d(2) = scaled(2)*Hmat(2,2)
1116 d(3) = scaled(3)*Hmat(3,3)
1117
1118 endif
1119
1120 endif
1121
1122 r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3)
1123
1124 end subroutine get_interatomic_vector
1125
1126 subroutine zero_work_arrays()
1127
1128 #ifdef IS_MPI
1129
1130 q_Row = 0.0_dp
1131 q_Col = 0.0_dp
1132
1133 q_group_Row = 0.0_dp
1134 q_group_Col = 0.0_dp
1135
1136 eFrame_Row = 0.0_dp
1137 eFrame_Col = 0.0_dp
1138
1139 A_Row = 0.0_dp
1140 A_Col = 0.0_dp
1141
1142 f_Row = 0.0_dp
1143 f_Col = 0.0_dp
1144 f_Temp = 0.0_dp
1145
1146 t_Row = 0.0_dp
1147 t_Col = 0.0_dp
1148 t_Temp = 0.0_dp
1149
1150 pot_Row = 0.0_dp
1151 pot_Col = 0.0_dp
1152 pot_Temp = 0.0_dp
1153 ppot_Temp = 0.0_dp
1154
1155 frho_row = 0.0_dp
1156 frho_col = 0.0_dp
1157 rho_row = 0.0_dp
1158 rho_col = 0.0_dp
1159 rho_tmp = 0.0_dp
1160 dfrhodrho_row = 0.0_dp
1161 dfrhodrho_col = 0.0_dp
1162
1163 #endif
1164 rho = 0.0_dp
1165 frho = 0.0_dp
1166 dfrhodrho = 0.0_dp
1167
1168 end subroutine zero_work_arrays
1169
1170 function skipThisPair(atom1, atom2) result(skip_it)
1171 integer, intent(in) :: atom1
1172 integer, intent(in), optional :: atom2
1173 logical :: skip_it
1174 integer :: unique_id_1, unique_id_2
1175 integer :: i
1176
1177 skip_it = .false.
1178
1179 !! there are a number of reasons to skip a pair or a particle
1180 !! mostly we do this to exclude atoms who are involved in short
1181 !! range interactions (bonds, bends, torsions), but we also need
1182 !! to exclude some overcounted interactions that result from
1183 !! the parallel decomposition
1184
1185 #ifdef IS_MPI
1186 !! in MPI, we have to look up the unique IDs for each atom
1187 unique_id_1 = AtomRowToGlobal(atom1)
1188 unique_id_2 = AtomColToGlobal(atom2)
1189 !! this situation should only arise in MPI simulations
1190 if (unique_id_1 == unique_id_2) then
1191 skip_it = .true.
1192 return
1193 end if
1194
1195 !! this prevents us from doing the pair on multiple processors
1196 if (unique_id_1 < unique_id_2) then
1197 if (mod(unique_id_1 + unique_id_2,2) == 0) then
1198 skip_it = .true.
1199 return
1200 endif
1201 else
1202 if (mod(unique_id_1 + unique_id_2,2) == 1) then
1203 skip_it = .true.
1204 return
1205 endif
1206 endif
1207 #else
1208 !! in the normal loop, the atom numbers are unique
1209 unique_id_1 = atom1
1210 unique_id_2 = atom2
1211 #endif
1212
1213 #ifdef IS_MPI
1214 do i = 1, nSkipsForRowAtom(atom1)
1215 if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then
1216 skip_it = .true.
1217 return
1218 endif
1219 end do
1220 #else
1221 do i = 1, nSkipsForLocalAtom(atom1)
1222 if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then
1223 skip_it = .true.
1224 return
1225 endif
1226 end do
1227 #endif
1228
1229 return
1230 end function skipThisPair
1231
1232 function getTopoDistance(atom1, atom2) result(topoDist)
1233 integer, intent(in) :: atom1
1234 integer, intent(in) :: atom2
1235 integer :: topoDist
1236 integer :: unique_id_2
1237 integer :: i
1238
1239 #ifdef IS_MPI
1240 unique_id_2 = AtomColToGlobal(atom2)
1241 #else
1242 unique_id_2 = atom2
1243 #endif
1244
1245 ! zero is default for unconnected (i.e. normal) pair interactions
1246
1247 topoDist = 0
1248
1249 do i = 1, nTopoPairsForAtom(atom1)
1250 if (toposForAtom(atom1, i) .eq. unique_id_2) then
1251 topoDist = topoDistance(atom1, i)
1252 return
1253 endif
1254 end do
1255
1256 return
1257 end function getTopoDistance
1258
1259 #ifdef PROFILE
1260 function getforcetime() result(totalforcetime)
1261 real(kind=dp) :: totalforcetime
1262 totalforcetime = forcetime
1263 end function getforcetime
1264 #endif
1265
1266 !! This cleans componets of force arrays belonging only to fortran
1267
1268 subroutine add_stress_tensor(dpair, fpair, tau)
1269
1270 real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1271 real( kind = dp ), dimension(9), intent(inout) :: tau
1272
1273 ! because the d vector is the rj - ri vector, and
1274 ! because fx, fy, fz are the force on atom i, we need a
1275 ! negative sign here:
1276
1277 tau(1) = tau(1) - dpair(1) * fpair(1)
1278 tau(2) = tau(2) - dpair(1) * fpair(2)
1279 tau(3) = tau(3) - dpair(1) * fpair(3)
1280 tau(4) = tau(4) - dpair(2) * fpair(1)
1281 tau(5) = tau(5) - dpair(2) * fpair(2)
1282 tau(6) = tau(6) - dpair(2) * fpair(3)
1283 tau(7) = tau(7) - dpair(3) * fpair(1)
1284 tau(8) = tau(8) - dpair(3) * fpair(2)
1285 tau(9) = tau(9) - dpair(3) * fpair(3)
1286
1287 end subroutine add_stress_tensor
1288
1289 end module doForces

Properties

Name Value
svn:keywords Author Id Revision Date