ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1168
Committed: Fri Jul 13 14:21:07 2007 UTC (18 years ago) by chuckv
Original Path: trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
File size: 9569 byte(s)
Log Message:
One more to add.

File Contents

# User Rev Content
1 chuckv 1168 !!
2     !! Copyright (c) 2007 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    
43     !! Calculates Metal-Non Metal interactions.
44     !! @author Charles F. Vardeman II
45     !! @version $Id: MetalNonMetal.F90,v 1.1 2007-07-13 14:21:07 chuckv Exp $, $Date: 2007-07-13 14:21:07 $, $Name: not supported by cvs2svn $, $Revision: 1.1 $
46    
47    
48     module MetalNonMetal
49     use definitions
50     use atype_module
51     use vector_class
52     use simulation
53     use status
54     use fForceOptions
55     #ifdef IS_MPI
56     use mpiSimulation
57     #endif
58     use force_globals
59    
60     implicit none
61     PRIVATE
62     #define __FORTRAN90
63     #include "UseTheForce/DarkSide/fInteractionMap.h"
64     #include "UseTheForce/DarkSide/fMnMInteractions.h"
65    
66     logical, save :: useGeometricDistanceMixing = .false.
67     logical, save :: haveInteractionLookup = .false.
68    
69     real(kind=DP), save :: defaultCutoff = 0.0_DP
70     logical, save :: defaultShiftPot = .false.
71     logical, save :: defaultShiftFrc = .false.
72     logical, save :: haveDefaultCutoff = .false.
73    
74     type :: MnMinteraction
75     integer :: metal_atid
76     integer :: nonmetal_atid
77     integer :: interaction_type
78     real(kind=dp) :: R0
79     real(kind=dp) :: D0
80     real(kind=dp) :: beta0
81     real(kind=dp) :: betaH
82     real(kind=dp) :: alpha
83     real(kind=dp) :: gamma
84     real(kind=dp) :: sigma
85     real(kind=dp) :: epsilon
86     real(kind=dp) :: rCut = 0.0_dp
87     logical :: rCutWasSet = .false.
88     logical :: shiftedPot = .false.
89     logical :: shiftedFrc = .false.
90     end type MnMinteraction
91    
92     type :: MnMinteractionMap
93     PRIVATE
94     integer :: initialCapacity = 10
95     integer :: capacityIncrement = 0
96     integer :: interactionCount = 0
97     type(MnMinteraction), pointer :: interactions(:) => null()
98     end type MnMinteractionMap
99    
100     type (MnMInteractionMap), pointer :: MnM_Map
101    
102     integer, allocatable, dimension(:,:) :: MnMinteractionLookup
103    
104     public :: setMnMDefaultCutoff
105     public :: addInteraction
106     public :: deleteInteractions
107     public :: MNMtype
108    
109    
110     contains
111    
112    
113    
114     subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
115     real(kind=dp), intent(in) :: thisRcut
116     logical, intent(in) :: shiftedPot
117     logical, intent(in) :: shiftedFrc
118     integer i, nInteractions
119     defaultCutoff = thisRcut
120     defaultShiftPot = shiftedPot
121     defaultShiftFrc = shiftedFrc
122    
123     if(MnM_Map%interactionCount /= 0) then
124     nInteractions = MnM_Map%interactionCount
125    
126     do i = 1, nInteractions
127     MnM_Map%interactions(i)%shiftedPot = shiftedPot
128     MnM_Map%interactions(i)%shiftedFrc = shiftedFrc
129     MnM_Map%interactions(i)%rCut = thisRcut
130     MnM_Map%interactions(i)%rCutWasSet = .true.
131     enddo
132     end if
133    
134     end subroutine setMnMDefaultCutoff
135    
136     subroutine copyAllData(v1, v2)
137     type(MnMinteractionMap), pointer :: v1
138     type(MnMinteractionMap), pointer :: v2
139     integer :: i, j
140    
141     do i = 1, v1%interactionCount
142     v2%interactions(i) = v1%interactions(i)
143     enddo
144    
145     v2%interactionCount = v1%interactionCount
146     return
147     end subroutine copyAllData
148    
149     subroutine addInteraction(myInteraction)
150     type(MNMtype) :: myInteraction
151     type(MnMinteraction) :: nt
152     integer :: id
153    
154     nt%interaction_type = myInteraction%MNMInteractionType
155     nt%metal_atid = myInteraction%metal_atid
156     nt%nonmetal_atid = myInteraction%nonmetal_atid
157    
158     select case (nt%interaction_type)
159     case (MNM_LENNARDJONES)
160     nt%sigma = myInteraction%sigma
161     nt%epsilon = myInteraction%epsilon
162     case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
163     nt%R0 = myInteraction%R0
164     nt%D0 = myInteraction%D0
165     nt%beta0 = myInteraction%beta0
166     case(MNM_MAW)
167     nt%R0 = myInteraction%R0
168     nt%D0 = myInteraction%D0
169     nt%beta0 = myInteraction%beta0
170     nt%betaH = myInteraction%betaH
171     nt%alpha = myInteraction%alpha
172     nt%gamma = myInteraction%gamma
173     case default
174     write(*,*) 'unknown MnM interaction type!'
175     end select
176    
177     if (.not. associated(MnM_Map)) then
178     call ensureCapacityHelper(MnM_Map, 1)
179     else
180     call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1)
181     end if
182    
183     MnM_Map%interactionCount = MnM_Map%interactionCount + 1
184     id = MnM_Map%interactionCount
185     MnM_Map%interactions(id) = nt
186     end subroutine addInteraction
187    
188     subroutine ensureCapacityHelper(this, minCapacity)
189     type(MnMinteractionMap), pointer :: this, that
190     integer, intent(in) :: minCapacity
191     integer :: oldCapacity
192     integer :: newCapacity
193     logical :: resizeFlag
194    
195     resizeFlag = .false.
196    
197     ! first time: allocate a new vector with default size
198    
199     if (.not. associated(this)) then
200     this => MnMinitialize(minCapacity, 0)
201     endif
202    
203     oldCapacity = size(this%interactions)
204    
205     if (minCapacity > oldCapacity) then
206     if (this%capacityIncrement .gt. 0) then
207     newCapacity = oldCapacity + this%capacityIncrement
208     else
209     newCapacity = oldCapacity * 2
210     endif
211     if (newCapacity .lt. minCapacity) then
212     newCapacity = minCapacity
213     endif
214     resizeFlag = .true.
215     else
216     newCapacity = oldCapacity
217     endif
218    
219     if (resizeFlag) then
220     that => MnMinitialize(newCapacity, this%capacityIncrement)
221     call copyAllData(this, that)
222     this => MnMdestroy(this)
223     this => that
224     endif
225     end subroutine ensureCapacityHelper
226    
227     function MnMinitialize(cap, capinc) result(this)
228     integer, intent(in) :: cap, capinc
229     integer :: error
230     type(MnMinteractionMap), pointer :: this
231    
232     nullify(this)
233    
234     if (cap < 0) then
235     write(*,*) 'Bogus Capacity:', cap
236     return
237     endif
238     allocate(this,stat=error)
239     if ( error /= 0 ) then
240     write(*,*) 'Could not allocate MnMinteractionMap!'
241     return
242     end if
243    
244     this%initialCapacity = cap
245     this%capacityIncrement = capinc
246    
247     allocate(this%interactions(this%initialCapacity), stat=error)
248     if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!'
249    
250     end function MnMinitialize
251    
252     subroutine createInteractionLookup(this)
253     type(MnMinteractionMap), pointer :: this
254     integer :: biggestAtid, i, metal_atid, nonmetal_atid, error
255    
256     biggestAtid=-1
257     do i = 1, this%interactionCount
258     metal_atid = this%interactions(i)%metal_atid
259     nonmetal_atid = this%interactions(i)%nonmetal_atid
260    
261     if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid
262     if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid
263     enddo
264    
265     allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error)
266     if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup'
267    
268     do i = 1, this%interactionCount
269     metal_atid = this%interactions(i)%metal_atid
270     nonmetal_atid = this%interactions(i)%nonmetal_atid
271    
272     MnMinteractionLookup(metal_atid, nonmetal_atid) = i
273     MnMinteractionLookup(nonmetal_atid, metal_atid) = i
274     enddo
275     end subroutine createInteractionLookup
276    
277    
278     function MnMdestroy(this) result(null_this)
279     logical :: done
280     type(MnMinteractionMap), pointer :: this
281     type(MnMinteractionMap), pointer :: null_this
282    
283     if (.not. associated(this)) then
284     null_this => null()
285     return
286     end if
287    
288     !! Walk down the list and deallocate each of the map's components
289     if(associated(this%interactions)) then
290     deallocate(this%interactions)
291     this%interactions=>null()
292     endif
293     deallocate(this)
294     this => null()
295     null_this => null()
296     end function MnMdestroy
297    
298    
299     subroutine deleteInteractions()
300     MnM_Map => MnMdestroy(MnM_Map)
301     return
302     end subroutine deleteInteractions
303    
304     end module MetalNonMetal