ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/simulation_module.F90
Revision: 295
Committed: Thu Mar 6 19:57:03 2003 UTC (21 years, 6 months ago) by chuckv
File size: 3486 byte(s)
Log Message:
Split force loop into subroutines. Moved wrap out of simulation module.

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 chuckv 295
21 mmeineke 270 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     interface getBox
33     module procedure getBox_3d
34     module procedure getBox_dim
35     end interface
36    
37    
38    
39    
40     contains
41    
42 chuckv 290 subroutine setSimulation(setThisSim,error)
43     type (simtype) :: setThisSim
44     integer :: error
45 mmeineke 270 integer :: alloc_stat
46 chuckv 290
47     error = 0
48 mmeineke 270 setSim = .true.
49    
50 chuckv 290 ! copy C struct into fortran type
51     thisSim = setThisSim
52 mmeineke 270 end subroutine setSimulation
53    
54     function getNparticles() result(nparticles)
55     integer :: nparticles
56     nparticles = thisSim%nLRparticles
57     end function getNparticles
58    
59    
60     subroutine change_box_size(new_box_size)
61     real(kind=dp), dimension(3) :: new_box_size
62    
63     thisSim%box = new_box_size
64    
65     end subroutine change_box_size
66    
67    
68     function getBox_3d() result(thisBox)
69     real( kind = dp ), dimension(3) :: thisBox
70     thisBox = thisSim%box
71     end function getBox_3d
72    
73     function getBox_dim(dim) result(thisBox)
74     integer, intent(in) :: dim
75     real( kind = dp ) :: thisBox
76    
77     thisBox = thisSim%box(dim)
78     end function getBox_dim
79    
80    
81    
82     subroutine getRcut(thisrcut,rcut2,rcut6,status)
83     real( kind = dp ), intent(out) :: thisrcut
84     real( kind = dp ), intent(out), optional :: rcut2
85     real( kind = dp ), intent(out), optional :: rcut6
86     integer, optional :: status
87    
88     if (present(status)) status = 0
89    
90     if (.not.setSim ) then
91     if (present(status)) status = -1
92     return
93     end if
94    
95     thisrcut = thisSim%rcut
96     if(present(rcut2)) rcut2 = thisSim%rcutsq
97     if(present(rcut6)) rcut6 = thisSim%rcut6
98    
99     end subroutine getRcut
100    
101    
102    
103    
104     subroutine getRlist(thisrlist,rlist2,status)
105     real( kind = dp ), intent(out) :: thisrlist
106     real( kind = dp ), intent(out), optional :: rlist2
107    
108     integer, optional :: status
109    
110     if (present(status)) status = 0
111    
112     if (.not.setSim ) then
113     if (present(status)) status = -1
114     return
115     end if
116    
117     thisrlist = thisSim%rlist
118     if(present(rlist2)) rlist2 = thisSim%rlistsq
119    
120    
121     end subroutine getRlist
122    
123    
124    
125     pure function getNlocal() result(nlocal)
126     integer :: nlocal
127     nlocal = thisSim%nLRparticles
128     end function getNlocal
129    
130    
131 mmeineke 285 function isEnsemble(this_ensemble) result(is_this_ensemble)
132     character(len = *) :: this_ensemble
133     logical :: is_this_enemble
134     is_this_ensemble = .false.
135     if (this_ensemble == thisSim%ensemble) is_this_ensemble = .true.
136     end function isEnsemble
137 mmeineke 270
138 mmeineke 285 function returnEnsemble() result(thisEnsemble)
139     character (len = len(thisSim%ensemble)) :: thisEnsemble
140     thisEnsemble = thisSim%ensemble
141     end function returnEnsemble
142    
143     function returnMixingRules() result(thisMixingRule)
144     character (len = len(thisSim%ensemble)) :: thisMixingRule
145     thisMixingRule = thisSim%MixingRule
146     end function returnMixingRules
147    
148     function isPBC() result(PBCset)
149     logical :: PBCset
150     PBCset = .false.
151     if (thisSim%use_pbc) PBCset = .true.
152     end function isPBC
153    
154     pure function getStringLen() result (thislen)
155     integer :: thislen
156     thislen = string_len
157     end function setStringLen
158    
159 mmeineke 270 end module simulation