ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/charge.F90
Revision: 1948
Committed: Fri Jan 14 20:31:16 2005 UTC (19 years, 7 months ago) by gezelter
File size: 6072 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 charge_charge
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 = 332.06508_DP
58     logical, save :: haveChargeMap = .false.
59    
60 gezelter 1633 public :: newChargeType
61     public :: do_charge_pair
62     public :: getCharge
63    
64 gezelter 1608 type :: ChargeList
65 gezelter 1930 integer :: c_ident
66 gezelter 1608 real(kind=DP) :: charge = 0.0_DP
67     end type ChargeList
68    
69     type(ChargeList), dimension(:), allocatable :: ChargeMap
70    
71     contains
72 gezelter 1633
73 gezelter 1930 subroutine newChargeType(c_ident, charge, status)
74     integer,intent(in) :: c_ident
75 gezelter 1633 real(kind=dp),intent(in) :: charge
76     integer,intent(out) :: status
77 gezelter 1930 integer :: nAtypes, myATID
78 gezelter 1633
79 gezelter 1608 status = 0
80    
81 gezelter 1930 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
82    
83 gezelter 1633 !! Be simple-minded and assume that we need a ChargeMap that
84     !! is the same size as the total number of atom types
85    
86     if (.not.allocated(ChargeMap)) then
87    
88     nAtypes = getSize(atypes)
89 gezelter 1608
90 gezelter 1633 if (nAtypes == 0) then
91     status = -1
92     return
93     end if
94    
95     if (.not. allocated(ChargeMap)) then
96     allocate(ChargeMap(nAtypes))
97     endif
98    
99     end if
100    
101 gezelter 1930 if (myATID .gt. size(ChargeMap)) then
102 gezelter 1608 status = -1
103     return
104 gezelter 1633 endif
105 gezelter 1608
106 gezelter 1633 ! set the values for ChargeMap for this atom type:
107 gezelter 1608
108 gezelter 1930 ChargeMap(myATID)%c_ident = c_ident
109     ChargeMap(myATID)%charge = charge
110 gezelter 1608
111 gezelter 1633 end subroutine newChargeType
112    
113     function getCharge(atid) result (c)
114     integer, intent(in) :: atid
115     integer :: localError
116     real(kind=dp) :: c
117 gezelter 1608
118 gezelter 1633 if (.not.allocated(ChargeMap)) then
119     call handleError("charge_charge", "no ChargeMap was present before first call of getCharge!")
120     return
121     end if
122 gezelter 1608
123 gezelter 1633 c = ChargeMap(atid)%charge
124     end function getCharge
125    
126 gezelter 1608 subroutine do_charge_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
127     pot, f, do_pot)
128    
129     logical :: do_pot
130    
131     integer atom1, atom2, me1, me2, id1, id2
132     integer :: localError
133     real(kind=dp) :: rij, r2, q1, q2, sw, vpair
134     real(kind=dp) :: drdx, drdy, drdz, dudr, fx, fy, fz
135     real(kind=dp) :: vterm
136    
137     real( kind = dp ) :: pot
138     real( kind = dp ), dimension(3) :: d, fpair
139     real( kind = dp ), dimension(3,nLocal) :: f
140    
141    
142 gezelter 1633 if (.not.allocated(ChargeMap)) then
143     call handleError("charge_charge", "no ChargeMap was present before first call of do_charge_pair!")
144     return
145     end if
146    
147 gezelter 1608 #ifdef IS_MPI
148     me1 = atid_Row(atom1)
149     me2 = atid_Col(atom2)
150     #else
151     me1 = atid(atom1)
152     me2 = atid(atom2)
153     #endif
154    
155     q1 = ChargeMap(me1)%charge
156     q2 = ChargeMap(me2)%charge
157    
158     vterm = pre * q1 * q2 / rij
159    
160     dudr = -sw * vterm / rij
161    
162     drdx = d(1) / rij
163     drdy = d(2) / rij
164     drdz = d(3) / rij
165    
166     fx = dudr * drdx
167     fy = dudr * drdy
168     fz = dudr * drdz
169    
170     vpair = vpair + vterm
171    
172     #ifdef IS_MPI
173     if (do_pot) then
174     pot_Row(atom1) = pot_Row(atom1) + sw*vterm*0.5
175     pot_Col(atom2) = pot_Col(atom2) + sw*vterm*0.5
176     endif
177    
178     f_Row(1,atom1) = f_Row(1,atom1) + fx
179     f_Row(2,atom1) = f_Row(2,atom1) + fy
180     f_Row(3,atom1) = f_Row(3,atom1) + fz
181    
182     f_Col(1,atom2) = f_Col(1,atom2) - fx
183     f_Col(2,atom2) = f_Col(2,atom2) - fy
184     f_Col(3,atom2) = f_Col(3,atom2) - fz
185    
186     #else
187    
188     if (do_pot) pot = pot + sw*vterm
189    
190     f(1,atom1) = f(1,atom1) + fx
191     f(2,atom1) = f(2,atom1) + fy
192     f(3,atom1) = f(3,atom1) + fz
193    
194     f(1,atom2) = f(1,atom2) - fx
195     f(2,atom2) = f(2,atom2) - fy
196     f(3,atom2) = f(3,atom2) - fz
197    
198     #endif
199    
200     #ifdef IS_MPI
201     id1 = AtomRowToGlobal(atom1)
202     id2 = AtomColToGlobal(atom2)
203     #else
204     id1 = atom1
205     id2 = atom2
206     #endif
207    
208     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
209    
210     fpair(1) = fpair(1) + fx
211     fpair(2) = fpair(2) + fy
212     fpair(3) = fpair(3) + fz
213    
214     endif
215     return
216     end subroutine do_charge_pair
217    
218     end module charge_charge