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

Comparing trunk/OOPSE_old/src/mdtools/libmdCode/simulation_module.F90 (file contents):
Revision 285 by mmeineke, Wed Feb 26 18:45:57 2003 UTC vs.
Revision 312 by gezelter, Tue Mar 11 17:46:18 2003 UTC

# Line 1 | Line 1
1 + !! Fortran interface to C entry plug.
2 +
3   module simulation
4    use definitions, ONLY :dp
5   #ifdef IS_MPI
# Line 8 | Line 10 | module simulation
10    PRIVATE
11  
12   #define __FORTRAN90
13 < #include "../headers/fsimulation.h"
13 > #include "fSimulation.h"
14  
15    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
17 >  logical :: setSim = .false.
18  
19 < !! WARNING: use_pbc hardcoded, fixme
23 <   logical :: setSim = .false.
19 >  integer,public :: natoms
20  
25 !! array for saving previous positions for neighbor lists.  
26  real( kind = dp ), allocatable,dimension(:,:),save :: q0
27
28
29  public :: wrap
21    public :: getBox
22    public :: getRcut
23    public :: getRlist
24 +  public :: getRrf
25 +  public :: getRt
26    public :: getNlocal
27    public :: setSimulation
28    public :: isEnsemble
29    public :: isPBC
30    public :: getStringLen
31    public :: returnMixingRules
32 +  public :: doStress
33  
34   !  public :: setRcut
41
42  interface wrap
43     module procedure wrap_1d
44     module procedure wrap_3d
45  end interface
46
35    interface getBox
36       module procedure getBox_3d
37       module procedure getBox_dim
38    end interface
39  
52
53
54
55
56
57
58
59
60
40   contains
41  
42 <  subroutine setSimulation(nLRParticles,box,rlist,rcut,ensemble,mixingRule,use_pbc)
43 <    integer, intent(in) :: nLRParticles
44 <    real(kind = dp ), intent(in), dimension(3) :: box
66 <    real(kind = dp ), intent(in) :: rlist
67 <    real(kind = dp ), intent(in) :: rcut
68 <    character( len = stringLen), intent(in)  :: ensemble
69 <    character( len = stringLen), intent(in)  :: mixingRule
70 <    logical, intent(in) :: use_pbc
42 >  subroutine setSimulation(setThisSim,error)
43 >    type (simtype) :: setThisSim
44 >    integer :: error
45      integer :: alloc_stat
72    if( setsim ) return  ! simulation is already initialized
73    setSim = .true.
46  
47 <    thisSim%nLRParticles = nLRParticles
48 <    thisSim%box          = box
49 <    thisSim%rlist        = rlist
50 <    thisSIm%rlistsq      = rlist * rlist
79 <    thisSim%rcut         = rcut
80 <    thisSim%rcutsq       = rcut * rcut
81 <    thisSim%rcut6        = thisSim%rcutsq * thisSim%rcutsq * thisSim%rcutsq
82 <    
83 <    thisSim%ensemble = ensemble
84 <    thisSim%mixingRule = mixingRule
85 <    thisSim%use_pbc = use_pbc
86 <
87 <    if (.not. allocated(q0)) then
88 <       allocate(q0(3,nLRParticles),stat=alloc_stat)
89 <    endif
47 >    error = 0
48 >    setSim = .true.
49 >    ! copy C struct into fortran type
50 >    thisSim = setThisSim
51    end subroutine setSimulation
52  
53    function getNparticles() result(nparticles)
# Line 94 | Line 55 | contains
55      nparticles = thisSim%nLRparticles
56    end function getNparticles
57  
97
58    subroutine change_box_size(new_box_size)
59      real(kind=dp), dimension(3) :: new_box_size
100
60      thisSim%box = new_box_size
102
61    end subroutine change_box_size
62  
105
63    function getBox_3d() result(thisBox)
64      real( kind = dp ), dimension(3) :: thisBox
65      thisBox = thisSim%box
# Line 111 | Line 68 | contains
68    function getBox_dim(dim) result(thisBox)
69      integer, intent(in) :: dim
70      real( kind = dp ) :: thisBox
71 <
71 >    
72      thisBox = thisSim%box(dim)
73    end function getBox_dim
74 <  
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 <    if (this_sim%use_pbc) then
142 <       !     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 <
74 >    
75    subroutine getRcut(thisrcut,rcut2,rcut6,status)
76      real( kind = dp ), intent(out) :: thisrcut
77      real( kind = dp ), intent(out), optional :: rcut2
# Line 164 | Line 88 | contains
88      thisrcut = thisSim%rcut
89      if(present(rcut2)) rcut2 = thisSim%rcutsq
90      if(present(rcut6)) rcut6 = thisSim%rcut6
167
91    end subroutine getRcut
92    
170  
171  
172
93    subroutine getRlist(thisrlist,rlist2,status)
94      real( kind = dp ), intent(out) :: thisrlist
95      real( kind = dp ), intent(out), optional :: rlist2
# Line 185 | Line 105 | contains
105      
106      thisrlist = thisSim%rlist
107      if(present(rlist2)) rlist2 = thisSim%rlistsq
188
189
108    end subroutine getRlist
191  
192  
109  
110 < pure function getNlocal() result(nlocal)
110 >  function getRrf() result(rrf)
111 >    real( kind = dp ) :: rrf
112 >    rrf = thisSim%rrf
113 >  end function getRrf
114 >  
115 >  function getRt() result(rt)
116 >    real( kind = dp ) :: rt
117 >    rt = thisSim%rt
118 >  end function getRt
119 >  
120 >  pure function getNlocal() result(nlocal)
121      integer :: nlocal
122      nlocal = thisSim%nLRparticles
123    end function getNlocal
124 <
125 <
124 >  
125 >  function doStress() result(do_stress)
126 >    logical :: do_stress
127 >    do_stress = thisSim%do_stress
128 >  end function doStress
129 >  
130    function isEnsemble(this_ensemble) result(is_this_ensemble)
131      character(len = *) :: this_ensemble
132 <    logical :: is_this_enemble
132 >    logical :: is_this_ensemble
133      is_this_ensemble = .false.
134      if (this_ensemble == thisSim%ensemble) is_this_ensemble = .true.
135    end function isEnsemble
136 <
136 >  
137    function returnEnsemble() result(thisEnsemble)
138      character (len = len(thisSim%ensemble)) :: thisEnsemble
139      thisEnsemble = thisSim%ensemble
140    end function returnEnsemble
141 <
141 >  
142    function returnMixingRules() result(thisMixingRule)
143      character (len = len(thisSim%ensemble)) :: thisMixingRule
144      thisMixingRule = thisSim%MixingRule
145    end function returnMixingRules
146 <
146 >  
147    function isPBC() result(PBCset)
148      logical :: PBCset
149      PBCset = .false.
150      if (thisSim%use_pbc) PBCset = .true.
151    end function isPBC
152 <
152 >  
153    pure function getStringLen() result (thislen)
154      integer :: thislen    
155      thislen = string_len
156 <  end function setStringLen
157 <
156 >  end function getStringLen
157 >  
158   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines