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

Comparing trunk/OOPSE/libmdtools/simulation_module.F90 (file contents):
Revision 378 by mmeineke, Fri Mar 21 17:42:12 2003 UTC vs.
Revision 747 by gezelter, Fri Sep 5 21:28:52 2003 UTC

# Line 6 | Line 6 | module simulation
6    use force_globals
7    use vector_class
8    use atype_module
9  use lj
9   #ifdef IS_MPI
10    use mpiSimulation
11   #endif
# Line 17 | Line 16 | module simulation
16   #define __FORTRAN90
17   #include "fSimulation.h"
18  
19 <  type (simtype), public :: thisSim
19 >  type (simtype), public, save :: thisSim
20  
21    logical, save :: simulation_setup_complete = .false.
22  
23 <  integer, public, save :: natoms
23 >  integer, public, save :: nLocal, nGlobal
24    integer, public, save :: nExcludes_Global = 0
25    integer, public, save :: nExcludes_Local = 0
26    integer, allocatable, dimension(:,:), public :: excludesLocal
27    integer, allocatable, dimension(:), public :: excludesGlobal
28 +  integer, allocatable, dimension(:), public :: molMembershipList
29  
30 <  real(kind=dp), save :: rcut2 = 0.0_DP
31 <  real(kind=dp), save :: rcut6 = 0.0_DP
32 <  real(kind=dp), save :: rlist2 = 0.0_DP
33 <  real(kind=dp), public, dimension(3), save :: box
34 <
30 >  real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
31 >  logical, public, save :: boxIsOrthorhombic
32    
33    public :: SimulationSetup
34    public :: getNlocal
35    public :: setBox
39  public :: setBox_3d
40  public :: getBox
41  public :: setRcut
42  public :: getRcut
43  public :: getRlist
44  public :: getRrf
45  public :: getRt
36    public :: getDielect
37    public :: SimUsesPBC
38    public :: SimUsesLJ
# Line 54 | Line 44 | module simulation
44    public :: SimRequiresPrepairCalc
45    public :: SimRequiresPostpairCalc
46    public :: SimUsesDirectionalAtoms
57
58  interface getBox
59     module procedure getBox_3d
60     module procedure getBox_1d
61  end interface
47    
63  interface setBox
64     module procedure setBox_3d
65     module procedure setBox_1d
66  end interface
67  
48   contains
49    
50 <  subroutine SimulationSetup(setThisSim, nComponents, c_idents, &
50 >  subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
51         CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, &
52 +       CmolMembership, &
53         status)    
54  
55      type (simtype) :: setThisSim
56 <    integer, intent(inout) :: nComponents
57 <    integer, dimension(nComponents),intent(inout) :: c_idents
56 >    integer, intent(inout) :: CnGlobal, CnLocal
57 >    integer, dimension(CnLocal),intent(inout) :: c_idents
58  
59      integer :: CnLocalExcludes
60      integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal
61      integer :: CnGlobalExcludes
62      integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal
63 +    integer, dimension(CnGlobal),intent(in) :: CmolMembership
64      !!  Result status, success = 0, status = -1
65      integer, intent(out) :: status
66      integer :: i, me, thisStat, alloc_stat, myNode
# Line 93 | Line 75 | contains
75      status = 0
76  
77      ! copy C struct into fortran type
78 +
79 +    nLocal = CnLocal
80 +    nGlobal = CnGlobal
81 +
82      thisSim = setThisSim
97    natoms = nComponents
98    rcut2 = thisSim%rcut * thisSim%rcut
99    rcut6 = rcut2 * rcut2 * rcut2
100    rlist2 = thisSim%rlist * thisSim%rlist
101    box = thisSim%box
83  
84      nExcludes_Global = CnGlobalExcludes
85      nExcludes_Local = CnLocalExcludes
86  
87 <    call InitializeForceGlobals(natoms, thisStat)
87 >    call InitializeForceGlobals(nLocal, thisStat)
88      if (thisStat /= 0) then
89 +       write(default_error,*) "SimSetup: InitializeForceGlobals error"
90         status = -1
91         return
92      endif
93  
94      call InitializeSimGlobals(thisStat)
95      if (thisStat /= 0) then
96 +       write(default_error,*) "SimSetup: InitializeSimGlobals error"
97         status = -1
98         return
99      endif
# Line 159 | Line 142 | contains
142         deallocate(c_idents_Row)
143      endif
144      
145 < #else
146 <    do i = 1, nComponents
145 > #endif
146 >
147 > ! We build the local atid's for both mpi and nonmpi
148 >    do i = 1, nLocal
149        
150         me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
151         atid(i) = me
152    
153      enddo
169 #endif
154  
171    !! Create neighbor lists
172    call expandNeighborList(nComponents, thisStat)
173    if (thisStat /= 0) then
174       status = -1
175       return
176    endif
155  
156 +
157 +
158      do i = 1, nExcludes_Local
159         excludesLocal(1,i) = CexcludesLocal(1,i)
160         excludesLocal(2,i) = CexcludesLocal(2,i)
# Line 183 | Line 163 | contains
163      do i = 1, nExcludes_Global
164         excludesGlobal(i) = CexcludesGlobal(i)
165      enddo
166 <    
166 >
167 >    do i = 1, nGlobal
168 >       molMemberShipList(i) = CmolMembership(i)
169 >     enddo
170 >
171      if (status == 0) simulation_setup_complete = .true.
172      
173    end subroutine SimulationSetup
174    
175 <  subroutine setBox_3d(new_box_size)
176 <    real(kind=dp), dimension(3) :: new_box_size
175 >  subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
176 >    real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
177 >    integer :: cBoxIsOrthorhombic
178      integer :: smallest, status, i
179 <
180 <    thisSim%box = new_box_size
181 <    box = thisSim%box
182 <
183 <    smallest = 1
184 <    do i = 2, 3
185 <       if (new_box_size(i) .lt. new_box_size(smallest)) smallest = i
186 <    end do
187 <    if (thisSim%rcut .gt. 0.5_dp * new_box_size(smallest)) &
203 <         call setRcut(0.5_dp * new_box_size(smallest), status)
179 >    
180 >    Hmat = cHmat
181 >    HmatInv = cHmatInv
182 >    if (cBoxIsOrthorhombic .eq. 0 ) then
183 >       boxIsOrthorhombic = .false.
184 >    else
185 >       boxIsOrthorhombic = .true.
186 >    endif
187 >    
188      return    
189 <  end subroutine setBox_3d
206 <
207 <  subroutine setBox_1d(dim, new_box_size)
208 <    integer :: dim, status
209 <    real(kind=dp) :: new_box_size
210 <    thisSim%box(dim) = new_box_size
211 <    box(dim) = thisSim%box(dim)
212 <    if (thisSim%rcut .gt. 0.5_dp * new_box_size) &
213 <         call setRcut(0.5_dp * new_box_size, status)
214 <  end subroutine setBox_1d
189 >  end subroutine setBox
190  
216  subroutine setRcut(new_rcut, status)
217    real(kind = dp) :: new_rcut
218    integer :: myStatus, status
219    thisSim%rcut = new_rcut
220    rcut2 = thisSim%rcut * thisSim%rcut
221    rcut6 = rcut2 * rcut2 * rcut2
222    myStatus = 0
223    call LJ_new_rcut(new_rcut, myStatus)
224    if (myStatus .ne. 0) then
225       write(default_error, *) 'LJ module refused our rcut!'
226       status = -1
227       return
228    endif
229    status = 0
230    return
231  end subroutine setRcut
232
233  function getBox_3d() result(thisBox)
234    real( kind = dp ), dimension(3) :: thisBox
235    thisBox = thisSim%box
236  end function getBox_3d
237  
238  function getBox_1d(dim) result(thisBox)
239    integer, intent(in) :: dim
240    real( kind = dp ) :: thisBox
241    
242    thisBox = thisSim%box(dim)
243  end function getBox_1d
244    
245  subroutine getRcut(thisrcut,rc2,rc6,status)
246    real( kind = dp ), intent(out) :: thisrcut
247    real( kind = dp ), intent(out), optional :: rc2
248    real( kind = dp ), intent(out), optional :: rc6
249    integer, optional :: status
250
251    if (present(status)) status = 0
252    
253    if (.not.simulation_setup_complete ) then
254       if (present(status)) status = -1
255       return
256    end if
257    
258    thisrcut = thisSim%rcut
259    if(present(rc2)) rc2 = rcut2
260    if(present(rc6)) rc6 = rcut6
261  end subroutine getRcut
262  
263  subroutine getRlist(thisrlist,rl2,status)
264    real( kind = dp ), intent(out) :: thisrlist
265    real( kind = dp ), intent(out), optional :: rl2
266
267    integer, optional :: status
268
269    if (present(status)) status = 0
270
271    if (.not.simulation_setup_complete ) then
272       if (present(status)) status = -1
273       return
274    end if
275    
276    thisrlist = thisSim%rlist
277    if(present(rl2)) rl2 = rlist2
278  end subroutine getRlist
279
280  function getRrf() result(rrf)
281    real( kind = dp ) :: rrf
282    rrf = thisSim%rrf
283  end function getRrf
284  
285  function getRt() result(rt)
286    real( kind = dp ) :: rt
287    rt = thisSim%rt
288  end function getRt
289
191    function getDielect() result(dielect)
192      real( kind = dp ) :: dielect
193      dielect = thisSim%dielect
# Line 362 | Line 263 | contains
263         thisStat = -1
264         return
265      endif
266 +
267 +    allocate(molMembershipList(nGlobal), stat=alloc_stat)
268 +    if (alloc_stat /= 0 ) then
269 +       thisStat = -1
270 +       return
271 +    endif
272      
273    end subroutine InitializeSimGlobals
274    
275    subroutine FreeSimGlobals()
276      
277      !We free in the opposite order in which we allocate in.
278 <    
278 >
279 >    if (allocated(molMembershipList)) deallocate(molMembershipList)    
280      if (allocated(excludesGlobal)) deallocate(excludesGlobal)
281      if (allocated(excludesLocal)) deallocate(excludesLocal)
282  
283    end subroutine FreeSimGlobals
284  
285 <  pure function getNlocal() result(nlocal)
286 <    integer :: nlocal
287 <    nlocal = natoms
285 >  pure function getNlocal() result(n)
286 >    integer :: n
287 >    n = nLocal
288    end function getNlocal
289  
290    

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines