1 |
module charge_charge |
2 |
|
3 |
use force_globals |
4 |
use definitions |
5 |
use atype_module |
6 |
use vector_class |
7 |
use simulation |
8 |
use status |
9 |
#ifdef IS_MPI |
10 |
use mpiSimulation |
11 |
#endif |
12 |
implicit none |
13 |
|
14 |
PRIVATE |
15 |
|
16 |
real(kind=dp), parameter :: pre = 332.06508_DP |
17 |
logical, save :: haveChargeMap = .false. |
18 |
|
19 |
public :: newChargeType |
20 |
public :: do_charge_pair |
21 |
public :: getCharge |
22 |
|
23 |
type :: ChargeList |
24 |
integer :: c_ident |
25 |
real(kind=DP) :: charge = 0.0_DP |
26 |
end type ChargeList |
27 |
|
28 |
type(ChargeList), dimension(:), allocatable :: ChargeMap |
29 |
|
30 |
contains |
31 |
|
32 |
subroutine newChargeType(c_ident, charge, status) |
33 |
integer,intent(in) :: c_ident |
34 |
real(kind=dp),intent(in) :: charge |
35 |
integer,intent(out) :: status |
36 |
integer :: nAtypes, myATID |
37 |
|
38 |
status = 0 |
39 |
|
40 |
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
41 |
|
42 |
!! Be simple-minded and assume that we need a ChargeMap that |
43 |
!! is the same size as the total number of atom types |
44 |
|
45 |
if (.not.allocated(ChargeMap)) then |
46 |
|
47 |
nAtypes = getSize(atypes) |
48 |
|
49 |
if (nAtypes == 0) then |
50 |
status = -1 |
51 |
return |
52 |
end if |
53 |
|
54 |
if (.not. allocated(ChargeMap)) then |
55 |
allocate(ChargeMap(nAtypes)) |
56 |
endif |
57 |
|
58 |
end if |
59 |
|
60 |
if (myATID .gt. size(ChargeMap)) then |
61 |
status = -1 |
62 |
return |
63 |
endif |
64 |
|
65 |
! set the values for ChargeMap for this atom type: |
66 |
|
67 |
ChargeMap(myATID)%c_ident = c_ident |
68 |
ChargeMap(myATID)%charge = charge |
69 |
|
70 |
end subroutine newChargeType |
71 |
|
72 |
function getCharge(atid) result (c) |
73 |
integer, intent(in) :: atid |
74 |
integer :: localError |
75 |
real(kind=dp) :: c |
76 |
|
77 |
if (.not.allocated(ChargeMap)) then |
78 |
call handleError("charge_charge", "no ChargeMap was present before first call of getCharge!") |
79 |
return |
80 |
end if |
81 |
|
82 |
c = ChargeMap(atid)%charge |
83 |
end function getCharge |
84 |
|
85 |
subroutine do_charge_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & |
86 |
pot, f, do_pot) |
87 |
|
88 |
logical :: do_pot |
89 |
|
90 |
integer atom1, atom2, me1, me2, id1, id2 |
91 |
integer :: localError |
92 |
real(kind=dp) :: rij, r2, q1, q2, sw, vpair |
93 |
real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz |
94 |
real(kind=dp) :: vterm |
95 |
|
96 |
real( kind = dp ) :: pot |
97 |
real( kind = dp ), dimension(3) :: d, fpair |
98 |
real( kind = dp ), dimension(3,nLocal) :: f |
99 |
|
100 |
|
101 |
if (.not.allocated(ChargeMap)) then |
102 |
call handleError("charge_charge", "no ChargeMap was present before first call of do_charge_pair!") |
103 |
return |
104 |
end if |
105 |
|
106 |
#ifdef IS_MPI |
107 |
me1 = atid_Row(atom1) |
108 |
me2 = atid_Col(atom2) |
109 |
#else |
110 |
me1 = atid(atom1) |
111 |
me2 = atid(atom2) |
112 |
#endif |
113 |
|
114 |
q1 = ChargeMap(me1)%charge |
115 |
q2 = ChargeMap(me2)%charge |
116 |
|
117 |
vterm = pre * q1 * q2 / rij |
118 |
|
119 |
dudr = -sw * vterm / rij |
120 |
|
121 |
drdx = d(1) / rij |
122 |
drdy = d(2) / rij |
123 |
drdz = d(3) / rij |
124 |
|
125 |
fx = dudr * drdx |
126 |
fy = dudr * drdy |
127 |
fz = dudr * drdz |
128 |
|
129 |
vpair = vpair + vterm |
130 |
|
131 |
#ifdef IS_MPI |
132 |
if (do_pot) then |
133 |
pot_Row(atom1) = pot_Row(atom1) + sw*vterm*0.5 |
134 |
pot_Col(atom2) = pot_Col(atom2) + sw*vterm*0.5 |
135 |
endif |
136 |
|
137 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
138 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
139 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
140 |
|
141 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
142 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
143 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
144 |
|
145 |
#else |
146 |
|
147 |
if (do_pot) pot = pot + sw*vterm |
148 |
|
149 |
f(1,atom1) = f(1,atom1) + fx |
150 |
f(2,atom1) = f(2,atom1) + fy |
151 |
f(3,atom1) = f(3,atom1) + fz |
152 |
|
153 |
f(1,atom2) = f(1,atom2) - fx |
154 |
f(2,atom2) = f(2,atom2) - fy |
155 |
f(3,atom2) = f(3,atom2) - fz |
156 |
|
157 |
#endif |
158 |
|
159 |
#ifdef IS_MPI |
160 |
id1 = AtomRowToGlobal(atom1) |
161 |
id2 = AtomColToGlobal(atom2) |
162 |
#else |
163 |
id1 = atom1 |
164 |
id2 = atom2 |
165 |
#endif |
166 |
|
167 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
168 |
|
169 |
fpair(1) = fpair(1) + fx |
170 |
fpair(2) = fpair(2) + fy |
171 |
fpair(3) = fpair(3) + fz |
172 |
|
173 |
endif |
174 |
return |
175 |
end subroutine do_charge_pair |
176 |
|
177 |
end module charge_charge |
178 |
|
179 |
subroutine newChargeType(ident, charge, status) |
180 |
|
181 |
use charge_charge, ONLY : module_newChargeType => newChargeType |
182 |
|
183 |
integer, parameter :: DP = selected_real_kind(15) |
184 |
integer,intent(inout) :: ident |
185 |
real(kind=dp),intent(inout) :: charge |
186 |
integer,intent(inout) :: status |
187 |
|
188 |
call module_newChargeType(ident, charge, status) |
189 |
|
190 |
end subroutine newChargeType |