ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90
(Generate patch)

Comparing trunk/OOPSE-4/src/UseTheForce/DarkSide/LJ.F90 (file contents):
Revision 1628 by gezelter, Thu Oct 21 20:15:31 2004 UTC vs.
Revision 1633 by gezelter, Fri Oct 22 20:22:48 2004 UTC

# Line 1 | Line 1
1   !! Calculates Long Range forces Lennard-Jones interactions.
2   !! @author Charles F. Vardeman II
3   !! @author Matthew Meineke
4 < !! @version $Id: LJ.F90,v 1.3 2004-10-21 20:15:25 gezelter Exp $, $Date: 2004-10-21 20:15:25 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
4 > !! @version $Id: LJ.F90,v 1.4 2004-10-22 20:22:47 gezelter Exp $, $Date: 2004-10-22 20:22:47 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
5  
6   module lj
7    use atype_module
# Line 51 | Line 51 | module lj
51    public :: useGeometricMixing
52    public :: do_lj_pair
53    public :: newLJtype  
54 +  public :: getSigma
55 +  public :: getEpsilon
56    
57   contains
58  
# Line 93 | Line 95 | contains
95      ParameterMap(ident)%sigma = sigma
96      
97    end subroutine newLJtype
98 +
99 +  function getSigma(atid) result (s)
100 +    integer, intent(in) :: atid
101 +    integer :: localError
102 +    real(kind=dp) :: s
103      
104 +    if (.not.allocated(ParameterMap)) then
105 +       call handleError("LJ", "no ParameterMap was present before first call of getSigma!")
106 +       return
107 +    end if
108 +    
109 +    s = ParameterMap(atid)%sigma
110 +  end function getSigma
111 +
112 +  function getEpsilon(atid) result (e)
113 +    integer, intent(in) :: atid
114 +    integer :: localError
115 +    real(kind=dp) :: e
116 +    
117 +    if (.not.allocated(ParameterMap)) then
118 +       call handleError("dipole-dipole", "no ParameterMap was present before first call of getEpsilon!")
119 +       return
120 +    end if
121 +    
122 +    e = ParameterMap(atid)%epsilon
123 +  end function getEpsilon
124 +
125 +
126    subroutine setCutoffLJ(rcut, do_shift, status)
127      logical, intent(in):: do_shift
128      integer :: status, myStatus

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines