ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/UseTheForce/DarkSide/dipole.F90
Revision: 1859
Committed: Mon Dec 6 23:29:20 2004 UTC (19 years, 6 months ago) by tim
File size: 6861 byte(s)
Log Message:
fix atom type ident in Charge and Dipole Module

File Contents

# User Rev Content
1 gezelter 1608 module dipole_dipole
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 = 14.38362_dp
17    
18 gezelter 1633 public :: newDipoleType
19     public :: do_dipole_pair
20     public :: getDipoleMoment
21 gezelter 1608
22     type :: MomentList
23 tim 1859 integer :: c_ident
24 gezelter 1608 real(kind=DP) :: dipole_moment = 0.0_DP
25     end type MomentList
26    
27     type(MomentList), dimension(:),allocatable :: MomentMap
28    
29     contains
30    
31 tim 1859 subroutine newDipoleType(c_ident, dipole_moment, status)
32     integer,intent(in) :: c_ident
33 gezelter 1633 real(kind=dp),intent(in) :: dipole_moment
34     integer,intent(out) :: status
35 tim 1859 integer :: nAtypes, myATID
36 gezelter 1608
37     status = 0
38 tim 1859 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
39 gezelter 1633
40     !! Be simple-minded and assume that we need a MomentMap that
41     !! is the same size as the total number of atom types
42 gezelter 1608
43 gezelter 1633 if (.not.allocated(MomentMap)) then
44    
45     nAtypes = getSize(atypes)
46 gezelter 1608
47 gezelter 1633 if (nAtypes == 0) then
48     status = -1
49     return
50     end if
51    
52     if (.not. allocated(MomentMap)) then
53     allocate(MomentMap(nAtypes))
54     endif
55    
56     end if
57    
58 tim 1859 if (myATID .gt. size(MomentMap)) then
59 gezelter 1608 status = -1
60     return
61 gezelter 1633 endif
62 gezelter 1608
63 gezelter 1633 ! set the values for MomentMap for this atom type:
64 gezelter 1608
65 tim 1859 MomentMap(myATID)%c_ident = c_ident
66     MomentMap(myATID)%dipole_moment = dipole_moment
67 gezelter 1633
68     end subroutine newDipoleType
69 gezelter 1608
70 gezelter 1633 function getDipoleMoment(atid) result (dm)
71     integer, intent(in) :: atid
72     integer :: localError
73     real(kind=dp) :: dm
74 gezelter 1608
75 gezelter 1633 if (.not.allocated(MomentMap)) then
76     call handleError("dipole-dipole", "no MomentMap was present before first call of getDipoleMoment!")
77     return
78     end if
79 gezelter 1608
80 gezelter 1633 dm = MomentMap(atid)%dipole_moment
81     end function getDipoleMoment
82    
83 gezelter 1608 subroutine do_dipole_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
84     pot, u_l, f, t, do_pot)
85    
86     logical :: do_pot
87    
88     integer atom1, atom2, me1, me2, id1, id2
89     integer :: localError
90     real(kind=dp) :: rij, mu1, mu2
91     real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5
92     real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
93     real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
94     real(kind=dp) :: sw, vpair, vterm
95    
96     real( kind = dp ) :: pot
97     real( kind = dp ), dimension(3) :: d, fpair
98     real( kind = dp ), dimension(3,nLocal) :: u_l
99     real( kind = dp ), dimension(3,nLocal) :: f
100     real( kind = dp ), dimension(3,nLocal) :: t
101    
102     real (kind = dp), dimension(3) :: ul1
103     real (kind = dp), dimension(3) :: ul2
104    
105    
106 gezelter 1633 if (.not.allocated(MomentMap)) then
107     call handleError("dipole-dipole", "no MomentMap was present before first call of do_dipole_pair!")
108     return
109     end if
110    
111    
112 gezelter 1608 #ifdef IS_MPI
113     me1 = atid_Row(atom1)
114     ul1(1) = u_l_Row(1,atom1)
115     ul1(2) = u_l_Row(2,atom1)
116     ul1(3) = u_l_Row(3,atom1)
117    
118     me2 = atid_Col(atom2)
119     ul2(1) = u_l_Col(1,atom2)
120     ul2(2) = u_l_Col(2,atom2)
121     ul2(3) = u_l_Col(3,atom2)
122     #else
123     me1 = atid(atom1)
124     ul1(1) = u_l(1,atom1)
125     ul1(2) = u_l(2,atom1)
126     ul1(3) = u_l(3,atom1)
127    
128     me2 = atid(atom2)
129     ul2(1) = u_l(1,atom2)
130     ul2(2) = u_l(2,atom2)
131     ul2(3) = u_l(3,atom2)
132     #endif
133    
134     mu1 = MomentMap(me1)%dipole_moment
135     mu2 = MomentMap(me2)%dipole_moment
136    
137     r3 = r2*rij
138     r5 = r3*r2
139    
140     rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
141     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
142     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
143    
144     dip2 = pre * mu1 * mu2
145     dfact1 = 3.0d0*dip2 / r2
146     dfact2 = 3.0d0*dip2 / r5
147    
148     vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
149    
150     vpair = vpair + vterm
151    
152     if (do_pot) then
153     #ifdef IS_MPI
154     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*sw
155     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*sw
156     #else
157     pot = pot + vterm*sw
158     #endif
159     endif
160    
161     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
162     (5.0d0*(rdotu1*rdotu2)/r5)) - &
163     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*sw
164    
165     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
166     (5.0d0*(rdotu1*rdotu2)/r5)) - &
167     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*sw
168    
169     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
170     (5.0d0*(rdotu1*rdotu2)/r5)) - &
171     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*sw
172    
173     dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*sw
174     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*sw
175     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*sw
176    
177     dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*sw
178     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*sw
179     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*sw
180    
181    
182     #ifdef IS_MPI
183     f_Row(1,atom1) = f_Row(1,atom1) + dudx
184     f_Row(2,atom1) = f_Row(2,atom1) + dudy
185     f_Row(3,atom1) = f_Row(3,atom1) + dudz
186    
187     f_Col(1,atom2) = f_Col(1,atom2) - dudx
188     f_Col(2,atom2) = f_Col(2,atom2) - dudy
189     f_Col(3,atom2) = f_Col(3,atom2) - dudz
190    
191     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
192     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
193     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
194    
195     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
196     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
197     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
198     #else
199     f(1,atom1) = f(1,atom1) + dudx
200     f(2,atom1) = f(2,atom1) + dudy
201     f(3,atom1) = f(3,atom1) + dudz
202    
203     f(1,atom2) = f(1,atom2) - dudx
204     f(2,atom2) = f(2,atom2) - dudy
205     f(3,atom2) = f(3,atom2) - dudz
206    
207     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
208     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
209     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
210    
211     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
212     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
213     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
214     #endif
215    
216     #ifdef IS_MPI
217     id1 = AtomRowToGlobal(atom1)
218     id2 = AtomColToGlobal(atom2)
219     #else
220     id1 = atom1
221     id2 = atom2
222     #endif
223    
224     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
225    
226     fpair(1) = fpair(1) + dudx
227     fpair(2) = fpair(2) + dudy
228     fpair(3) = fpair(3) + dudz
229    
230     endif
231    
232     return
233     end subroutine do_dipole_pair
234    
235     end module dipole_dipole
236 gezelter 1633
237     subroutine newDipoleType(ident, dipole_moment, status)
238    
239     use dipole_dipole, ONLY : module_newDipoleType => newDipoleType
240    
241     integer, parameter :: DP = selected_real_kind(15)
242     integer,intent(inout) :: ident
243     real(kind=dp),intent(inout) :: dipole_moment
244     integer,intent(inout) :: status
245    
246     call module_newDipoleType(ident, dipole_moment, status)
247    
248     end subroutine newDipoleType