ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/fForceOptions.F90
Revision: 1442
Committed: Mon May 10 17:28:26 2010 UTC (15 years, 2 months ago) by gezelter
File size: 4935 byte(s)
Log Message:
Adding property set to svn entries

File Contents

# Content
1 !!
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. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41 !! Created by Charles F. Vardeman II on 12/6/05.
42 !!
43 !! PURPOSE:
44 !!
45 !! @author Charles F. Vardeman II
46 !! @version $Id$
47
48 !! Handles Mixing options for Fortran.
49
50 module fForceOptions
51 use definitions
52 implicit none
53 PRIVATE
54
55 #define __FORTRAN90
56 #include "UseTheForce/fForceOptions.h"
57
58 type(ForceOptions), save :: fortranForceOptions
59 logical, save :: haveForceOptions = .false.
60 real(kind=dp), save, dimension(0:3) :: vdwScale
61 real(kind=dp), save, dimension(0:3) :: electrostaticScale
62
63
64 public :: ForceOptions
65 public :: vdwScale
66 public :: electrostaticScale
67 public :: getGayBerneMu
68 public :: getGayBerneNu
69 public :: getEnergyMixingRule
70 public :: getDistanceMixingRule
71 public :: usesGeometricDistanceMixing
72 public :: usesGeometricEnergyMixing
73 public :: setForceOptions
74
75 contains
76
77 subroutine setForceOptions(theseOptions)
78 type(ForceOptions),intent(in) :: theseOptions
79 fortranForceOptions = theseOptions
80 haveForceOptions = .true.
81
82 vdwScale(0) = 1.0_dp ! default is topologically unconnected
83 vdwScale(1) = theseOptions%vdw12scale
84 vdwScale(2) = theseOptions%vdw13scale
85 vdwScale(3) = theseOptions%vdw14scale
86
87 electrostaticScale(0) = 1.0_dp ! default is topologically unconnected
88 electrostaticScale(1) = theseOptions%electrostatic12scale
89 electrostaticScale(2) = theseOptions%electrostatic13scale
90 electrostaticScale(3) = theseOptions%electrostatic14scale
91 end subroutine setForceOptions
92
93 function getGayBerneMu() result(thisMu)
94 real(kind=dp) :: thisMu
95 thisMu = fortranForceOptions%GayBerneMu
96 end function getGayBerneMu
97
98 function getGayBerneNu() result(thisNu)
99 real(kind=dp) :: thisNu
100 thisNu = fortranForceOptions%GayBerneNu
101 end function getGayBerneNu
102
103 function usesGeometricDistanceMixing() result(doesit)
104 logical :: doesit
105 doesit = .false.
106 if (.not.haveForceOptions) return
107 if (fortranForceOptions%DistanceMixingRule == GEOMETRIC_MIXING_RULE) then
108 doesit = .true.
109 endif
110 end function usesGeometricDistanceMixing
111
112 function usesGeometricEnergyMixing() result(doesit)
113 logical :: doesit
114 doesit = .false.
115 if (.not.haveForceOptions) return
116 if (fortranForceOptions%EnergyMixingRule == GEOMETRIC_MIXING_RULE) then
117 doesit = .true.
118 endif
119 end function usesGeometricEnergyMixing
120
121 function getEnergyMixingRule() result(MixingRule)
122 integer :: MixingRule
123 MixingRule = 0
124 if (.not.haveForceOptions) return
125 MixingRule = fortranForceOptions%EnergyMixingRule
126 end function getEnergyMixingRule
127
128 function getDistanceMixingRule() result(MixingRule)
129 integer :: MixingRule
130 MixingRule = 0
131 if (.not.haveForceOptions) return
132 MixingRule = fortranForceOptions%DistanceMixingRule
133 end function getDistanceMixingRule
134
135 end module fForceOptions

Properties

Name Value
svn:keywords Author Id Revision Date