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: 1683
Committed: Thu Oct 28 22:34:02 2004 UTC (19 years, 8 months ago)
File size: 6775 byte(s)
Log Message:
This commit was manufactured by cvs2svn to create branch 'new_design'.

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