ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_charge_charge.F90
Revision: 1132
Committed: Sat Apr 24 04:31:36 2004 UTC (20 years, 2 months ago) by tim
File size: 6289 byte(s)
Log Message:
add reaction field correction to charge-charge interaction

File Contents

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