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 221 by chuckv, Thu Dec 19 21:59:51 2002 UTC vs.
Revision 222 by chuckv, Thu Jan 2 21:45:45 2003 UTC

# Line 6 | Line 6 | module simulation
6    PRIVATE
7  
8  
9 <  
10 <  real ( kind = dp ), public, dimension(3) :: box
9 >  integer :: nLR
10  
11  
12 +  real ( kind = dp ), dimension(3) :: box
13 +
14 +
15    real ( kind = dp ), public :: rlist
16    real ( kind = dp ), public :: rcut
17    real ( kind = dp ), public :: rlistsq
# Line 19 | Line 21 | module simulation
21    integer,public, allocatable, dimension(:)   :: list
22  
23  
24 +
25 +
26 +
27 +  public :: check
28 +  public :: save_nlist
29 +  public :: wrap
30 +  public :: getBox
31 +
32 +  interface wrap
33 +     module procedure wrap_1d
34 +     module procedure wrap_3d
35 +  end interface
36 +
37 +  interface getBox
38 +     module procedure getBox_3d
39 +     module procedure getBox_dim
40 +  end interface
41 +
42 +
43   !MPI dependent routines
44  
45   #ifdef IS_MPI
# Line 70 | Line 91 | contains
91  
92    end subroutine change_box_size
93  
94 +
95 +  elemental function getBox_3d() result(thisBox)
96 +    real( kind = dp ), dimension(3) :: thisBox
97 +    thisBox = box
98 +  end function getBox_3d
99 +
100 +  function getBox_dim(dim) result(thisBox)
101 +    integer, intent(in) :: dim
102 +    real( kind = dp ) :: thisBox
103 +
104 +    thisBox = box(dim)
105 +  end function getBox_dim
106 +    
107 +  subroutine check(update_nlist)
108    
109 +    integer :: i
110 +    real( kind = DP ) :: dispmx
111 +    logical, intent(out) :: update_nlist
112 +    real( kind = DP ) :: dispmx_tmp
113 +    
114 +    dispmx = 0.0E0_DP
115 +    
116 +    !! calculate the largest displacement of any atom in any direction
117 +    
118 + #ifdef MPI
119 +    dispmx_tmp = 0.0E0_DP
120 +    do i = 1, nlocal
121 +       dispmx_tmp = max( abs ( q(1,i) - q0(1,i) ), dispmx )
122 +       dispmx_tmp = max( abs ( q(2,i) - q0(2,i) ), dispmx )
123 +       dispmx_tmp = max( abs ( q(3,i) - q0(3,i) ), dispmx )
124 +    end do
125 +    call mpi_allreduce(dispmx_tmp,dispmx,1,mpi_double_precision, &
126 +       mpi_max,mpi_comm_world,mpi_err)
127 + #else
128 +    do i = 1, natoms
129 +       dispmx = max( abs ( q(1,i) - q0(1,i) ), dispmx )
130 +       dispmx = max( abs ( q(2,i) - q0(2,i) ), dispmx )
131 +       dispmx = max( abs ( q(3,i) - q0(3,i) ), dispmx )
132 +    end do
133 + #endif
134 +  
135 +    !! a conservative test of list skin crossings
136 +    dispmx = 2.0E0_DP * sqrt (3.0E0_DP * dispmx * dispmx)
137 +    
138 +    update_nlist = (dispmx.gt.(skin_thickness))
139 +    
140 +  end subroutine check
141 +  
142 +  subroutine save_nlist()
143 +    integer :: i
144 + #ifdef MPI
145 +    do i = 1, nlocal
146 + #else
147 +       do i = 1, natoms
148 + #endif
149 +          q0(1,i) = q(1,i)
150 +          q0(2,i) = q(2,i)
151 +          q0(3,i) = q(3,i)
152 +       end do
153 +
154 +  end subroutine save_nlist
155 +  
156 +
157 +  function wrap_1d(r,dim) result(this_wrap)
158 +    
159 +    
160 +    real( kind = DP ) :: r
161 +    real( kind = DP ) :: this_wrap
162 +    integer           :: dim
163 +    
164 +    if (use_pbc) then
165 +       !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
166 +       this_wrap = r - box(dim)*nint(r/box(dim))
167 +    else
168 +       this_wrap = r
169 +    endif
170 +    
171 +    return
172 +  end function wrap_1d
173 +
174 +  elemental function wrap_3d(r) result(this_wrap)
175 +    real( kind = dp ), dimension(3), intent(in) :: r
176 +    real( kind = dp ), dimension(3) :: this_wrap
177 +
178 +    
179 +    if (use_pbc) then
180 +       !     this_wrap = r - box(dim)*dsign(1.0E0_DP,r)*int(abs(r/box(dim)) + 0.5E0_DP)
181 +       this_wrap(1:3) = r(1:3) - box(1:3)*nint(r(1:3)/box(1:3))
182 +    else
183 +       this_wrap = r
184 +    endif
185 +  end function wrap_3d
186 +
187 +  
188   end module simulation

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines