ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/dipole.F90
Revision: 1948
Committed: Fri Jan 14 20:31:16 2005 UTC (19 years, 7 months ago) by gezelter
File size: 8595 byte(s)
Log Message:
separating modules and C/Fortran interface subroutines

File Contents

# User Rev Content
1 gezelter 1930 !!
2     !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     !!
4     !! The University of Notre Dame grants you ("Licensee") a
5     !! non-exclusive, royalty free, license to use, modify and
6     !! redistribute this software in source and binary code form, provided
7     !! that the following conditions are met:
8     !!
9     !! 1. Acknowledgement of the program authors must be made in any
10     !! publication of scientific results based in part on use of the
11     !! program. An acceptable form of acknowledgement is citation of
12     !! the article in which the program was described (Matthew
13     !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42 gezelter 1608 module dipole_dipole
43    
44     use force_globals
45     use definitions
46     use atype_module
47     use vector_class
48     use simulation
49     use status
50     #ifdef IS_MPI
51     use mpiSimulation
52     #endif
53     implicit none
54    
55     PRIVATE
56    
57     real(kind=dp), parameter :: pre = 14.38362_dp
58    
59 gezelter 1633 public :: newDipoleType
60     public :: do_dipole_pair
61     public :: getDipoleMoment
62 gezelter 1608
63     type :: MomentList
64 gezelter 1930 integer :: c_ident
65 gezelter 1608 real(kind=DP) :: dipole_moment = 0.0_DP
66     end type MomentList
67    
68     type(MomentList), dimension(:),allocatable :: MomentMap
69    
70     contains
71    
72 gezelter 1930 subroutine newDipoleType(c_ident, dipole_moment, status)
73     integer,intent(in) :: c_ident
74 gezelter 1633 real(kind=dp),intent(in) :: dipole_moment
75     integer,intent(out) :: status
76 gezelter 1930 integer :: nAtypes, myATID
77 gezelter 1608
78     status = 0
79 gezelter 1930 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
80 gezelter 1633
81     !! Be simple-minded and assume that we need a MomentMap that
82     !! is the same size as the total number of atom types
83 gezelter 1608
84 gezelter 1633 if (.not.allocated(MomentMap)) then
85    
86     nAtypes = getSize(atypes)
87 gezelter 1608
88 gezelter 1633 if (nAtypes == 0) then
89     status = -1
90     return
91     end if
92    
93     if (.not. allocated(MomentMap)) then
94     allocate(MomentMap(nAtypes))
95     endif
96    
97     end if
98    
99 gezelter 1930 if (myATID .gt. size(MomentMap)) then
100 gezelter 1608 status = -1
101     return
102 gezelter 1633 endif
103 gezelter 1608
104 gezelter 1633 ! set the values for MomentMap for this atom type:
105 gezelter 1608
106 gezelter 1930 MomentMap(myATID)%c_ident = c_ident
107     MomentMap(myATID)%dipole_moment = dipole_moment
108 gezelter 1633
109     end subroutine newDipoleType
110 gezelter 1608
111 gezelter 1633 function getDipoleMoment(atid) result (dm)
112     integer, intent(in) :: atid
113     integer :: localError
114     real(kind=dp) :: dm
115 gezelter 1608
116 gezelter 1633 if (.not.allocated(MomentMap)) then
117     call handleError("dipole-dipole", "no MomentMap was present before first call of getDipoleMoment!")
118     return
119     end if
120 gezelter 1608
121 gezelter 1633 dm = MomentMap(atid)%dipole_moment
122     end function getDipoleMoment
123    
124 gezelter 1608 subroutine do_dipole_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
125 gezelter 1930 pot, eFrame, f, t, do_pot)
126 gezelter 1608
127     logical :: do_pot
128    
129     integer atom1, atom2, me1, me2, id1, id2
130     integer :: localError
131     real(kind=dp) :: rij, mu1, mu2
132     real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5
133     real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
134     real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
135     real(kind=dp) :: sw, vpair, vterm
136    
137     real( kind = dp ) :: pot
138     real( kind = dp ), dimension(3) :: d, fpair
139 gezelter 1930 real( kind = dp ), dimension(9,nLocal) :: eFrame
140 gezelter 1608 real( kind = dp ), dimension(3,nLocal) :: f
141     real( kind = dp ), dimension(3,nLocal) :: t
142    
143     real (kind = dp), dimension(3) :: ul1
144     real (kind = dp), dimension(3) :: ul2
145    
146    
147 gezelter 1633 if (.not.allocated(MomentMap)) then
148     call handleError("dipole-dipole", "no MomentMap was present before first call of do_dipole_pair!")
149     return
150     end if
151    
152    
153 gezelter 1608 #ifdef IS_MPI
154     me1 = atid_Row(atom1)
155 gezelter 1930 ul1(1) = eFrame_Row(3,atom1)
156     ul1(2) = eFrame_Row(6,atom1)
157     ul1(3) = eFrame_Row(9,atom1)
158 gezelter 1608
159     me2 = atid_Col(atom2)
160 gezelter 1930 ul2(1) = eFrame_Col(3,atom2)
161     ul2(2) = eFrame_Col(6,atom2)
162     ul2(3) = eFrame_Col(9,atom2)
163 gezelter 1608 #else
164     me1 = atid(atom1)
165 gezelter 1930 ul1(1) = eFrame(3,atom1)
166     ul1(2) = eFrame(6,atom1)
167     ul1(3) = eFrame(9,atom1)
168 gezelter 1608
169     me2 = atid(atom2)
170 gezelter 1930 ul2(1) = eFrame(3,atom2)
171     ul2(2) = eFrame(6,atom2)
172     ul2(3) = eFrame(9,atom2)
173 gezelter 1608 #endif
174    
175     mu1 = MomentMap(me1)%dipole_moment
176     mu2 = MomentMap(me2)%dipole_moment
177    
178     r3 = r2*rij
179     r5 = r3*r2
180    
181     rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
182     rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
183     u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
184    
185     dip2 = pre * mu1 * mu2
186     dfact1 = 3.0d0*dip2 / r2
187     dfact2 = 3.0d0*dip2 / r5
188    
189     vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
190    
191     vpair = vpair + vterm
192    
193     if (do_pot) then
194     #ifdef IS_MPI
195     pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*sw
196     pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*sw
197     #else
198     pot = pot + vterm*sw
199     #endif
200     endif
201    
202     dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
203     (5.0d0*(rdotu1*rdotu2)/r5)) - &
204     dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*sw
205    
206     dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
207     (5.0d0*(rdotu1*rdotu2)/r5)) - &
208     dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*sw
209    
210     dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
211     (5.0d0*(rdotu1*rdotu2)/r5)) - &
212     dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*sw
213    
214     dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*sw
215     dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*sw
216     dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*sw
217    
218     dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*sw
219     dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*sw
220     dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*sw
221    
222    
223     #ifdef IS_MPI
224     f_Row(1,atom1) = f_Row(1,atom1) + dudx
225     f_Row(2,atom1) = f_Row(2,atom1) + dudy
226     f_Row(3,atom1) = f_Row(3,atom1) + dudz
227    
228     f_Col(1,atom2) = f_Col(1,atom2) - dudx
229     f_Col(2,atom2) = f_Col(2,atom2) - dudy
230     f_Col(3,atom2) = f_Col(3,atom2) - dudz
231    
232     t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
233     t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
234     t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
235    
236     t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
237     t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
238     t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
239     #else
240     f(1,atom1) = f(1,atom1) + dudx
241     f(2,atom1) = f(2,atom1) + dudy
242     f(3,atom1) = f(3,atom1) + dudz
243    
244     f(1,atom2) = f(1,atom2) - dudx
245     f(2,atom2) = f(2,atom2) - dudy
246     f(3,atom2) = f(3,atom2) - dudz
247    
248     t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
249     t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
250     t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
251    
252     t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
253     t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
254     t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
255     #endif
256    
257     #ifdef IS_MPI
258     id1 = AtomRowToGlobal(atom1)
259     id2 = AtomColToGlobal(atom2)
260     #else
261     id1 = atom1
262     id2 = atom2
263     #endif
264    
265     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
266    
267     fpair(1) = fpair(1) + dudx
268     fpair(2) = fpair(2) + dudy
269     fpair(3) = fpair(3) + dudz
270    
271     endif
272    
273     return
274     end subroutine do_dipole_pair
275    
276     end module dipole_dipole
277 gezelter 1633