ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/simulation_module.F90
Revision: 285
Committed: Wed Feb 26 18:45:57 2003 UTC (21 years, 6 months ago) by mmeineke
File size: 5505 byte(s)
Log Message:
OOPSE has been braought up to date with the mdtools devfelopment tree. mdtools should no longer be used. All development should take place in the OOPSE tree.

File Contents

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