ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/fForceOptions.F90
Revision: 1286
Committed: Wed Sep 10 17:57:55 2008 UTC (16 years, 10 months ago) by gezelter
File size: 5019 byte(s)
Log Message:
Added support for scaled 1-2, 1-3 and 1-4 interactions.

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. 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 !! fForceOptions.F90
43 !! OOPSE-2.0
44 !!
45 !! Created by Charles F. Vardeman II on 12/6/05.
46 !!
47 !! PURPOSE:
48 !!
49 !! @author Charles F. Vardeman II
50 !! @version $Id: fForceOptions.F90,v 1.4 2008-09-10 17:57:55 gezelter Exp $
51
52 !! Handles Mixing options for Fortran.
53
54 module fForceOptions
55 use definitions
56 implicit none
57 PRIVATE
58
59 #define __FORTRAN90
60 #include "UseTheForce/fForceOptions.h"
61
62 type(ForceOptions), save :: fortranForceOptions
63 logical, save :: haveForceOptions = .false.
64 real(kind=dp), save, dimension(0:3) :: vdwScale
65 real(kind=dp), save, dimension(0:3) :: electrostaticScale
66
67
68 public :: ForceOptions
69 public :: vdwScale
70 public :: electrostaticScale
71 public :: getGayBerneMu
72 public :: getGayBerneNu
73 public :: getEnergyMixingRule
74 public :: getDistanceMixingRule
75 public :: usesGeometricDistanceMixing
76 public :: usesGeometricEnergyMixing
77 public :: setForceOptions
78
79 contains
80
81 subroutine setForceOptions(theseOptions)
82 type(ForceOptions),intent(in) :: theseOptions
83 fortranForceOptions = theseOptions
84 haveForceOptions = .true.
85
86 vdwScale(0) = 1.0_dp ! default is topologically unconnected
87 vdwScale(1) = theseOptions%vdw12scale
88 vdwScale(2) = theseOptions%vdw13scale
89 vdwScale(3) = theseOptions%vdw14scale
90
91 electrostaticScale(0) = 1.0_dp ! default is topologically unconnected
92 electrostaticScale(1) = theseOptions%electrostatic12scale
93 electrostaticScale(2) = theseOptions%electrostatic13scale
94 electrostaticScale(3) = theseOptions%electrostatic14scale
95 end subroutine setForceOptions
96
97 function getGayBerneMu() result(thisMu)
98 real(kind=dp) :: thisMu
99 thisMu = fortranForceOptions%GayBerneMu
100 end function getGayBerneMu
101
102 function getGayBerneNu() result(thisNu)
103 real(kind=dp) :: thisNu
104 thisNu = fortranForceOptions%GayBerneNu
105 end function getGayBerneNu
106
107 function usesGeometricDistanceMixing() result(doesit)
108 logical :: doesit
109 doesit = .false.
110 if (.not.haveForceOptions) return
111 if (fortranForceOptions%DistanceMixingRule == GEOMETRIC_MIXING_RULE) then
112 doesit = .true.
113 endif
114 end function usesGeometricDistanceMixing
115
116 function usesGeometricEnergyMixing() result(doesit)
117 logical :: doesit
118 doesit = .false.
119 if (.not.haveForceOptions) return
120 if (fortranForceOptions%EnergyMixingRule == GEOMETRIC_MIXING_RULE) then
121 doesit = .true.
122 endif
123 end function usesGeometricEnergyMixing
124
125 function getEnergyMixingRule() result(MixingRule)
126 integer :: MixingRule
127 MixingRule = 0
128 if (.not.haveForceOptions) return
129 MixingRule = fortranForceOptions%EnergyMixingRule
130 end function getEnergyMixingRule
131
132 function getDistanceMixingRule() result(MixingRule)
133 integer :: MixingRule
134 MixingRule = 0
135 if (.not.haveForceOptions) return
136 MixingRule = fortranForceOptions%DistanceMixingRule
137 end function getDistanceMixingRule
138
139 end module fForceOptions