83
84
87
88
89
90 integer, intent(in) :: ng, tile
91 integer, intent(in) :: LBi, UBi, LBj, UBj
92 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
93 integer, intent(in) :: nrhs, nnew
94
95# ifdef ASSUMED_SHAPE
96# ifdef MASKING
97 real(r8), intent(in) :: umask(LBi:,LBj:)
98 real(r8), intent(in) :: vmask(LBi:,LBj:)
99# endif
100 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
101 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
102 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
103 real(r8), intent(in) :: pm(LBi:,LBj:)
104 real(r8), intent(in) :: pn(LBi:,LBj:)
105 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
106
107 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
108# else
109# ifdef MASKING
110 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
112# endif
113 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
114 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
115 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
116 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
117 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
118 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
119
120 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
121# endif
122
123
124
125 integer :: i, itrc, j, k
126 real(r8) :: adfac
127
128 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FE
129 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FX
130
131# include "set_bounds.h"
132
133
134
135
136
137
138
139 ad_fx=0.0_r8
140 ad_fe=0.0_r8
141
143 IF (
tl_tdiff(itrc,ng).gt.0.0_r8)
THEN
145
146
147
148 DO j=jstr,jend
149 DO i=istr,iend
150
151
152
153
154
155 adfac=
dt(ng)*pm(i,j)*pn(i,j)*ad_t(i,j,k,nnew,itrc)
156 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac
157 ad_fx(i,j)=ad_fx(i,j)-adfac
158 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac
159 ad_fe(i,j)=ad_fe(i,j)+adfac
160 END DO
161 END DO
162
163
164
165
166 DO j=jstr,jend+1
167 DO i=istr,iend
168# ifdef MASKING
169
170
171 ad_fe(i,j)=ad_fe(i,j)*vmask(i,j)
172# endif
173
174
175
176
177
178 adfac=0.5_r8*
tl_tdiff(itrc,ng)*pnom_v(i,j)* &
179 & (hz(i,j,k)+hz(i,j-1,k))*ad_fe(i,j)
180 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac
181 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac
182 ad_fe(i,j)=0.0_r8
183 END DO
184 END DO
185
186 DO j=jstr,jend
187 DO i=istr,iend+1
188# ifdef MASKING
189
190
191 ad_fx(i,j)=ad_fx(i,j)*umask(i,j)
192# endif
193
194
195
196
197
198 adfac=0.5_r8*
tl_tdiff(itrc,ng)*pmon_u(i,j)* &
199 & (hz(i,j,k)+hz(i-1,j,k))*ad_fx(i,j)
200 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac
201 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac
202 ad_fx(i,j)=0.0_r8
203 END DO
204 END DO
205 END DO
206 END IF
207 END DO
208
209 RETURN
integer, dimension(:), allocatable n
integer, dimension(:), allocatable nt
real(r8), dimension(:,:), allocatable tl_tdiff
real(dp), dimension(:), allocatable dt