ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/fForceOptions.F90
Revision: 2496
Committed: Wed Dec 7 19:46:56 2005 UTC (18 years, 7 months ago) by chuckv
File size: 4320 byte(s)
Log Message:
Further fixes for compiler issues.

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.2 2005-12-07 19:46:56 chuckv Exp $
51
52 !! Handles Mixing options for Fortran.
53
54
55
56 module fForceOptions
57 use definitions
58 implicit none
59 PRIVATE
60
61 #define __FORTRAN90
62 #include "UseTheForce/fForceOptions.h"
63
64 type(ForceOptions), save :: fortranForceOptions
65 logical, save :: haveForceOptions = .false.
66
67
68 public :: ForceOptions
69 public :: getVDW14Scale
70 public :: getElectrostatic14Scale
71 public :: getEnergyMixingRule
72 public :: getDistanceMixingRule
73 public :: usesGeometricDistanceMixing
74 public :: usesGeometricEnergyMixing
75 public :: setForceOptions
76
77
78 contains
79
80 subroutine setForceOptions(theseOptions)
81 type(ForceOptions),intent(in) :: theseOptions
82 fortranForceOptions = theseOptions
83 haveForceOptions = .true.
84 end subroutine setForceOptions
85
86
87 function getVDW14Scale() result(thisScale)
88 real(kind=dp) :: thisScale
89 thisScale = fortranForceOptions%vdw14scale
90 end function getVDW14Scale
91
92 function getElectrostatic14Scale() result(thisScale)
93 real(kind=dp) :: thisScale
94 thisScale = fortranForceOptions%electrostatic14scale
95 end function getElectrostatic14Scale
96
97
98 function usesGeometricDistanceMixing() result(doesit)
99 logical :: doesit
100 doesit = .false.
101 if (.not.haveForceOptions) return
102 if (fortranForceOptions%DistanceMixingRule == GEOMETRIC_MIXING_RULE) then
103 doesit = .true.
104 endif
105
106 end function usesGeometricDistanceMixing
107
108
109 function usesGeometricEnergyMixing() result(doesit)
110 logical :: doesit
111 doesit = .false.
112 if (.not.haveForceOptions) return
113 if (fortranForceOptions%EnergyMixingRule == GEOMETRIC_MIXING_RULE) then
114 doesit = .true.
115 endif
116 end function usesGeometricEnergyMixing
117
118
119 function getEnergyMixingRule() result(MixingRule)
120 integer :: MixingRule
121 MixingRule = 0
122 if (.not.haveForceOptions) return
123 MixingRule = fortranForceOptions%EnergyMixingRule
124 end function getEnergyMixingRule
125
126 function getDistanceMixingRule() result(MixingRule)
127 integer :: MixingRule
128 MixingRule = 0
129 if (.not.haveForceOptions) return
130 MixingRule = fortranForceOptions%DistanceMixingRule
131 end function getDistanceMixingRule
132
133
134
135 end module fForceOptions