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