1 |
gezelter |
940 |
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 |
gezelter |
945 |
real(kind=dp), save :: ecr = 0.0_DP |
16 |
gezelter |
940 |
real(kind=dp), save :: rt = 0.0_DP |
17 |
tim |
1132 |
real(kind=dp), save :: dielect = 0.0_DP |
18 |
gezelter |
940 |
real(kind=dp), save :: pre = 0.0_DP |
19 |
tim |
1132 |
|
20 |
gezelter |
940 |
logical, save :: haveCutoffs = .false. |
21 |
|
|
logical, save :: haveChargeMap = .false. |
22 |
tim |
1132 |
logical, save :: haveDielectric = .false. |
23 |
gezelter |
940 |
|
24 |
gezelter |
945 |
public::setCutoffsCharge |
25 |
gezelter |
940 |
public::do_charge_pair |
26 |
tim |
1132 |
public::initialize_charge |
27 |
gezelter |
940 |
|
28 |
|
|
type :: ChargeList |
29 |
|
|
real(kind=DP) :: charge = 0.0_DP |
30 |
|
|
end type ChargeList |
31 |
|
|
|
32 |
|
|
type(ChargeList), dimension(:), allocatable :: ChargeMap |
33 |
|
|
|
34 |
|
|
contains |
35 |
tim |
1132 |
|
36 |
|
|
subroutine initialize_charge(this_dielect) |
37 |
|
|
real(kind=dp), intent(in) :: this_dielect |
38 |
gezelter |
940 |
|
39 |
tim |
1132 |
dielect = this_dielect |
40 |
|
|
haveDielectric = .true. |
41 |
|
|
|
42 |
|
|
! because setCutoffsCharge is called before initialize_charge |
43 |
|
|
! we need to call it agin to make sure all of the precalculation |
44 |
|
|
! value is correct. for the time being, just a quick hack |
45 |
|
|
call setCutoffsCharge(ecr, rt) |
46 |
|
|
return |
47 |
|
|
end subroutine initialize_charge |
48 |
|
|
|
49 |
|
|
|
50 |
gezelter |
945 |
subroutine setCutoffsCharge(this_ecr, this_rt) |
51 |
|
|
real(kind=dp), intent(in) :: this_ecr, this_rt |
52 |
|
|
ecr = this_ecr |
53 |
gezelter |
940 |
rt = this_rt |
54 |
|
|
|
55 |
|
|
! pre converts from fundamental charge to kcal/mol |
56 |
|
|
pre = 332.06508_DP |
57 |
tim |
1132 |
|
58 |
|
|
!if (.not.haveDielectric)then |
59 |
|
|
! write(default_error,*) 'Dielect constant in charge module is not set!' |
60 |
|
|
!endif |
61 |
|
|
|
62 |
gezelter |
940 |
haveCutoffs = .true. |
63 |
|
|
|
64 |
|
|
return |
65 |
gezelter |
945 |
end subroutine setCutoffsCharge |
66 |
gezelter |
940 |
|
67 |
|
|
subroutine createChargeMap(status) |
68 |
|
|
integer :: nAtypes |
69 |
|
|
integer :: status |
70 |
|
|
integer :: i |
71 |
|
|
real (kind=DP) :: thisCharge |
72 |
|
|
logical :: thisProperty |
73 |
|
|
|
74 |
|
|
status = 0 |
75 |
|
|
|
76 |
|
|
nAtypes = getSize(atypes) |
77 |
|
|
|
78 |
|
|
if (nAtypes == 0) then |
79 |
|
|
status = -1 |
80 |
|
|
return |
81 |
|
|
end if |
82 |
|
|
|
83 |
|
|
if (.not. allocated(ChargeMap)) then |
84 |
|
|
allocate(ChargeMap(nAtypes)) |
85 |
|
|
endif |
86 |
|
|
|
87 |
|
|
do i = 1, nAtypes |
88 |
|
|
|
89 |
|
|
call getElementProperty(atypes, i, "is_Charge", thisProperty) |
90 |
|
|
|
91 |
|
|
if (thisProperty) then |
92 |
|
|
call getElementProperty(atypes, i, "charge", thisCharge) |
93 |
|
|
ChargeMap(i)%charge = thisCharge |
94 |
|
|
endif |
95 |
|
|
|
96 |
|
|
end do |
97 |
|
|
|
98 |
|
|
haveChargeMap = .true. |
99 |
|
|
|
100 |
|
|
end subroutine createChargeMap |
101 |
|
|
|
102 |
|
|
|
103 |
gezelter |
1150 |
subroutine do_charge_pair(atom1, atom2, d, rij, r2, sw, vpair, & |
104 |
|
|
pot, f, do_pot, do_stress) |
105 |
gezelter |
940 |
|
106 |
gezelter |
1150 |
logical :: do_pot, do_stress |
107 |
gezelter |
940 |
|
108 |
|
|
integer atom1, atom2, me1, me2, id1, id2 |
109 |
|
|
integer :: localError |
110 |
gezelter |
1150 |
real(kind=dp) :: rij, r2, q1, q2, sw, vpair |
111 |
gezelter |
940 |
real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz |
112 |
gezelter |
1150 |
real(kind=dp) :: vterm |
113 |
gezelter |
940 |
|
114 |
|
|
real( kind = dp ) :: pot |
115 |
gezelter |
1150 |
real( kind = dp ), dimension(3) :: d |
116 |
gezelter |
940 |
real( kind = dp ), dimension(3,nLocal) :: f |
117 |
gezelter |
1150 |
|
118 |
gezelter |
940 |
|
119 |
|
|
if (.not. haveCutoffs) then |
120 |
|
|
write(default_error,*) 'charge-charge does not have cutoffs set!' |
121 |
|
|
return |
122 |
|
|
endif |
123 |
|
|
|
124 |
|
|
if (.not.haveChargeMap) then |
125 |
|
|
localError = 0 |
126 |
|
|
call createChargeMap(localError) |
127 |
|
|
if ( localError .ne. 0 ) then |
128 |
|
|
call handleError("charge-charge", "ChargeMap creation failed!") |
129 |
|
|
return |
130 |
|
|
end if |
131 |
|
|
endif |
132 |
|
|
|
133 |
|
|
#ifdef IS_MPI |
134 |
|
|
me1 = atid_Row(atom1) |
135 |
|
|
me2 = atid_Col(atom2) |
136 |
|
|
#else |
137 |
|
|
me1 = atid(atom1) |
138 |
|
|
me2 = atid(atom2) |
139 |
|
|
#endif |
140 |
|
|
|
141 |
|
|
q1 = ChargeMap(me1)%charge |
142 |
|
|
q2 = ChargeMap(me2)%charge |
143 |
gezelter |
1138 |
|
144 |
gezelter |
1150 |
vterm = sw * pre * q1 * q2 / rij |
145 |
gezelter |
1138 |
|
146 |
gezelter |
1150 |
dudr = -vterm / rij |
147 |
|
|
drdx = d(1) / rij |
148 |
|
|
drdy = d(2) / rij |
149 |
|
|
drdz = d(3) / rij |
150 |
|
|
|
151 |
|
|
fx = dudr * drdx |
152 |
|
|
fy = dudr * drdy |
153 |
|
|
fz = dudr * drdz |
154 |
gezelter |
940 |
|
155 |
gezelter |
1150 |
vpair = vpair + vterm |
156 |
|
|
|
157 |
gezelter |
940 |
#ifdef IS_MPI |
158 |
gezelter |
1150 |
if (do_pot) then |
159 |
|
|
pot_Row(atom1) = pot_Row(atom1) + vterm*0.5 |
160 |
|
|
pot_Col(atom2) = pot_Col(atom2) + vterm*0.5 |
161 |
|
|
endif |
162 |
|
|
|
163 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + fx |
164 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + fy |
165 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + fz |
166 |
|
|
|
167 |
|
|
f_Col(1,atom2) = f_Col(1,atom2) - fx |
168 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - fy |
169 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - fz |
170 |
|
|
|
171 |
gezelter |
940 |
#else |
172 |
|
|
|
173 |
gezelter |
1150 |
if (do_pot) pot = pot + vterm |
174 |
|
|
|
175 |
|
|
f(1,atom1) = f(1,atom1) + fx |
176 |
|
|
f(2,atom1) = f(2,atom1) + fy |
177 |
|
|
f(3,atom1) = f(3,atom1) + fz |
178 |
|
|
|
179 |
|
|
f(1,atom2) = f(1,atom2) - fx |
180 |
|
|
f(2,atom2) = f(2,atom2) - fy |
181 |
|
|
f(3,atom2) = f(3,atom2) - fz |
182 |
|
|
|
183 |
gezelter |
940 |
#endif |
184 |
gezelter |
1150 |
|
185 |
|
|
if (do_stress) then |
186 |
|
|
|
187 |
gezelter |
940 |
#ifdef IS_MPI |
188 |
gezelter |
1150 |
id1 = tagRow(atom1) |
189 |
|
|
id2 = tagColumn(atom2) |
190 |
gezelter |
940 |
#else |
191 |
gezelter |
1150 |
id1 = atom1 |
192 |
|
|
id2 = atom2 |
193 |
gezelter |
940 |
#endif |
194 |
gezelter |
1150 |
|
195 |
|
|
|
196 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
197 |
|
|
|
198 |
|
|
! because the d vector is the rj - ri vector, and |
199 |
|
|
! because fx, fy, fz are the force on atom i, we need a |
200 |
|
|
! negative sign here: |
201 |
|
|
|
202 |
|
|
tau_Temp(1) = tau_Temp(1) - d(1) * fx |
203 |
|
|
tau_Temp(2) = tau_Temp(2) - d(1) * fy |
204 |
|
|
tau_Temp(3) = tau_Temp(3) - d(1) * fz |
205 |
|
|
tau_Temp(4) = tau_Temp(4) - d(2) * fx |
206 |
|
|
tau_Temp(5) = tau_Temp(5) - d(2) * fy |
207 |
|
|
tau_Temp(6) = tau_Temp(6) - d(2) * fz |
208 |
|
|
tau_Temp(8) = tau_Temp(8) - d(3) * fy |
209 |
|
|
tau_Temp(9) = tau_Temp(9) - d(3) * fz |
210 |
|
|
|
211 |
|
|
|
212 |
|
|
virial_Temp = virial_Temp + & |
213 |
|
|
(tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) |
214 |
|
|
|
215 |
gezelter |
940 |
endif |
216 |
|
|
endif |
217 |
gezelter |
1150 |
|
218 |
|
|
|
219 |
gezelter |
940 |
return |
220 |
|
|
end subroutine do_charge_pair |
221 |
|
|
|
222 |
|
|
end module charge_charge |