ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-4/src/UseTheForce/DarkSide/charge.F90
Revision: 1683
Committed: Thu Oct 28 22:34:02 2004 UTC (19 years, 8 months ago)
File size: 4269 byte(s)
Log Message:
This commit was manufactured by cvs2svn to create branch 'new_design'.

File Contents

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