ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/atype_module.F90
Revision: 378
Committed: Fri Mar 21 17:42:12 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 3176 byte(s)
Log Message:
This commit was generated by cvs2svn to compensate for changes in r377,
which included commits to RCS files with non-trunk default branches.

File Contents

# User Rev Content
1 mmeineke 377 !! module defines atypes available to simulation
2    
3     module atype_module
4     use definitions, only: dp
5     use vector_class
6     use sticky_pair
7     use gb_pair
8     implicit none
9     private
10    
11     type (Vector), pointer, public :: atypes => null()
12    
13    
14     public :: new_atype
15    
16     contains
17    
18     subroutine new_atype(c_ident, is_LJ, is_Sticky, is_DP, is_GB, &
19     lj_epsilon, lj_sigma, &
20     dipole_moment, &
21     sticky_w0, sticky_v0, &
22     GB_sigma, GB_l2b_ratio, GB_eps, GB_eps_ratio, GB_mu, GB_nu, &
23     status)
24    
25     real( kind = dp ), intent(in) :: lj_epsilon
26     real( kind = dp ), intent(in) :: lj_sigma
27     real( kind = dp ), intent(in) :: sticky_w0
28     real( kind = dp ), intent(in) :: sticky_v0
29     real( kind = dp ), intent(in) :: dipole_moment
30     real( kind = dp ), intent(in) :: GB_sigma, GB_l2b_ratio, GB_eps
31     real( kind = dp ), intent(in) :: GB_eps_ratio, GB_mu, GB_nu
32    
33     integer, intent(in) :: c_ident
34     integer, intent(out) :: status
35     integer, intent(in) :: is_Sticky
36     integer, intent(in) :: is_DP
37     integer, intent(in) :: is_GB
38     integer, intent(in) :: is_LJ
39     integer :: me
40     logical :: l_is_LJ, l_is_DP, l_is_Sticky, l_is_GB
41     integer :: FFcheckStatus
42     status = 0
43    
44     if (.not. associated(atypes)) then
45     !! There are 16 properties to worry about for now.
46     !! Fix this if needed for more atomic properties
47     atypes => initialize(16)
48     if (.not.associated(atypes)) then
49     status = -1
50     return
51     endif
52     endif
53    
54     me = addElement(atypes)
55     call setElementProperty(atypes, me, "c_ident", c_ident)
56    
57     l_is_LJ = (is_LJ .ne. 0)
58     l_is_DP = (is_DP .ne. 0)
59     l_is_Sticky = (is_Sticky .ne. 0)
60     l_is_GB = (is_GB .ne. 0)
61    
62     call setElementProperty(atypes, me, "is_LJ", l_is_LJ)
63     call setElementProperty(atypes, me, "is_DP", l_is_DP)
64     call setElementProperty(atypes, me, "is_Sticky", l_is_Sticky)
65     call setElementProperty(atypes, me, "is_GB", l_is_GB)
66    
67     if (l_is_LJ) then
68     call setElementProperty(atypes, me, "lj_sigma", lj_sigma)
69     call setElementProperty(atypes, me, "lj_epsilon", lj_epsilon)
70     endif
71     if (l_is_DP) then
72     call setElementProperty(atypes, me, "dipole_moment", dipole_moment)
73     endif
74     if (l_is_Sticky) then
75     call setElementProperty(atypes, me, "sticky_w0", sticky_w0)
76     call setElementProperty(atypes, me, "sticky_v0", sticky_v0)
77     call check_sticky_FF(FFcheckStatus)
78     if (FFcheckStatus .ne. 0) call set_sticky_params(sticky_w0, sticky_v0)
79     endif
80     if (l_is_GB) then
81     call setElementProperty(atypes, me, "GB_sigma", GB_sigma)
82     call setElementProperty(atypes, me, "GB_l2b_ratio", GB_l2b_ratio)
83     call setElementProperty(atypes, me, "GB_eps", GB_eps)
84     call setElementProperty(atypes, me, "GB_eps_ratio", GB_eps_ratio)
85     call setElementProperty(atypes, me, "GB_mu", GB_mu)
86     call setElementProperty(atypes, me, "GB_nu", GB_nu)
87     call check_gb_pair_FF(FFcheckStatus)
88     if (FFcheckStatus .ne. 0) call set_gb_pair_params(GB_sigma, GB_l2b_ratio, &
89     GB_eps, GB_eps_ratio, GB_mu, GB_nu)
90     endif
91    
92     end subroutine new_atype
93    
94     end module atype_module