ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/simulation_module.F90
(Generate patch)

Comparing trunk/mdtools/md_code/simulation_module.F90 (file contents):
Revision 264 by chuckv, Tue Feb 4 20:16:08 2003 UTC vs.
Revision 281 by chuckv, Mon Feb 24 21:25:15 2003 UTC

# Line 7 | Line 7 | module simulation
7    implicit none
8    PRIVATE
9  
10 + #define __FORTRAN90
11 + #include "../headers/fsimulation.h"
12  
11
12  type, public :: simtype
13     PRIVATE
14 !     SEQUENCE
15 !! Number of particles on this processor
16     integer :: nLRparticles
17 !! Periodic Box    
18     real ( kind = dp ), dimension(3) :: box
19 !! List Cutoff    
20     real ( kind = dp ) :: rlist = 0.0_dp
21 !! Radial cutoff
22     real ( kind = dp ) :: rcut  = 0.0_dp
23 !! List cutoff squared
24     real ( kind = dp ) :: rlistsq = 0.0_dp
25 !! Radial Cutoff squared
26     real ( kind = dp ) :: rcutsq  = 0.0_dp
27 !! Radial Cutoff^6
28     real ( kind = dp ) :: rcut6  = 0.0_dp
29
30  end type simtype
31
13    type (simtype), public :: thisSim
14   !! Tag for MPI calculations  
15    integer, allocatable, dimension(:) :: tag
# Line 39 | Line 20 | module simulation
20   #endif
21  
22   !! WARNING: use_pbc hardcoded, fixme
23 <  logical :: use_pbc = .true.
43 <  logical :: setSim = .false.
23 >   logical :: setSim = .false.
24  
25   !! array for saving previous positions for neighbor lists.  
26    real( kind = dp ), allocatable,dimension(:,:),save :: q0
# Line 54 | Line 34 | module simulation
34    public :: getRlist
35    public :: getNlocal
36    public :: setSimulation
37 +  public :: isEnsemble
38 +  public :: isPBC
39 +  public :: getStringLen
40 +  public :: returnMixingRules
41 +
42   !  public :: setRcut
43  
44    interface wrap
# Line 77 | Line 62 | contains
62  
63   contains
64  
65 <  subroutine setSimulation(nLRParticles,box,rlist,rcut)
65 >  subroutine setSimulation(nLRParticles,box,rlist,rcut,ensemble,mixingRule,use_pbc)
66      integer, intent(in) :: nLRParticles
67      real(kind = dp ), intent(in), dimension(3) :: box
68      real(kind = dp ), intent(in) :: rlist
69      real(kind = dp ), intent(in) :: rcut
70 +    character( len = stringLen), intent(in)  :: ensemble
71 +    character( len = stringLen), intent(in)  :: mixingRule
72 +    logical, intent(in) :: use_pbc
73      integer :: alloc_stat
74      if( setsim ) return  ! simulation is already initialized
75      setSim = .true.
# Line 94 | Line 82 | contains
82      thisSim%rcutsq       = rcut * rcut
83      thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
84      
85 +    thisSim%ensemble = ensemble
86 +    thisSim%mixingRule = mixingRule
87 +    thisSim%use_pbc = use_pbc
88 +
89      if (.not. allocated(q0)) then
90         allocate(q0(3,nLRParticles),stat=alloc_stat)
91      endif
# Line 205 | Line 197 | contains
197      real( kind = dp ), dimension(3) :: this_wrap
198  
199      
200 <    if (use_pbc) then
200 >    if (this_sim%use_pbc) then
201         !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
202         this_wrap = r - thisSim%box*nint(r/thisSim%box)
203      else
# Line 264 | Line 256 | contains
256    end function getNlocal
257  
258  
259 +  function isEnsemble(this_ensemble) result(is_this_ensemble)
260 +    character(len = *) :: this_ensemble
261 +    logical :: is_this_enemble
262 +    is_this_ensemble = .false.
263 +    if (this_ensemble == thisSim%ensemble) is_this_ensemble = .true.
264 +  end function isEnsemble
265  
266 +  function returnEnsemble() result(thisEnsemble)
267 +    character (len = len(thisSim%ensemble)) :: thisEnsemble
268 +    thisEnsemble = thisSim%ensemble
269 +  end function returnEnsemble
270 +
271 +  function returnMixingRules() result(thisMixingRule)
272 +    character (len = len(thisSim%ensemble)) :: thisMixingRule
273 +    thisMixingRule = thisSim%MixingRule
274 +  end function returnMixingRules
275 +
276 +  function isPBC() result(PBCset)
277 +    logical :: PBCset
278 +    PBCset = .false.
279 +    if (thisSim%use_pbc) PBCset = .true.
280 +  end function isPBC
281 +
282 +  pure function getStringLen() result (thislen)
283 +    integer :: thislen    
284 +    thislen = string_len
285 +  end function setStringLen
286 +
287   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines