ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/simulation_module.F90
Revision: 290
Committed: Thu Feb 27 21:25:47 2003 UTC (21 years, 6 months ago) by chuckv
File size: 4391 byte(s)
Log Message:
Made changes to lj_FF.cpp to pass stress tensor.

File Contents

# User Rev Content
1 chuckv 290 !! Fortran interface to C entry plug.
2    
3 mmeineke 270 module simulation
4     use definitions, ONLY :dp
5     #ifdef IS_MPI
6     use mpiSimulation
7     #endif
8    
9     implicit none
10     PRIVATE
11    
12 mmeineke 285 #define __FORTRAN90
13     #include "../headers/fsimulation.h"
14 mmeineke 270
15     type (simtype), public :: thisSim
16    
17 chuckv 290 logical :: setSim = .false.
18 mmeineke 270
19    
20     public :: wrap
21     public :: getBox
22     public :: getRcut
23     public :: getRlist
24     public :: getNlocal
25     public :: setSimulation
26 mmeineke 285 public :: isEnsemble
27     public :: isPBC
28     public :: getStringLen
29     public :: returnMixingRules
30    
31 mmeineke 270 ! public :: setRcut
32    
33     interface wrap
34     module procedure wrap_1d
35     module procedure wrap_3d
36     end interface
37    
38     interface getBox
39     module procedure getBox_3d
40     module procedure getBox_dim
41     end interface
42    
43    
44    
45    
46     contains
47    
48 chuckv 290 subroutine setSimulation(setThisSim,error)
49     type (simtype) :: setThisSim
50     integer :: error
51 mmeineke 270 integer :: alloc_stat
52 chuckv 290
53     error = 0
54 mmeineke 270 setSim = .true.
55    
56 chuckv 290 ! copy C struct into fortran type
57     thisSim = setThisSim
58 mmeineke 270 end subroutine setSimulation
59    
60     function getNparticles() result(nparticles)
61     integer :: nparticles
62     nparticles = thisSim%nLRparticles
63     end function getNparticles
64    
65    
66     subroutine change_box_size(new_box_size)
67     real(kind=dp), dimension(3) :: new_box_size
68    
69     thisSim%box = new_box_size
70    
71     end subroutine change_box_size
72    
73    
74     function getBox_3d() result(thisBox)
75     real( kind = dp ), dimension(3) :: thisBox
76     thisBox = thisSim%box
77     end function getBox_3d
78    
79     function getBox_dim(dim) result(thisBox)
80     integer, intent(in) :: dim
81     real( kind = dp ) :: thisBox
82    
83     thisBox = thisSim%box(dim)
84     end function getBox_dim
85    
86    
87     function wrap_1d(r,dim) result(this_wrap)
88    
89    
90     real( kind = DP ) :: r
91     real( kind = DP ) :: this_wrap
92     integer :: dim
93    
94     if (use_pbc) then
95     ! this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
96     this_wrap = r - thisSim%box(dim)*nint(r/thisSim%box(dim))
97     else
98     this_wrap = r
99     endif
100    
101     return
102     end function wrap_1d
103    
104     function wrap_3d(r) result(this_wrap)
105     real( kind = dp ), dimension(3), intent(in) :: r
106     real( kind = dp ), dimension(3) :: this_wrap
107    
108    
109 mmeineke 285 if (this_sim%use_pbc) then
110 mmeineke 270 ! this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
111     this_wrap = r - thisSim%box*nint(r/thisSim%box)
112     else
113     this_wrap = r
114     endif
115     end function wrap_3d
116    
117    
118    
119     subroutine getRcut(thisrcut,rcut2,rcut6,status)
120     real( kind = dp ), intent(out) :: thisrcut
121     real( kind = dp ), intent(out), optional :: rcut2
122     real( kind = dp ), intent(out), optional :: rcut6
123     integer, optional :: status
124    
125     if (present(status)) status = 0
126    
127     if (.not.setSim ) then
128     if (present(status)) status = -1
129     return
130     end if
131    
132     thisrcut = thisSim%rcut
133     if(present(rcut2)) rcut2 = thisSim%rcutsq
134     if(present(rcut6)) rcut6 = thisSim%rcut6
135    
136     end subroutine getRcut
137    
138    
139    
140    
141     subroutine getRlist(thisrlist,rlist2,status)
142     real( kind = dp ), intent(out) :: thisrlist
143     real( kind = dp ), intent(out), optional :: rlist2
144    
145     integer, optional :: status
146    
147     if (present(status)) status = 0
148    
149     if (.not.setSim ) then
150     if (present(status)) status = -1
151     return
152     end if
153    
154     thisrlist = thisSim%rlist
155     if(present(rlist2)) rlist2 = thisSim%rlistsq
156    
157    
158     end subroutine getRlist
159    
160    
161    
162     pure function getNlocal() result(nlocal)
163     integer :: nlocal
164     nlocal = thisSim%nLRparticles
165     end function getNlocal
166    
167    
168 mmeineke 285 function isEnsemble(this_ensemble) result(is_this_ensemble)
169     character(len = *) :: this_ensemble
170     logical :: is_this_enemble
171     is_this_ensemble = .false.
172     if (this_ensemble == thisSim%ensemble) is_this_ensemble = .true.
173     end function isEnsemble
174 mmeineke 270
175 mmeineke 285 function returnEnsemble() result(thisEnsemble)
176     character (len = len(thisSim%ensemble)) :: thisEnsemble
177     thisEnsemble = thisSim%ensemble
178     end function returnEnsemble
179    
180     function returnMixingRules() result(thisMixingRule)
181     character (len = len(thisSim%ensemble)) :: thisMixingRule
182     thisMixingRule = thisSim%MixingRule
183     end function returnMixingRules
184    
185     function isPBC() result(PBCset)
186     logical :: PBCset
187     PBCset = .false.
188     if (thisSim%use_pbc) PBCset = .true.
189     end function isPBC
190    
191     pure function getStringLen() result (thislen)
192     integer :: thislen
193     thislen = string_len
194     end function setStringLen
195    
196 mmeineke 270 end module simulation