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

Comparing trunk/OOPSE_old/src/mdtools/libmdCode/calc_reaction_field.F90 (file contents):
Revision 360 by gezelter, Tue Mar 18 16:46:47 2003 UTC vs.
Revision 361 by gezelter, Tue Mar 18 19:09:12 2003 UTC

# Line 36 | Line 36 | contains
36  
37      integer, intent(in) :: atom1, atom2
38      real (kind = dp), intent(in) :: rij
39 <    real (kind = dp), dimension(3, getNlocal()) :: u_l    
39 >    real (kind = dp), dimension(:,:) :: u_l    
40  
41      integer :: me1, me2
42      real (kind = dp) :: taper, mu1, mu2
# Line 108 | Line 108 | contains
108      
109      integer, intent(in) :: atom1
110      real(kind=dp), intent(in) :: mu1
111 <    real(kind=dp), dimension(3,getNlocal()) :: u_l
111 >    real(kind=dp), dimension(:,:) :: u_l
112      
113      !! should work for both MPI and non-MPI version since this is not pairwise.
114      rf(1,atom1) = rf(1,atom1) + u_l(1,atom1)*mu1
# Line 124 | Line 124 | contains
124      real (kind=dp), intent(in) :: mu1
125      real (kind=dp), intent(inout) :: rfpot
126      logical, intent(in) :: do_pot
127 <    real (kind = dp), dimension(3, getNlocal()) :: u_l    
128 <    real (kind = dp), dimension(3, getNlocal()) :: t
127 >    real (kind = dp), dimension(:,:) :: u_l    
128 >    real (kind = dp), dimension(:,:) :: t
129  
130      if (.not.rf_initialized) then
131         write(default_error,*) 'Reaction field not initialized!'
# Line 156 | Line 156 | contains
156      integer, intent(in) :: atom1, atom2
157      real(kind=dp), dimension(3), intent(in) :: d
158      real(kind=dp), intent(in) :: rij
159 <    real( kind = dp ), dimension(3,getNlocal()) :: u_l
160 <    real( kind = dp ), dimension(3,getNlocal()) :: f
159 >    real( kind = dp ), dimension(:,:) :: u_l
160 >    real( kind = dp ), dimension(:,:) :: f
161      logical, intent(in) :: do_stress
162      
163      real (kind = dp), dimension(3) :: ul1

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines